PolyFEM
Loading...
Searching...
No Matches
GenericElastic.cpp
Go to the documentation of this file.
1#include "GenericElastic.hpp"
2
8
10
11namespace polyfem::assembler
12{
13 template <typename Derived>
17
18 template <typename Derived>
20 const OutputData &data,
21 const int all_size,
22 const ElasticityTensorType &type,
23 Eigen::MatrixXd &all,
24 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const
25 {
26 Eigen::MatrixXd deformation_grad(size(), size());
27 Eigen::MatrixXd stress_tensor(size(), size());
30
31 const auto &displacement = data.fun;
32 const auto &local_pts = data.local_pts;
33 const auto &bs = data.bs;
34 const auto &gbs = data.gbs;
35 const auto el_id = data.el_id;
36
37 assert(displacement.cols() == 1);
38
39 all.resize(local_pts.rows(), all_size);
40 DiffScalarBase::setVariableCount(deformation_grad.size());
41
42 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
45 vals.compute(el_id, size() == 3, local_pts, bs, gbs);
46
47 for (long p = 0; p < local_pts.rows(); ++p)
48 {
49 compute_diplacement_grad(size(), bs, vals, local_pts, p, displacement, deformation_grad);
50
51 // Id + grad d
52 for (int d = 0; d < size(); ++d)
53 deformation_grad(d, d) += 1;
54
55 if (type == ElasticityTensorType::F)
56 {
57 all.row(p) = fun(deformation_grad);
58 continue;
59 }
60
61 for (int d1 = 0; d1 < size(); ++d1)
62 {
63 for (int d2 = 0; d2 < size(); ++d2)
64 def_grad(d1, d2) = Diff(d1 * size() + d2, deformation_grad(d1, d2));
65 }
66
67 const auto val = derived().elastic_energy(local_pts.row(p), data.t, vals.element_id, def_grad);
68
69 for (int d1 = 0; d1 < size(); ++d1)
70 {
71 for (int d2 = 0; d2 < size(); ++d2)
72 stress_tensor(d1, d2) = val.getGradient()(d1 * size() + d2);
73 }
74
75 stress_tensor = 1.0 / deformation_grad.determinant() * stress_tensor * deformation_grad.transpose();
76
77 if (type == ElasticityTensorType::PK1)
78 stress_tensor = pk1_from_cauchy(stress_tensor, deformation_grad);
79 else if (type == ElasticityTensorType::PK2)
80 stress_tensor = pk2_from_cauchy(stress_tensor, deformation_grad);
81
82 all.row(p) = fun(stress_tensor);
83 }
84 }
85
86 template <typename Derived>
88 {
89 return compute_energy_aux<double>(data);
90 }
91
92 template <typename Derived>
94 {
95 const int n_bases = data.vals.basis_values.size();
97 size(), n_bases, data,
98 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 6, 1>>>(data); },
99 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 8, 1>>>(data); },
100 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 12, 1>>>(data); },
101 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 18, 1>>>(data); },
102 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 24, 1>>>(data); },
103 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 30, 1>>>(data); },
104 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 60, 1>>>(data); },
105 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 81, 1>>>(data); },
106 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, SMALL_N, 1>>>(data); },
107 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, BIG_N, 1>>>(data); },
108 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::VectorXd>>(data); });
109 }
110
111 template <typename Derived>
113 {
114 const int n_bases = data.vals.basis_values.size();
116 size(), n_bases, data,
117 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 6, 1>, Eigen::Matrix<double, 6, 6>>>(data); },
118 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 8, 1>, Eigen::Matrix<double, 8, 8>>>(data); },
119 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 12, 1>, Eigen::Matrix<double, 12, 12>>>(data); },
120 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 18, 1>, Eigen::Matrix<double, 18, 18>>>(data); },
121 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 24, 1>, Eigen::Matrix<double, 24, 24>>>(data); },
122 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 30, 1>, Eigen::Matrix<double, 30, 30>>>(data); },
123 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 60, 1>, Eigen::Matrix<double, 60, 60>>>(data); },
124 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 81, 1>, Eigen::Matrix<double, 81, 81>>>(data); },
125 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, SMALL_N, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, SMALL_N, SMALL_N>>>(data); },
126 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::VectorXd, Eigen::MatrixXd>>(data); });
127 }
128
129 template <typename Derived>
131 const OptAssemblerData &data,
132 const Eigen::MatrixXd &mat,
133 Eigen::MatrixXd &stress,
134 Eigen::MatrixXd &result) const
135 {
136 typedef DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9>> Diff;
137
138 const double t = data.t;
139 const int el_id = data.el_id;
140 const Eigen::MatrixXd &local_pts = data.local_pts;
141 const Eigen::MatrixXd &global_pts = data.global_pts;
142 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
143
144 DiffScalarBase::setVariableCount(size() * size());
145 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
146
147 Eigen::MatrixXd F = grad_u_i;
148 for (int d = 0; d < size(); ++d)
149 F(d, d) += 1.;
150
151 assert(local_pts.rows() == 1);
152 for (int i = 0; i < size(); ++i)
153 for (int j = 0; j < size(); ++j)
154 def_grad(i, j) = Diff(i + j * size(), F(i, j));
155
156 auto energy = derived().elastic_energy(global_pts, t, el_id, def_grad);
157
158 // Grad is ∂W(F)/∂F_ij
159 Eigen::MatrixXd grad = energy.getGradient().reshaped(size(), size());
160 // Hessian is ∂W(F)/(∂F_ij*∂F_kl)
161 Eigen::MatrixXd hess = energy.getHessian();
162
163 // Stress is S_ij = ∂W(F)/∂F_ij
164 stress = grad;
165 // Compute ∂S_ij/∂F_kl * M_kl, same as M_ij * ∂S_ij/∂F_kl since the hessian is symmetric
166 result = (hess * mat.reshaped(size() * size(), 1)).reshaped(size(), size());
167 }
168
169 template <typename Derived>
171 const OptAssemblerData &data,
172 Eigen::MatrixXd &stress,
173 Eigen::MatrixXd &result) const
174 {
175 typedef DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9>> Diff;
176
177 const double t = data.t;
178 const int el_id = data.el_id;
179 const Eigen::MatrixXd &local_pts = data.local_pts;
180 const Eigen::MatrixXd &global_pts = data.global_pts;
181 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
182
183 DiffScalarBase::setVariableCount(size() * size());
184 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
185
186 Eigen::MatrixXd F = grad_u_i;
187 for (int d = 0; d < size(); ++d)
188 F(d, d) += 1.;
189
190 assert(local_pts.rows() == 1);
191 for (int i = 0; i < size(); ++i)
192 for (int j = 0; j < size(); ++j)
193 def_grad(i, j) = Diff(i + j * size(), F(i, j));
194
195 auto energy = derived().elastic_energy(global_pts, t, el_id, def_grad);
196
197 // Grad is ∂W(F)/∂F_ij
198 Eigen::MatrixXd grad = energy.getGradient().reshaped(size(), size());
199 // Hessian is ∂W(F)/(∂F_ij*∂F_kl)
200 Eigen::MatrixXd hess = energy.getHessian();
201
202 // Stress is S_ij = ∂W(F)/∂F_ij
203 stress = grad;
204 // Compute ∂S_ij/∂F_kl * S_kl, same as S_ij * ∂S_ij/∂F_kl since the hessian is symmetric
205 result = (hess * stress.reshaped(size() * size(), 1)).reshaped(size(), size());
206 }
207
208 template <typename Derived>
210 const OptAssemblerData &data,
211 const Eigen::MatrixXd &vect,
212 Eigen::MatrixXd &stress,
213 Eigen::MatrixXd &result) const
214 {
215 typedef DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9>> Diff;
216
217 const double t = data.t;
218 const int el_id = data.el_id;
219 const Eigen::MatrixXd &local_pts = data.local_pts;
220 const Eigen::MatrixXd &global_pts = data.global_pts;
221 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
222
223 DiffScalarBase::setVariableCount(size() * size());
224 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
225
226 Eigen::MatrixXd F = grad_u_i;
227 for (int d = 0; d < size(); ++d)
228 F(d, d) += 1.;
229
230 assert(local_pts.rows() == 1);
231 for (int i = 0; i < size(); ++i)
232 for (int j = 0; j < size(); ++j)
233 def_grad(i, j) = Diff(i + j * size(), F(i, j));
234
235 auto energy = derived().elastic_energy(global_pts, t, el_id, def_grad);
236
237 // Grad is ∂W(F)/∂F_ij
238 Eigen::MatrixXd grad = energy.getGradient().reshaped(size(), size());
239 // Hessian is ∂W(F)/(∂F_ij*∂F_kl)
240 Eigen::MatrixXd hess = energy.getHessian();
241
242 // Stress is S_ij = ∂W(F)/∂F_ij
243 stress = grad;
244 result.resize(hess.rows(), vect.size());
245 for (int i = 0; i < hess.rows(); ++i)
246 if (vect.rows() == 1)
247 // Compute ∂S_ij/∂F_kl * v_k, same as ∂S_ij/∂F_kl * v_i since the hessian is symmetric
248 result.row(i) = vect * hess.row(i).reshaped(size(), size());
249 else
250 // Compute ∂S_ij/∂F_kl * v_l, same as ∂S_ij/∂F_kl * v_j since the hessian is symmetric
251 result.row(i) = hess.row(i).reshaped(size(), size()) * vect;
252 }
259} // namespace polyfem::assembler
double val
Definition Assembler.cpp:86
ElementAssemblyValues vals
Definition Assembler.cpp:22
stores per element basis values at given quadrature points and geometric mapping
void compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
computes the per element values at the local (ref el) points (pts) sets basis_values,...
void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void assign_stress_tensor(const OutputData &data, const int all_size, const ElasticityTensorType &type, Eigen::MatrixXd &all, const std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
void compute_stress_grad_multiply_vect(const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
double compute_energy(const NonLinearAssemblerData &data) const override
void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
const basis::ElementBases & bs
const Eigen::MatrixXd & fun
const Eigen::MatrixXd & local_pts
const basis::ElementBases & gbs
Used for test only.
DScalar1< double, Eigen::Matrix< double, Eigen::Dynamic, 1 > > Diff
Eigen::MatrixXd pk2_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
Eigen::MatrixXd hessian_from_energy(const int size, const int n_bases, const assembler::NonLinearAssemblerData &data, const std::function< DScalar2< double, Eigen::Matrix< double, 6, 1 >, Eigen::Matrix< double, 6, 6 > >(const assembler::NonLinearAssemblerData &)> &fun6, const std::function< DScalar2< double, Eigen::Matrix< double, 8, 1 >, Eigen::Matrix< double, 8, 8 > >(const assembler::NonLinearAssemblerData &)> &fun8, const std::function< DScalar2< double, Eigen::Matrix< double, 12, 1 >, Eigen::Matrix< double, 12, 12 > >(const assembler::NonLinearAssemblerData &)> &fun12, const std::function< DScalar2< double, Eigen::Matrix< double, 18, 1 >, Eigen::Matrix< double, 18, 18 > >(const assembler::NonLinearAssemblerData &)> &fun18, const std::function< DScalar2< double, Eigen::Matrix< double, 24, 1 >, Eigen::Matrix< double, 24, 24 > >(const assembler::NonLinearAssemblerData &)> &fun24, const std::function< DScalar2< double, Eigen::Matrix< double, 30, 1 >, Eigen::Matrix< double, 30, 30 > >(const assembler::NonLinearAssemblerData &)> &fun30, const std::function< DScalar2< double, Eigen::Matrix< double, 60, 1 >, Eigen::Matrix< double, 60, 60 > >(const assembler::NonLinearAssemblerData &)> &fun60, const std::function< DScalar2< double, Eigen::Matrix< double, 81, 1 >, Eigen::Matrix< double, 81, 81 > >(const assembler::NonLinearAssemblerData &)> &fun81, const std::function< DScalar2< double, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, SMALL_N, 1 >, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, SMALL_N, SMALL_N > >(const assembler::NonLinearAssemblerData &)> &funN, const std::function< DScalar2< double, Eigen::VectorXd, Eigen::MatrixXd >(const assembler::NonLinearAssemblerData &)> &funn)
void compute_diplacement_grad(const int size, const ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const int p, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &displacement_grad)
Eigen::MatrixXd pk1_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
Eigen::VectorXd gradient_from_energy(const int size, const int n_bases, const assembler::NonLinearAssemblerData &data, const std::function< DScalar1< double, Eigen::Matrix< double, 6, 1 > >(const assembler::NonLinearAssemblerData &)> &fun6, const std::function< DScalar1< double, Eigen::Matrix< double, 8, 1 > >(const assembler::NonLinearAssemblerData &)> &fun8, const std::function< DScalar1< double, Eigen::Matrix< double, 12, 1 > >(const assembler::NonLinearAssemblerData &)> &fun12, const std::function< DScalar1< double, Eigen::Matrix< double, 18, 1 > >(const assembler::NonLinearAssemblerData &)> &fun18, const std::function< DScalar1< double, Eigen::Matrix< double, 24, 1 > >(const assembler::NonLinearAssemblerData &)> &fun24, const std::function< DScalar1< double, Eigen::Matrix< double, 30, 1 > >(const assembler::NonLinearAssemblerData &)> &fun30, const std::function< DScalar1< double, Eigen::Matrix< double, 60, 1 > >(const assembler::NonLinearAssemblerData &)> &fun60, const std::function< DScalar1< double, Eigen::Matrix< double, 81, 1 > >(const assembler::NonLinearAssemblerData &)> &fun81, const std::function< DScalar1< double, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, SMALL_N, 1 > >(const assembler::NonLinearAssemblerData &)> &funN, const std::function< DScalar1< double, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, BIG_N, 1 > >(const assembler::NonLinearAssemblerData &)> &funBigN, const std::function< DScalar1< double, Eigen::VectorXd >(const assembler::NonLinearAssemblerData &)> &funn)
Automatic differentiation scalar with first-order derivatives.
Definition autodiff.h:112
Automatic differentiation scalar with first- and second-order derivatives.
Definition autodiff.h:493
static void setVariableCount(size_t value)
Set the independent variable count used by the automatic differentiation layer.
Definition autodiff.h:54