PolyFEM
Loading...
Searching...
No Matches
ElasticityUtils.cpp
Go to the documentation of this file.
1#include "ElasticityUtils.hpp"
3
4#include <tinyexpr.h>
5
6namespace polyfem
7{
8 using namespace assembler;
9 using namespace basis;
10 using namespace utils;
11
12 double convert_to_lambda(const bool is_volume, const double E, const double nu)
13 {
14 if (is_volume)
15 return (E * nu) / ((1.0 + nu) * (1.0 - 2.0 * nu));
16
17 return (nu * E) / (1.0 - nu * nu);
18 }
19
20 double convert_to_mu(const double E, const double nu)
21 {
22 return E / (2.0 * (1.0 + nu));
23 }
24
25 Eigen::Matrix2d d_lambda_mu_d_E_nu(const bool is_volume, const double E, const double nu)
26 {
27 Eigen::Matrix2d A;
28 if (is_volume)
29 {
30 A(0, 0) = nu / ((1.0 + nu) * (1.0 - 2.0 * nu));
31 A(0, 1) = (E*(0.5*nu*nu + 0.25))/(pow(nu + 1., 2)*pow(0.5 - nu, 2));
32 }
33 else
34 {
35 A(0, 0) = nu / (1.0 - nu * nu);
36 A(0, 1) = (E*(1. + nu*nu))/pow(-1. + nu*nu, 2);
37 }
38 A(1, 0) = 1 / (2 * (1 + nu));
39 A(1, 1) = -E / 2 * pow(1 + nu, -2);
40
41 return A;
42 }
43
44 Eigen::Matrix2d d_E_nu_d_lambda_mu(const bool is_volume, const double lambda, const double mu)
45 {
46 Eigen::Matrix2d A;
47 if (is_volume)
48 {
49 A(0, 0) = pow(mu / (lambda + mu), 2);
50 A(0, 1) = (3 * lambda * lambda + 4 * lambda * mu + 2 * mu * mu) / pow(lambda + mu, 2);
51 A(1, 0) = mu / 2 / pow(lambda + mu, 2);
52 A(1, 1) = -lambda / 2 / pow(lambda + mu, 2);
53 }
54 else
55 {
56 A(0, 0) = pow(2 * mu / (lambda + 2 * mu), 2);
57 A(0, 1) = (4 * lambda * lambda + 8 * lambda * mu + 8 * mu * mu) / pow(lambda + 2 * mu, 2);
58 A(1, 0) = mu * 2 / pow(lambda + 2 * mu, 2);
59 A(1, 1) = -lambda * 2 / pow(lambda + 2 * mu, 2);
60 }
61
62 return A;
63 }
64
65 double convert_to_E(const bool is_volume, const double lambda, const double mu)
66 {
67 if (is_volume)
68 return mu * (3.0 * lambda + 2.0 * mu) / (lambda + mu);
69
70 return 2 * mu * (2.0 * lambda + 2.0 * mu) / (lambda + 2.0 * mu);
71 }
72
73 double convert_to_nu(const bool is_volume, const double lambda, const double mu)
74 {
75 if (is_volume)
76 return lambda / (2.0 * (lambda + mu));
77
78 return lambda / (lambda + 2.0 * mu);
79 }
80
81 Eigen::VectorXd gradient_from_energy(const int size, const int n_bases, const assembler::NonLinearAssemblerData &data,
82 const std::function<DScalar1<double, Eigen::Matrix<double, 6, 1>>(const assembler::NonLinearAssemblerData &)> &fun6,
83 const std::function<DScalar1<double, Eigen::Matrix<double, 8, 1>>(const assembler::NonLinearAssemblerData &)> &fun8,
84 const std::function<DScalar1<double, Eigen::Matrix<double, 12, 1>>(const assembler::NonLinearAssemblerData &)> &fun12,
85 const std::function<DScalar1<double, Eigen::Matrix<double, 18, 1>>(const assembler::NonLinearAssemblerData &)> &fun18,
86 const std::function<DScalar1<double, Eigen::Matrix<double, 24, 1>>(const assembler::NonLinearAssemblerData &)> &fun24,
87 const std::function<DScalar1<double, Eigen::Matrix<double, 30, 1>>(const assembler::NonLinearAssemblerData &)> &fun30,
88 const std::function<DScalar1<double, Eigen::Matrix<double, 60, 1>>(const assembler::NonLinearAssemblerData &)> &fun60,
89 const std::function<DScalar1<double, Eigen::Matrix<double, 81, 1>>(const assembler::NonLinearAssemblerData &)> &fun81,
90 const std::function<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, SMALL_N, 1>>(const assembler::NonLinearAssemblerData &)> &funN,
91 const std::function<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, BIG_N, 1>>(const assembler::NonLinearAssemblerData &)> &funBigN,
93 {
94 Eigen::VectorXd grad;
95
96 switch (size * n_bases)
97 {
98 case 6:
99 {
100 auto auto_diff_energy = fun6(data);
101 grad = auto_diff_energy.getGradient();
102 break;
103 }
104 case 8:
105 {
106 auto auto_diff_energy = fun8(data);
107 grad = auto_diff_energy.getGradient();
108 break;
109 }
110 case 12:
111 {
112 auto auto_diff_energy = fun12(data);
113 grad = auto_diff_energy.getGradient();
114 break;
115 }
116 case 18:
117 {
118 auto auto_diff_energy = fun18(data);
119 grad = auto_diff_energy.getGradient();
120 break;
121 }
122 case 24:
123 {
124 auto auto_diff_energy = fun24(data);
125 grad = auto_diff_energy.getGradient();
126 break;
127 }
128 case 30:
129 {
130 auto auto_diff_energy = fun30(data);
131 grad = auto_diff_energy.getGradient();
132 break;
133 }
134 case 60:
135 {
136 auto auto_diff_energy = fun60(data);
137 grad = auto_diff_energy.getGradient();
138 break;
139 }
140 case 81:
141 {
142 auto auto_diff_energy = fun81(data);
143 grad = auto_diff_energy.getGradient();
144 break;
145 }
146 default: // default handled after with grad size
147 break;
148 }
149
150 if (grad.size() <= 0)
151 {
152 if (n_bases * size <= SMALL_N)
153 {
154 auto auto_diff_energy = funN(data);
155 grad = auto_diff_energy.getGradient();
156 }
157 else if (n_bases * size <= BIG_N)
158 {
159 auto auto_diff_energy = funBigN(data);
160 grad = auto_diff_energy.getGradient();
161 }
162 else
163 {
164 static bool show_message = true;
165
166 if (show_message)
167 {
168 logger().debug("[Warning] grad {}^{} not using static sizes", n_bases, size);
169 show_message = false;
170 }
171
172 auto auto_diff_energy = funn(data);
173 grad = auto_diff_energy.getGradient();
174 }
175 }
176
177 return grad;
178 }
179
180 Eigen::MatrixXd hessian_from_energy(const int size, const int n_bases, const assembler::NonLinearAssemblerData &data,
181 const std::function<DScalar2<double, Eigen::Matrix<double, 6, 1>, Eigen::Matrix<double, 6, 6>>(const assembler::NonLinearAssemblerData &)> &fun6,
182 const std::function<DScalar2<double, Eigen::Matrix<double, 8, 1>, Eigen::Matrix<double, 8, 8>>(const assembler::NonLinearAssemblerData &)> &fun8,
183 const std::function<DScalar2<double, Eigen::Matrix<double, 12, 1>, Eigen::Matrix<double, 12, 12>>(const assembler::NonLinearAssemblerData &)> &fun12,
184 const std::function<DScalar2<double, Eigen::Matrix<double, 18, 1>, Eigen::Matrix<double, 18, 18>>(const assembler::NonLinearAssemblerData &)> &fun18,
185 const std::function<DScalar2<double, Eigen::Matrix<double, 24, 1>, Eigen::Matrix<double, 24, 24>>(const assembler::NonLinearAssemblerData &)> &fun24,
186 const std::function<DScalar2<double, Eigen::Matrix<double, 30, 1>, Eigen::Matrix<double, 30, 30>>(const assembler::NonLinearAssemblerData &)> &fun30,
187 const std::function<DScalar2<double, Eigen::Matrix<double, 60, 1>, Eigen::Matrix<double, 60, 60>>(const assembler::NonLinearAssemblerData &)> &fun60,
188 const std::function<DScalar2<double, Eigen::Matrix<double, 81, 1>, Eigen::Matrix<double, 81, 81>>(const assembler::NonLinearAssemblerData &)> &fun81,
189 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,
191 {
192 Eigen::MatrixXd hessian;
193
194 switch (size * n_bases)
195 {
196 case 6:
197 {
198 auto auto_diff_energy = fun6(data);
199 hessian = auto_diff_energy.getHessian();
200 break;
201 }
202 case 8:
203 {
204 auto auto_diff_energy = fun8(data);
205 hessian = auto_diff_energy.getHessian();
206 break;
207 }
208 case 12:
209 {
210 auto auto_diff_energy = fun12(data);
211 hessian = auto_diff_energy.getHessian();
212 break;
213 }
214 case 18:
215 {
216 auto auto_diff_energy = fun18(data);
217 hessian = auto_diff_energy.getHessian();
218 break;
219 }
220 case 24:
221 {
222 auto auto_diff_energy = fun24(data);
223 hessian = auto_diff_energy.getHessian();
224 break;
225 }
226 case 30:
227 {
228 auto auto_diff_energy = fun30(data);
229 hessian = auto_diff_energy.getHessian();
230 break;
231 }
232 case 60:
233 {
234 auto auto_diff_energy = fun60(data);
235 hessian = auto_diff_energy.getHessian();
236 break;
237 }
238 case 81:
239 {
240 auto auto_diff_energy = fun81(data);
241 hessian = auto_diff_energy.getHessian();
242 break;
243 }
244 default: // default handled after with hessian size
245 break;
246 }
247
248 if (hessian.size() <= 0)
249 {
250 if (n_bases * size <= SMALL_N)
251 {
252 auto auto_diff_energy = funN(data);
253 hessian = auto_diff_energy.getHessian();
254 }
255 else
256 {
257 static bool show_message = true;
258
259 if (show_message)
260 {
261 logger().debug("[Warning] hessian {}*{} not using static sizes", n_bases, size);
262 show_message = false;
263 }
264
265 auto auto_diff_energy = funn(data);
266 hessian = auto_diff_energy.getHessian();
267 }
268 }
269
270 // time.stop();
271 // std::cout << "-- hessian: " << time.getElapsedTime() << std::endl;
272
273 return hessian;
274 }
275
276 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)
277 {
278 assert(displacement.cols() == 1);
279
280 displacement_grad.setZero();
281
282 for (std::size_t j = 0; j < vals.basis_values.size(); ++j)
283 {
284 const auto &loc_val = vals.basis_values[j];
285
286 assert(loc_val.grad.rows() == local_pts.rows());
287 assert(loc_val.grad.cols() == size);
288
289 for (int d = 0; d < size; ++d)
290 {
291 for (std::size_t ii = 0; ii < loc_val.global.size(); ++ii)
292 {
293 displacement_grad.row(d) += loc_val.global[ii].val * loc_val.grad.row(p) * displacement(loc_val.global[ii].index * size + d);
294 }
295 }
296 }
297
298 displacement_grad = (displacement_grad * vals.jac_it[p]).eval();
299 }
300
301 void compute_diplacement_grad(const int size, const ElementBases &bs, const ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const int p, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &displacement_grad)
302 {
303 assert(displacement.cols() == 1);
304
305 displacement_grad.setZero();
306
307 for (std::size_t j = 0; j < bs.bases.size(); ++j)
308 {
309 const Basis &b = bs.bases[j];
310 const auto &loc_val = vals.basis_values[j];
311
312 assert(bs.bases.size() == vals.basis_values.size());
313 assert(loc_val.grad.rows() == local_pts.rows());
314 assert(loc_val.grad.cols() == size);
315
316 for (int d = 0; d < size; ++d)
317 {
318 for (std::size_t ii = 0; ii < b.global().size(); ++ii)
319 {
320 displacement_grad.row(d) += b.global()[ii].val * loc_val.grad.row(p) * displacement(b.global()[ii].index * size + d);
321 }
322 }
323 }
324
325 displacement_grad = (displacement_grad * vals.jac_it[p]).eval();
326 }
327
328 double von_mises_stress_for_stress_tensor(const Eigen::MatrixXd &stress)
329 {
330 double von_mises_stress;
331
332 if (stress.rows() == 3)
333 {
334 von_mises_stress = 0.5 * (stress(0, 0) - stress(1, 1)) * (stress(0, 0) - stress(1, 1)) + 3.0 * stress(0, 1) * stress(1, 0);
335 von_mises_stress += 0.5 * (stress(2, 2) - stress(1, 1)) * (stress(2, 2) - stress(1, 1)) + 3.0 * stress(2, 1) * stress(2, 1);
336 von_mises_stress += 0.5 * (stress(2, 2) - stress(0, 0)) * (stress(2, 2) - stress(0, 0)) + 3.0 * stress(2, 0) * stress(2, 0);
337 }
338 else
339 {
340 // von_mises_stress = ( stress(0, 0) - stress(1, 1) ) * ( stress(0, 0) - stress(1, 1) ) + 3.0 * stress(0, 1) * stress(1, 0);
341 von_mises_stress = stress(0, 0) * stress(0, 0) - stress(0, 0) * stress(1, 1) + stress(1, 1) * stress(1, 1) + 3.0 * stress(0, 1) * stress(1, 0);
342 }
343
344 von_mises_stress = sqrt(fabs(von_mises_stress));
345
346 return von_mises_stress;
347 }
348
349 Eigen::MatrixXd pk1_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
350 {
351 return F.determinant() * stress * F.inverse().transpose();
352 }
353
354 Eigen::MatrixXd pk2_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
355 {
356 return F.determinant() * F.inverse() * stress * F.inverse().transpose();
357 }
358} // namespace polyfem
ElementAssemblyValues vals
Definition Assembler.cpp:22
stores per element basis values at given quadrature points and geometric mapping
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
Represents one basis function and its gradient.
Definition Basis.hpp:44
const std::vector< Local2Global > & global() const
Definition Basis.hpp:104
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
std::vector< Basis > bases
one basis function per node in the element
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
constexpr int SMALL_N
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)
double convert_to_E(const bool is_volume, const double lambda, const double mu)
constexpr int BIG_N
double convert_to_mu(const double E, const double nu)
double von_mises_stress_for_stress_tensor(const Eigen::MatrixXd &stress)
double convert_to_nu(const bool is_volume, const double lambda, const double mu)
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)
double convert_to_lambda(const bool is_volume, const double E, const double nu)
Eigen::Matrix2d d_E_nu_d_lambda_mu(const bool is_volume, const double lambda, const double mu)
Eigen::Matrix2d d_lambda_mu_d_E_nu(const bool is_volume, const double E, const double nu)
Automatic differentiation scalar with first-order derivatives.
Definition autodiff.h:112
Automatic differentiation scalar with first- and second-order derivatives.
Definition autodiff.h:493