PolyFEM
Loading...
Searching...
No Matches
NeoHookeanElasticity.cpp
Go to the documentation of this file.
2
5
6#include <type_traits>
7
8namespace polyfem::assembler
9{
10 namespace
11 {
12 bool delta(int i, int j)
13 {
14 return (i == j) ? true : false;
15 }
16 } // namespace
17
21
22 void NeoHookeanElasticity::add_multimaterial(const int index, const json &params, const Units &units)
23 {
24 assert(size() == 2 || size() == 3);
25
26 params_.add_multimaterial(index, params, size() == 3, units.stress());
27 }
28
29 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1>
31 {
32 assert(pt.size() == size());
33 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1> res;
34
35 double lambda, mu;
36 // TODO!
37 params_.lambda_mu(0, 0, 0, pt(0).getValue(), pt(1).getValue(), size() == 2 ? 0. : pt(2).getValue(), 0, 0, lambda, mu);
38
39 if (size() == 2)
40 autogen::neo_hookean_2d_function(pt, lambda, mu, res);
41 else if (size() == 3)
42 autogen::neo_hookean_3d_function(pt, lambda, mu, res);
43 else
44 assert(false);
45
46 return res;
47 }
48
49 Eigen::VectorXd
51 {
52 Eigen::Matrix<double, Eigen::Dynamic, 1> gradient;
53
54 if (size() == 2)
55 {
56 switch (data.vals.basis_values.size())
57 {
58 case 3:
59 {
60 gradient.resize(6);
61 compute_energy_aux_gradient_fast<3, 2>(data, gradient);
62 break;
63 }
64 case 6:
65 {
66 gradient.resize(12);
67 compute_energy_aux_gradient_fast<6, 2>(data, gradient);
68 break;
69 }
70 case 10:
71 {
72 gradient.resize(20);
73 compute_energy_aux_gradient_fast<10, 2>(data, gradient);
74 break;
75 }
76 default:
77 {
78 gradient.resize(data.vals.basis_values.size() * 2);
79 compute_energy_aux_gradient_fast<Eigen::Dynamic, 2>(data, gradient);
80 break;
81 }
82 }
83 }
84 else // if (size() == 3)
85 {
86 assert(size() == 3);
87 switch (data.vals.basis_values.size())
88 {
89 case 4:
90 {
91 gradient.resize(12);
92 compute_energy_aux_gradient_fast<4, 3>(data, gradient);
93 break;
94 }
95 case 10:
96 {
97 gradient.resize(30);
98 compute_energy_aux_gradient_fast<10, 3>(data, gradient);
99 break;
100 }
101 case 20:
102 {
103 gradient.resize(60);
104 compute_energy_aux_gradient_fast<20, 3>(data, gradient);
105 break;
106 }
107 default:
108 {
109 gradient.resize(data.vals.basis_values.size() * 3);
110 compute_energy_aux_gradient_fast<Eigen::Dynamic, 3>(data, gradient);
111 break;
112 }
113 }
114 }
115
116 return gradient;
117 }
118
121 const Eigen::MatrixXd &local_pts,
122 const Eigen::MatrixXd &displacement,
123 Eigen::MatrixXd &tensor) const
124 {
125 tensor.resize(local_pts.rows(), size() * size() * size() * size());
126 assert(displacement.cols() == 1);
127
128 Eigen::MatrixXd displacement_grad(size(), size());
129
130 for (long p = 0; p < local_pts.rows(); ++p)
131 {
132 double lambda, mu;
133 params_.lambda_mu(local_pts.row(p), vals.val.row(p), t, vals.element_id, lambda, mu);
134
135 compute_diplacement_grad(size(), vals, local_pts, p, displacement, displacement_grad);
136 const Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(size(), size()) + displacement_grad;
137 const Eigen::MatrixXd FmT = def_grad.inverse().transpose();
138 const Eigen::VectorXd FmT_vec = utils::flatten(FmT);
139 const double J = def_grad.determinant();
140 const double tmp1 = mu - lambda * std::log(J);
141 for (int i = 0, idx = 0; i < size(); i++)
142 for (int j = 0; j < size(); j++)
143 for (int k = 0; k < size(); k++)
144 for (int l = 0; l < size(); l++)
145 {
146 tensor(p, idx) = mu * delta(i, k) * delta(j, l) + tmp1 * FmT(i, l) * FmT(k, j);
147 idx++;
148 }
149
150 tensor.row(p) += lambda * utils::flatten(FmT_vec * FmT_vec.transpose());
151
152 // {
153 // Eigen::MatrixXd hess = utils::unflatten(tensor.row(p), size()*size());
154 // Eigen::MatrixXd fhess;
155 // Eigen::VectorXd x0 = utils::flatten(def_grad);
156 // fd::finite_jacobian(
157 // x0, [this, lambda, mu](const Eigen::VectorXd &x1) -> Eigen::VectorXd
158 // {
159 // Eigen::MatrixXd def_grad = utils::unflatten(x1, this->size());
160 // const Eigen::MatrixXd FmT = def_grad.inverse().transpose();
161 // const double J = def_grad.determinant();
162 // Eigen::MatrixXd stress_tensor = mu * (def_grad - FmT) + lambda * std::log(J) * FmT;
163 // return utils::flatten(stress_tensor);
164 // }, fhess);
165
166 // if (!fd::compare_hessian(hess, fhess))
167 // {
168 // std::cout << "Hessian: " << hess << std::endl;
169 // std::cout << "Finite hessian: " << fhess << std::endl;
170 // log_and_throw_error("Hessian in Neohookean mismatch");
171 // }
172 // }
173 }
174 }
175
176 Eigen::MatrixXd
178 {
179 Eigen::MatrixXd hessian;
180
181 if (size() == 2)
182 {
183 switch (data.vals.basis_values.size())
184 {
185 case 3:
186 {
187 hessian.resize(6, 6);
188 hessian.setZero();
189 compute_energy_hessian_aux_fast<3, 2>(data, hessian);
190 break;
191 }
192 case 6:
193 {
194 hessian.resize(12, 12);
195 hessian.setZero();
196 compute_energy_hessian_aux_fast<6, 2>(data, hessian);
197 break;
198 }
199 case 10:
200 {
201 hessian.resize(20, 20);
202 hessian.setZero();
203 compute_energy_hessian_aux_fast<10, 2>(data, hessian);
204 break;
205 }
206 default:
207 {
208 hessian.resize(data.vals.basis_values.size() * 2, data.vals.basis_values.size() * 2);
209 hessian.setZero();
210 compute_energy_hessian_aux_fast<Eigen::Dynamic, 2>(data, hessian);
211 break;
212 }
213 }
214 }
215 else // if (size() == 3)
216 {
217 assert(size() == 3);
218 switch (data.vals.basis_values.size())
219 {
220 case 4:
221 {
222 hessian.resize(12, 12);
223 hessian.setZero();
224 compute_energy_hessian_aux_fast<4, 3>(data, hessian);
225 break;
226 }
227 case 10:
228 {
229 hessian.resize(30, 30);
230 hessian.setZero();
231 compute_energy_hessian_aux_fast<10, 3>(data, hessian);
232 break;
233 }
234 case 20:
235 {
236 hessian.resize(60, 60);
237 hessian.setZero();
238 compute_energy_hessian_aux_fast<20, 3>(data, hessian);
239 break;
240 }
241 default:
242 {
243 hessian.resize(data.vals.basis_values.size() * 3, data.vals.basis_values.size() * 3);
244 hessian.setZero();
245 compute_energy_hessian_aux_fast<Eigen::Dynamic, 3>(data, hessian);
246 break;
247 }
248 }
249 }
250
251 return hessian;
252 }
253
255 const int all_size,
256 const ElasticityTensorType &type,
257 Eigen::MatrixXd &all,
258 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const
259 {
260 const auto &displacement = data.fun;
261 const auto &local_pts = data.local_pts;
262 const auto &bs = data.bs;
263 const auto &gbs = data.gbs;
264 const auto el_id = data.el_id;
265 const auto t = data.t;
266
267 Eigen::MatrixXd displacement_grad(size(), size());
268
269 assert(displacement.cols() == 1);
270
271 all.resize(local_pts.rows(), all_size);
272
274 vals.compute(el_id, size() == 3, local_pts, bs, gbs);
275 const auto I = Eigen::MatrixXd::Identity(size(), size());
276
277 for (long p = 0; p < local_pts.rows(); ++p)
278 {
279 compute_diplacement_grad(size(), bs, vals, local_pts, p, displacement, displacement_grad);
280
281 const Eigen::MatrixXd def_grad = I + displacement_grad;
282 if (type == ElasticityTensorType::F)
283 {
284 all.row(p) = fun(def_grad);
285 continue;
286 }
287 const double J = def_grad.determinant();
288 const Eigen::MatrixXd b = def_grad * def_grad.transpose();
289
290 double lambda, mu;
291 params_.lambda_mu(local_pts.row(p), vals.val.row(p), t, vals.element_id, lambda, mu);
292
293 Eigen::MatrixXd stress_tensor = (lambda * std::log(J) * I + mu * (b - I)) / J;
294 if (type == ElasticityTensorType::PK1)
295 stress_tensor = pk1_from_cauchy(stress_tensor, def_grad);
296 else if (type == ElasticityTensorType::PK2)
297 stress_tensor = pk2_from_cauchy(stress_tensor, def_grad);
298
299 all.row(p) = fun(stress_tensor);
300 }
301 }
302
304 {
305 return compute_energy_aux<double>(data);
306 }
307
308 // Compute ∫ ½μ (tr(FᵀF) - 3 - 2ln(J)) + ½λ ln²(J) du
309 template <typename T>
311 {
312 if constexpr (std::is_same_v<T, double>)
313 {
314 Eigen::MatrixXd local_disp(data.vals.basis_values.size(), size());
315 local_disp.setZero();
316 for (size_t i = 0; i < data.vals.basis_values.size(); ++i)
317 {
318 const auto &bs = data.vals.basis_values[i];
319 for (size_t ii = 0; ii < bs.global.size(); ++ii)
320 {
321 for (int d = 0; d < size(); ++d)
322 {
323 local_disp(i, d) += bs.global[ii].val * data.x(bs.global[ii].index * size() + d);
324 }
325 }
326 }
327 Eigen::VectorXd jacs;
329 jacs = data.vals.eval_deformed_jacobian_determinant(data.x);
330 Eigen::MatrixXd def_grad(size(), size());
331
332 T energy = T(0.0);
333
334 const int n_pts = data.da.size();
335 for (long p = 0; p < n_pts; ++p)
336 {
337 Eigen::MatrixXd grad(data.vals.basis_values.size(), size());
338
339 for (size_t i = 0; i < data.vals.basis_values.size(); ++i)
340 {
341 grad.row(i) = data.vals.basis_values[i].grad.row(p);
342 }
343
344 const Eigen::MatrixXd jac_it = data.vals.jac_it[p];
345 // Id + grad d
346 def_grad = local_disp.transpose() * grad * jac_it + Eigen::MatrixXd::Identity(size(), size());
347 double lambda, mu;
348 params_.lambda_mu(data.vals.quadrature.points.row(p), data.vals.val.row(p), data.t, data.vals.element_id, lambda, mu);
349 const T J = use_robust_jacobian ? jacs(p) * jac_it.determinant() : def_grad.determinant();
350 const T log_det_j = log(J);
351 const T val = mu / 2 * ((def_grad.transpose() * def_grad).trace() - size() - 2 * log_det_j) +
352 lambda / 2 * log_det_j * log_det_j;
353 energy += val * data.da(p);
354 }
355 return energy;
356 }
357 else
358 {
359 typedef Eigen::Matrix<T, Eigen::Dynamic, 1> AutoDiffVect;
360 typedef Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> AutoDiffGradMat;
361
362 AutoDiffVect local_disp;
363 get_local_disp(data, size(), local_disp);
364
365 AutoDiffGradMat def_grad(size(), size());
366
367 T energy = T(0.0);
368
369 const int n_pts = data.da.size();
370 for (long p = 0; p < n_pts; ++p)
371 {
372 compute_disp_grad_at_quad(data, local_disp, p, size(), def_grad);
373
374 // Id + grad d
375 for (int d = 0; d < size(); ++d)
376 def_grad(d, d) += T(1);
377
378 double lambda, mu;
379 params_.lambda_mu(data.vals.quadrature.points.row(p), data.vals.val.row(p), data.t, data.vals.element_id, lambda, mu);
380
381 const T log_det_j = log(polyfem::utils::determinant(def_grad));
382 const T val = mu / 2 * ((def_grad.transpose() * def_grad).trace() - size() - 2 * log_det_j) + lambda / 2 * log_det_j * log_det_j;
383
384 energy += val * data.da(p);
385 }
386 return energy;
387 }
388 }
389
390 template <int dim>
391 Eigen::Matrix<double, dim, dim> hat(const Eigen::Matrix<double, dim, 1> &x)
392 {
393
394 Eigen::Matrix<double, dim, dim> prod;
395 prod.setZero();
396
397 prod(0, 1) = -x(2);
398 prod(0, 2) = x(1);
399 prod(1, 0) = x(2);
400 prod(1, 2) = -x(0);
401 prod(2, 0) = -x(1);
402 prod(2, 1) = x(0);
403
404 return prod;
405 }
406
407 template <int dim>
408 Eigen::Matrix<double, dim, 1> cross(const Eigen::Matrix<double, dim, 1> &x, const Eigen::Matrix<double, dim, 1> &y)
409 {
410
411 Eigen::Matrix<double, dim, 1> z;
412 z.setZero();
413
414 z(0) = x(1) * y(2) - x(2) * y(1);
415 z(1) = x(2) * y(0) - x(0) * y(2);
416 z(2) = x(0) * y(1) - x(1) * y(0);
417
418 return z;
419 }
420
421 template <int n_basis, int dim>
422 void NeoHookeanElasticity::compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::Matrix<double, Eigen::Dynamic, 1> &G_flattened) const
423 {
424 assert(data.x.cols() == 1);
425
426 const int n_pts = data.da.size();
427
428 Eigen::Matrix<double, n_basis, dim> local_disp(data.vals.basis_values.size(), size());
429 local_disp.setZero();
430 for (size_t i = 0; i < data.vals.basis_values.size(); ++i)
431 {
432 const auto &bs = data.vals.basis_values[i];
433 for (size_t ii = 0; ii < bs.global.size(); ++ii)
434 {
435 for (int d = 0; d < size(); ++d)
436 {
437 local_disp(i, d) += bs.global[ii].val * data.x(bs.global[ii].index * size() + d);
438 }
439 }
440 }
441
442 Eigen::VectorXd jacs;
444 jacs = data.vals.eval_deformed_jacobian_determinant(data.x);
445
446 Eigen::Matrix<double, dim, dim> def_grad(size(), size());
447
448 Eigen::Matrix<double, n_basis, dim> G(data.vals.basis_values.size(), size());
449 G.setZero();
450
451 for (long p = 0; p < n_pts; ++p)
452 {
453 Eigen::Matrix<double, n_basis, dim> grad(data.vals.basis_values.size(), size());
454
455 for (size_t i = 0; i < data.vals.basis_values.size(); ++i)
456 {
457 grad.row(i) = data.vals.basis_values[i].grad.row(p);
458 }
459
460 Eigen::Matrix<double, dim, dim> jac_it = data.vals.jac_it[p];
461
462 // Id + grad d
463 def_grad = local_disp.transpose() * grad * jac_it + Eigen::Matrix<double, dim, dim>::Identity(size(), size());
464
465 const double J = use_robust_jacobian ? jacs(p) * jac_it.determinant() : def_grad.determinant();
466 const double log_det_j = log(J);
467
468 Eigen::Matrix<double, dim, dim> delJ_delF(size(), size());
469 delJ_delF.setZero();
470
471 if (dim == 2)
472 {
473
474 delJ_delF(0, 0) = def_grad(1, 1);
475 delJ_delF(0, 1) = -def_grad(1, 0);
476 delJ_delF(1, 0) = -def_grad(0, 1);
477 delJ_delF(1, 1) = def_grad(0, 0);
478 }
479
480 else if (dim == 3)
481 {
482
483 Eigen::Matrix<double, dim, 1> u(def_grad.rows());
484 Eigen::Matrix<double, dim, 1> v(def_grad.rows());
485 Eigen::Matrix<double, dim, 1> w(def_grad.rows());
486
487 u = def_grad.col(0);
488 v = def_grad.col(1);
489 w = def_grad.col(2);
490
491 delJ_delF.col(0) = cross<dim>(v, w);
492 delJ_delF.col(1) = cross<dim>(w, u);
493 delJ_delF.col(2) = cross<dim>(u, v);
494 }
495
496 double lambda, mu;
497 params_.lambda_mu(data.vals.quadrature.points.row(p), data.vals.val.row(p), data.t, data.vals.element_id, lambda, mu);
498
499 Eigen::Matrix<double, n_basis, dim> delF_delU = grad * jac_it;
500
501 Eigen::Matrix<double, dim, dim> gradient_temp = mu * def_grad - mu * (1 / J) * delJ_delF + lambda * log_det_j * (1 / J) * delJ_delF;
502 Eigen::Matrix<double, n_basis, dim> gradient = delF_delU * gradient_temp.transpose();
503
504 double val = mu / 2 * ((def_grad.transpose() * def_grad).trace() - size() - 2 * log_det_j) + lambda / 2 * log_det_j * log_det_j;
505
506 G.noalias() += gradient * data.da(p);
507 }
508
509 Eigen::Matrix<double, dim, n_basis> G_T = G.transpose();
510
511 constexpr int N = (n_basis == Eigen::Dynamic) ? Eigen::Dynamic : n_basis * dim;
512 Eigen::Matrix<double, N, 1> temp(Eigen::Map<Eigen::Matrix<double, N, 1>>(G_T.data(), G_T.size()));
513 G_flattened = temp;
514 }
515
516 template <int n_basis, int dim>
518 {
519 assert(data.x.cols() == 1);
520
521 constexpr int N = (n_basis == Eigen::Dynamic) ? Eigen::Dynamic : n_basis * dim;
522 const int n_pts = data.da.size();
523
524 Eigen::Matrix<double, n_basis, dim> local_disp(data.vals.basis_values.size(), size());
525 local_disp.setZero();
526 for (size_t i = 0; i < data.vals.basis_values.size(); ++i)
527 {
528 const auto &bs = data.vals.basis_values[i];
529 for (size_t ii = 0; ii < bs.global.size(); ++ii)
530 {
531 for (int d = 0; d < size(); ++d)
532 {
533 local_disp(i, d) += bs.global[ii].val * data.x(bs.global[ii].index * size() + d);
534 }
535 }
536 }
537
538 Eigen::VectorXd jacs;
540 jacs = data.vals.eval_deformed_jacobian_determinant(data.x);
541
542 Eigen::Matrix<double, dim, dim> def_grad(size(), size());
543
544 for (long p = 0; p < n_pts; ++p)
545 {
546 Eigen::Matrix<double, n_basis, dim> grad(data.vals.basis_values.size(), size());
547
548 for (size_t i = 0; i < data.vals.basis_values.size(); ++i)
549 {
550 grad.row(i) = data.vals.basis_values[i].grad.row(p);
551 }
552
553 Eigen::Matrix<double, dim, dim> jac_it = data.vals.jac_it[p];
554
555 // Id + grad d
556 def_grad = local_disp.transpose() * grad * jac_it + Eigen::Matrix<double, dim, dim>::Identity(size(), size());
557
558 const double J = use_robust_jacobian ? jacs(p) * jac_it.determinant() : def_grad.determinant();
559 double log_det_j = log(J);
560
561 Eigen::Matrix<double, dim, dim> delJ_delF(size(), size());
562 delJ_delF.setZero();
563 Eigen::Matrix<double, dim * dim, dim * dim> del2J_delF2(size() * size(), size() * size());
564 del2J_delF2.setZero();
565
566 if (dim == 2)
567 {
568 delJ_delF(0, 0) = def_grad(1, 1);
569 delJ_delF(0, 1) = -def_grad(1, 0);
570 delJ_delF(1, 0) = -def_grad(0, 1);
571 delJ_delF(1, 1) = def_grad(0, 0);
572
573 del2J_delF2(0, 3) = 1;
574 del2J_delF2(1, 2) = -1;
575 del2J_delF2(2, 1) = -1;
576 del2J_delF2(3, 0) = 1;
577 }
578 else if (size() == 3)
579 {
580 Eigen::Matrix<double, dim, 1> u(def_grad.rows());
581 Eigen::Matrix<double, dim, 1> v(def_grad.rows());
582 Eigen::Matrix<double, dim, 1> w(def_grad.rows());
583
584 u = def_grad.col(0);
585 v = def_grad.col(1);
586 w = def_grad.col(2);
587
588 delJ_delF.col(0) = cross<dim>(v, w);
589 delJ_delF.col(1) = cross<dim>(w, u);
590 delJ_delF.col(2) = cross<dim>(u, v);
591
592 del2J_delF2.template block<dim, dim>(0, 6) = hat<dim>(v);
593 del2J_delF2.template block<dim, dim>(6, 0) = -hat<dim>(v);
594 del2J_delF2.template block<dim, dim>(0, 3) = -hat<dim>(w);
595 del2J_delF2.template block<dim, dim>(3, 0) = hat<dim>(w);
596 del2J_delF2.template block<dim, dim>(3, 6) = -hat<dim>(u);
597 del2J_delF2.template block<dim, dim>(6, 3) = hat<dim>(u);
598 }
599
600 double lambda, mu;
601 params_.lambda_mu(data.vals.quadrature.points.row(p), data.vals.val.row(p), data.t, data.vals.element_id, lambda, mu);
602
603 Eigen::Matrix<double, dim * dim, dim * dim> id = Eigen::Matrix<double, dim * dim, dim * dim>::Identity(size() * size(), size() * size());
604
605 Eigen::Matrix<double, dim * dim, 1> g_j = Eigen::Map<const Eigen::Matrix<double, dim * dim, 1>>(delJ_delF.data(), delJ_delF.size());
606
607 Eigen::Matrix<double, dim * dim, dim * dim> hessian_temp = (mu * id) + (((mu + lambda * (1 - log_det_j)) / (J * J)) * (g_j * g_j.transpose())) + (((lambda * log_det_j - mu) / (J)) * del2J_delF2);
608
609 Eigen::Matrix<double, dim * dim, N> delF_delU_tensor(jac_it.size(), grad.size());
610
611 for (size_t i = 0; i < local_disp.rows(); ++i)
612 {
613 for (size_t j = 0; j < local_disp.cols(); ++j)
614 {
615 Eigen::Matrix<double, dim, dim> temp(size(), size());
616 temp.setZero();
617 temp.row(j) = grad.row(i);
618 temp = temp * jac_it;
619 Eigen::Matrix<double, dim * dim, 1> temp_flattened(Eigen::Map<Eigen::Matrix<double, dim * dim, 1>>(temp.data(), temp.size()));
620 delF_delU_tensor.col(i * size() + j) = temp_flattened;
621 }
622 }
623
624 Eigen::Matrix<double, N, N> hessian = delF_delU_tensor.transpose() * hessian_temp * delF_delU_tensor;
625
626 double val = mu / 2 * ((def_grad.transpose() * def_grad).trace() - size() - 2 * log_det_j) + lambda / 2 * log_det_j * log_det_j;
627
628 H += hessian * data.da(p);
629 }
630 }
631
633 const OptAssemblerData &data,
634 const Eigen::MatrixXd &mat,
635 Eigen::MatrixXd &stress,
636 Eigen::MatrixXd &result) const
637 {
638 const double t = data.t;
639 const int el_id = data.el_id;
640 const Eigen::MatrixXd &local_pts = data.local_pts;
641 const Eigen::MatrixXd &global_pts = data.global_pts;
642 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
643
644 double lambda, mu;
645 params_.lambda_mu(local_pts, global_pts, t, el_id, lambda, mu);
646
647 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
648 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
649
650 stress = mu * (def_grad - FmT) + lambda * std::log(def_grad.determinant()) * FmT;
651 result = mu * mat + FmT * mat.transpose() * FmT * (mu - lambda * std::log(def_grad.determinant())) + lambda * (FmT.array() * mat.array()).sum() * FmT;
652 }
653
655 const OptAssemblerData &data,
656 Eigen::MatrixXd &stress,
657 Eigen::MatrixXd &result) const
658 {
659 const double t = data.t;
660 const int el_id = data.el_id;
661 const Eigen::MatrixXd &local_pts = data.local_pts;
662 const Eigen::MatrixXd &global_pts = data.global_pts;
663 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
664
665 double lambda, mu;
666 params_.lambda_mu(local_pts, global_pts, t, el_id, lambda, mu);
667
668 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
669 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
670
671 stress = mu * (def_grad - FmT) + lambda * std::log(def_grad.determinant()) * FmT;
672 result = mu * stress + FmT * stress.transpose() * FmT * (mu - lambda * std::log(def_grad.determinant())) + lambda * (FmT.array() * stress.array()).sum() * FmT;
673 }
674
676 const OptAssemblerData &data,
677 const Eigen::MatrixXd &vect,
678 Eigen::MatrixXd &stress,
679 Eigen::MatrixXd &result) const
680 {
681 const double t = data.t;
682 const int el_id = data.el_id;
683 const Eigen::MatrixXd &local_pts = data.local_pts;
684 const Eigen::MatrixXd &global_pts = data.global_pts;
685 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
686
687 double lambda, mu;
688 params_.lambda_mu(local_pts, global_pts, t, el_id, lambda, mu);
689
690 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
691 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
692
693 stress = mu * (def_grad - FmT) + lambda * std::log(def_grad.determinant()) * FmT;
694 result.setZero(size() * size(), size());
695 if (vect.rows() == 1)
696 for (int i = 0; i < size(); ++i)
697 for (int j = 0; j < size(); ++j)
698 for (int l = 0; l < size(); ++l)
699 {
700 result(i * size() + j, l) += mu * vect(i) * ((j == l) ? 1 : 0);
701 // For some reason, the second gives lower error. From the formula, though, it should be the following.
702 // result(i * size() + j, l) += (mu - lambda * std::log(def_grad.determinant())) * FmT(j, l) * (FmT.col(i).array() * vect.transpose().array()).sum();
703 result(i * size() + j, l) += (mu - lambda * std::log(def_grad.determinant())) * FmT(i, l) * (FmT.col(j).array() * vect.transpose().array()).sum();
704 result(i * size() + j, l) += lambda * FmT(i, j) * (FmT.col(l).array() * vect.transpose().array()).sum();
705 }
706 else
707 for (int i = 0; i < size(); ++i)
708 for (int j = 0; j < size(); ++j)
709 for (int k = 0; k < size(); ++k)
710 {
711 result(i * size() + j, k) += mu * vect(j) * ((i == k) ? 1 : 0);
712 // For some reason, the second gives lower error. From the formula, though, it should be the following.
713 // result(i * size() + j, k) += (mu - lambda * std::log(def_grad.determinant())) * FmT(k, j) * (FmT.row(i).array() * vect.transpose().array()).sum();
714 result(i * size() + j, k) += (mu - lambda * std::log(def_grad.determinant())) * FmT(k, i) * (FmT.row(j).array() * vect.transpose().array()).sum();
715 result(i * size() + j, k) += lambda * FmT(i, j) * (FmT.row(k).array() * vect.transpose().array()).sum();
716 }
717 }
718
720 const OptAssemblerData &data,
721 Eigen::MatrixXd &dstress_dmu,
722 Eigen::MatrixXd &dstress_dlambda) const
723 {
724 const double t = data.t;
725 const int el_id = data.el_id;
726 const Eigen::MatrixXd &local_pts = data.local_pts;
727 const Eigen::MatrixXd &global_pts = data.global_pts;
728 const Eigen::MatrixXd &grad_u_i = data.grad_u_i;
729
730 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
731 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
732 dstress_dmu = def_grad - FmT;
733 dstress_dlambda = std::log(def_grad.determinant()) * FmT;
734 }
735
736 std::map<std::string, Assembler::ParamFunc> NeoHookeanElasticity::parameters() const
737 {
738 std::map<std::string, ParamFunc> res;
739 const auto &params = params_;
740 const int size = this->size();
741
742 res["lambda"] = [&params](const RowVectorNd &uv, const RowVectorNd &p, double t, int e) {
743 double lambda, mu;
744
745 params.lambda_mu(uv, p, t, e, lambda, mu);
746 return lambda;
747 };
748
749 res["mu"] = [&params](const RowVectorNd &uv, const RowVectorNd &p, double t, int e) {
750 double lambda, mu;
751
752 params.lambda_mu(uv, p, t, e, lambda, mu);
753 return mu;
754 };
755
756 res["E"] = [&params, size](const RowVectorNd &uv, const RowVectorNd &p, double t, int e) {
757 double lambda, mu;
758 params.lambda_mu(uv, p, t, e, lambda, mu);
759
760 if (size == 3)
761 return mu * (3.0 * lambda + 2.0 * mu) / (lambda + mu);
762 else
763 return 2 * mu * (2.0 * lambda + 2.0 * mu) / (lambda + 2.0 * mu);
764 };
765
766 res["nu"] = [&params, size](const RowVectorNd &uv, const RowVectorNd &p, double t, int e) {
767 double lambda, mu;
768
769 params.lambda_mu(uv, p, t, e, lambda, mu);
770
771 if (size == 3)
772 return lambda / (2.0 * (lambda + mu));
773 else
774 return lambda / (lambda + 2.0 * mu);
775 };
776
777 return res;
778 }
779
780} // namespace polyfem::assembler
double val
Definition Assembler.cpp:86
ElementAssemblyValues vals
Definition Assembler.cpp:22
int y
int z
int x
std::string stress() const
Definition Units.hpp:25
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,...
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
Eigen::VectorXd eval_deformed_jacobian_determinant(const Eigen::VectorXd &disp) const
void lambda_mu(double px, double py, double pz, double x, double y, double z, double t, int el_id, double &lambda, double &mu) const
void add_multimaterial(const int index, const json &params, const bool is_volume, const std::string &stress_unit)
double compute_energy(const NonLinearAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() const override
void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const
void compute_dstress_dmu_dlambda(const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const override
void compute_stress_grad_multiply_vect(const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, 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
void add_multimaterial(const int index, const json &params, const Units &units) override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
T compute_energy_aux(const NonLinearAssemblerData &data) const
void compute_stiffness_value(const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const override
const basis::ElementBases & bs
const Eigen::MatrixXd & fun
const Eigen::MatrixXd & local_pts
const basis::ElementBases & gbs
Used for test only.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AutoDiffGradMat
Eigen::Matrix< double, dim, dim > hat(const Eigen::Matrix< double, dim, 1 > &x)
Eigen::Matrix< double, dim, 1 > cross(const Eigen::Matrix< double, dim, 1 > &x, const Eigen::Matrix< double, dim, 1 > &y)
Eigen::Matrix< T, Eigen::Dynamic, 1 > AutoDiffVect
void neo_hookean_3d_function(const AutodiffHessianPt &pt, const double lambda, const double mu, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > &res)
void neo_hookean_2d_function(const AutodiffHessianPt &pt, const double lambda, const double mu, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > &res)
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
T determinant(const Eigen::Matrix< T, rows, cols, option, maxRow, maxCol > &mat)
Eigen::MatrixXd pk2_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
void compute_disp_grad_at_quad(const assembler::NonLinearAssemblerData &data, const AutoDiffVect &local_disp, const int p, const int size, AutoDiffGradMat &def_grad)
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
nlohmann::json json
Definition Common.hpp:9
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)
void get_local_disp(const assembler::NonLinearAssemblerData &data, const int size, AutoDiffVect &local_disp)
Eigen::MatrixXd pk1_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13