PolyFEM
Loading...
Searching...
No Matches
GenericElastic.cpp
Go to the documentation of this file.
1#include "GenericElastic.hpp"
2
10
12
14
15namespace polyfem::assembler
16{
17 template <typename Derived>
21
22 template <typename Derived>
24 const OutputData &data,
25 const int all_size,
27 Eigen::MatrixXd &all,
28 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const
29 {
30 Eigen::MatrixXd deformation_grad(size(), size());
31 Eigen::MatrixXd stress_tensor(size(), size());
32
35 const auto &displacement = data.fun;
36 const auto &local_pts = data.local_pts;
37 const auto &bs = data.bs;
38 const auto &gbs = data.gbs;
39 const auto el_id = data.el_id;
40
41 assert(displacement.cols() == 1);
42
43 all.resize(local_pts.rows(), all_size);
44 DiffScalarBase::setVariableCount(deformation_grad.size());
45
46 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
47
49 vals.compute(el_id, size() == 3, local_pts, bs, gbs);
50
51 for (long p = 0; p < local_pts.rows(); ++p)
52 {
53 compute_diplacement_grad(size(), bs, vals, local_pts, p, displacement, deformation_grad);
54
55 // Id + grad d
56 for (int d = 0; d < size(); ++d)
57 deformation_grad(d, d) += 1;
58
59 if (type == ElasticityTensorType::F)
60 {
61 all.row(p) = fun(deformation_grad);
62 continue;
63 }
64
65 for (int d1 = 0; d1 < size(); ++d1)
66 {
67 for (int d2 = 0; d2 < size(); ++d2)
68 def_grad(d1, d2) = Diff(d1 * size() + d2, deformation_grad(d1, d2));
69 }
70
71 const auto val = derived().elastic_energy(local_pts.row(p), data.t, vals.element_id, def_grad);
72
73 for (int d1 = 0; d1 < size(); ++d1)
74 {
75 for (int d2 = 0; d2 < size(); ++d2)
76 stress_tensor(d1, d2) = val.getGradient()(d1 * size() + d2);
77 }
78
79 stress_tensor = 1.0 / deformation_grad.determinant() * stress_tensor * deformation_grad.transpose();
80
81 if (type == ElasticityTensorType::PK1)
82 stress_tensor = pk1_from_cauchy(stress_tensor, deformation_grad);
83 else if (type == ElasticityTensorType::PK2)
84 stress_tensor = pk2_from_cauchy(stress_tensor, deformation_grad);
85
86 all.row(p) = fun(stress_tensor);
87 }
88 }
89
90 template <typename Derived>
92 {
93 return compute_energy_aux<double>(data);
94 }
95
96 template <typename Derived>
98 {
99 const int n_bases = data.vals.basis_values.size();
101 size(), n_bases, data,
102 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 6, 1>>>(data); },
103 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 8, 1>>>(data); },
104 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 12, 1>>>(data); },
105 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 18, 1>>>(data); },
106 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 24, 1>>>(data); },
107 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 30, 1>>>(data); },
108 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 60, 1>>>(data); },
109 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 81, 1>>>(data); },
110 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, SMALL_N, 1>>>(data); },
111 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, BIG_N, 1>>>(data); },
112 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::VectorXd>>(data); });
113 }
114
115 template <typename Derived>
117 {
118 const int n_bases = data.vals.basis_values.size();
120 size(), n_bases, data,
121 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 6, 1>, Eigen::Matrix<double, 6, 6>>>(data); },
122 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 8, 1>, Eigen::Matrix<double, 8, 8>>>(data); },
123 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 12, 1>, Eigen::Matrix<double, 12, 12>>>(data); },
124 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 18, 1>, Eigen::Matrix<double, 18, 18>>>(data); },
125 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 24, 1>, Eigen::Matrix<double, 24, 24>>>(data); },
126 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 30, 1>, Eigen::Matrix<double, 30, 30>>>(data); },
127 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 60, 1>, Eigen::Matrix<double, 60, 60>>>(data); },
128 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 81, 1>, Eigen::Matrix<double, 81, 81>>>(data); },
129 [&](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); },
130 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::VectorXd, Eigen::MatrixXd>>(data); });
131 }
132
133 template <typename Derived>
135 const OptAssemblerData &data,
136 const Eigen::MatrixXd &mat,
137 Eigen::MatrixXd &stress,
138 Eigen::MatrixXd &result) const
139 {
140 typedef DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9>> Diff;
141
142 const double t = data.t;
143 const int el_id = data.el_id;
144 const Eigen::MatrixXd &local_pts = data.local_pts;
145 const Eigen::MatrixXd &global_pts = data.global_pts;
146 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
147
148 DiffScalarBase::setVariableCount(size() * size());
149 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
150
151 Eigen::MatrixXd F = grad_u_i;
152 for (int d = 0; d < size(); ++d)
153 F(d, d) += 1.;
154
155 assert(local_pts.rows() == 1);
156 for (int i = 0; i < size(); ++i)
157 for (int j = 0; j < size(); ++j)
158 def_grad(i, j) = Diff(i + j * size(), F(i, j));
159
160 auto energy = derived().elastic_energy(global_pts, t, el_id, def_grad);
161
162 // Grad is ∂W(F)/∂F_ij
163 Eigen::MatrixXd grad = energy.getGradient().reshaped(size(), size());
164 // Hessian is ∂W(F)/(∂F_ij*∂F_kl)
165 Eigen::MatrixXd hess = energy.getHessian();
166
167 // Stress is S_ij = ∂W(F)/∂F_ij
168 stress = grad;
169 // Compute ∂S_ij/∂F_kl * M_kl, same as M_ij * ∂S_ij/∂F_kl since the hessian is symmetric
170 result = (hess * mat.reshaped(size() * size(), 1)).reshaped(size(), size());
171 }
172
173 template <typename Derived>
175 const OptAssemblerData &data,
176 Eigen::MatrixXd &stress,
177 Eigen::MatrixXd &result) const
178 {
179 typedef DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9>> Diff;
180
181 const double t = data.t;
182 const int el_id = data.el_id;
183 const Eigen::MatrixXd &local_pts = data.local_pts;
184 const Eigen::MatrixXd &global_pts = data.global_pts;
185 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
186
187 DiffScalarBase::setVariableCount(size() * size());
188 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
189
190 Eigen::MatrixXd F = grad_u_i;
191 for (int d = 0; d < size(); ++d)
192 F(d, d) += 1.;
193
194 assert(local_pts.rows() == 1);
195 for (int i = 0; i < size(); ++i)
196 for (int j = 0; j < size(); ++j)
197 def_grad(i, j) = Diff(i + j * size(), F(i, j));
198
199 auto energy = derived().elastic_energy(global_pts, t, el_id, def_grad);
200
201 // Grad is ∂W(F)/∂F_ij
202 Eigen::MatrixXd grad = energy.getGradient().reshaped(size(), size());
203 // Hessian is ∂W(F)/(∂F_ij*∂F_kl)
204 Eigen::MatrixXd hess = energy.getHessian();
205
206 // Stress is S_ij = ∂W(F)/∂F_ij
207 stress = grad;
208 // Compute ∂S_ij/∂F_kl * S_kl, same as S_ij * ∂S_ij/∂F_kl since the hessian is symmetric
209 result = (hess * stress.reshaped(size() * size(), 1)).reshaped(size(), size());
210 }
211
212 template <typename Derived>
214 const OptAssemblerData &data,
215 const Eigen::MatrixXd &vect,
216 Eigen::MatrixXd &stress,
217 Eigen::MatrixXd &result) const
218 {
219 typedef DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9>> Diff;
220
221 const double t = data.t;
222 const int el_id = data.el_id;
223 const Eigen::MatrixXd &local_pts = data.local_pts;
224 const Eigen::MatrixXd &global_pts = data.global_pts;
225 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
226
227 DiffScalarBase::setVariableCount(size() * size());
228 Eigen::Matrix<Diff, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> def_grad(size(), size());
229
230 Eigen::MatrixXd F = grad_u_i;
231 for (int d = 0; d < size(); ++d)
232 F(d, d) += 1.;
233
234 assert(local_pts.rows() == 1);
235 for (int i = 0; i < size(); ++i)
236 for (int j = 0; j < size(); ++j)
237 def_grad(i, j) = Diff(i + j * size(), F(i, j));
238
239 auto energy = derived().elastic_energy(global_pts, t, el_id, def_grad);
240
241 // Grad is ∂W(F)/∂F_ij
242 Eigen::MatrixXd grad = energy.getGradient().reshaped(size(), size());
243 // Hessian is ∂W(F)/(∂F_ij*∂F_kl)
244 Eigen::MatrixXd hess = energy.getHessian();
245
246 // Stress is S_ij = ∂W(F)/∂F_ij
247 stress = grad;
248 result.resize(hess.rows(), vect.size());
249 for (int i = 0; i < hess.rows(); ++i)
250 if (vect.rows() == 1)
251 // Compute ∂S_ij/∂F_kl * v_k, same as ∂S_ij/∂F_kl * v_i since the hessian is symmetric
252 result.row(i) = vect * hess.row(i).reshaped(size(), size());
253 else
254 // Compute ∂S_ij/∂F_kl * v_l, same as ∂S_ij/∂F_kl * v_j since the hessian is symmetric
255 result.row(i) = hess.row(i).reshaped(size(), size()) * vect;
256 }
264 template class GenericElastic<VolumePenalty>;
265
266 template class GenericElastic<HGOFiber>;
267
268} // 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