36 static scalar_type frobenius_product_trans(
const base_matrix &A,
37 const base_matrix &B) {
39 scalar_type res = scalar_type(0);
42 res += A(i, j) * B(j, i);
46 struct compute_invariants {
51 scalar_type i1_, i2_, i3_, j1_, j2_;
52 bool i1_c, i2_c, i3_c, j1_c, j2_c;
54 base_matrix di1, di2, di3, dj1, dj2;
55 bool di1_c, di2_c, di3_c, dj1_c, dj2_c;
57 base_tensor ddi1, ddi2, ddi3, ddj1, ddj2;
58 bool ddi1_c, ddi2_c, ddi3_c, ddj1_c, ddj2_c;
69 gmm::copy(gmm::identity_matrix(), di1);
74 ddi1 = base_tensor(N, N, N, N);
78 inline scalar_type i1()
79 {
if (!i1_c) compute_i1();
return i1_; }
81 inline const base_matrix &grad_i1()
82 {
if (!di1_c) compute_di1();
return di1; }
84 inline const base_tensor &sym_grad_grad_i1()
85 {
if (!ddi1_c) compute_ddi1();
return ddi1; }
90 i2_ = (gmm::sqr(gmm::mat_trace(E))
91 - frobenius_product_trans(E, E)) / scalar_type(2);
97 gmm::copy(gmm::identity_matrix(), di2);
98 gmm::scale(di2, i1());
100 gmm::add(gmm::scaled(E, -scalar_type(1)), di2);
104 void compute_ddi2() {
105 ddi2 = base_tensor(N, N, N, N);
108 ddi2(i,i,k,k) += scalar_type(1);
111 ddi2(i,j,j,i) -= scalar_type(1)/scalar_type(2);
112 ddi2(j,i,j,i) -= scalar_type(1)/scalar_type(2);
117 inline scalar_type i2()
118 {
if (!i2_c) compute_i2();
return i2_; }
120 inline const base_matrix &grad_i2()
121 {
if (!di2_c) compute_di2();
return di2; }
123 inline const base_tensor &sym_grad_grad_i2()
124 {
if (!ddi2_c) compute_ddi2();
return ddi2; }
129 i3_ = bgeot::lu_inverse(&(*(Einv.begin())), gmm::mat_nrows(Einv));
134 scalar_type det = i3();
139 gmm::scale(di3, det);
143 void compute_ddi3() {
144 ddi3 = base_tensor(N, N, N, N);
145 scalar_type det = i3() / scalar_type(2);
150 ddi3(i,j,k,l) = det*(Einv(j,i)*Einv(l,k) - Einv(j,k)*Einv(l,i)
151 + Einv(i,j)*Einv(l,k) - Einv(i,k)*Einv(l,j));
155 inline scalar_type i3()
156 {
if (!i3_c) compute_i3();
return i3_; }
158 inline const base_matrix &grad_i3()
159 {
if (!di3_c) compute_di3();
return di3; }
161 inline const base_tensor &sym_grad_grad_i3()
162 {
if (!ddi3_c) compute_ddi3();
return ddi3; }
166 j1_ = i1() * ::pow(gmm::abs(i3()), -scalar_type(1) / scalar_type(3));
172 gmm::add(gmm::scaled(grad_i3(), -i1() / (scalar_type(3) * i3())), dj1);
173 gmm::scale(dj1, ::pow(gmm::abs(i3()), -scalar_type(1) / scalar_type(3)));
177 void compute_ddj1() {
178 const base_matrix &di1_ = grad_i1();
179 const base_matrix &di3_ = grad_i3();
180 scalar_type coeff1 = scalar_type(1) / (scalar_type(3)*i3());
181 scalar_type coeff2 = scalar_type(4) * coeff1 * coeff1 * i1();
182 ddj1 = sym_grad_grad_i3();
183 gmm::scale(ddj1.as_vector(), -i1() * coeff1);
190 (di3_(i, j) * di3_(k, l)) * coeff2
191 - (di1_(i, j) * di3_(k, l) + di1_(k, l) * di3_(i, j)) * coeff1;
193 gmm::scale(ddj1.as_vector(),
194 ::pow(gmm::abs(i3()), -scalar_type(1)/scalar_type(3)));
198 inline scalar_type j1()
199 {
if (!j1_c) compute_j1();
return j1_; }
201 inline const base_matrix &grad_j1()
202 {
if (!dj1_c) compute_dj1();
return dj1; }
204 inline const base_tensor &sym_grad_grad_j1()
205 {
if (!ddj1_c) compute_ddj1();
return ddj1; }
209 j2_ = i2() * ::pow(gmm::abs(i3()), -scalar_type(2) / scalar_type(3));
215 gmm::add(gmm::scaled(grad_i3(), -scalar_type(2) * i2() / (scalar_type(3) * i3())), dj2);
216 gmm::scale(dj2, ::pow(gmm::abs(i3()), -scalar_type(2) / scalar_type(3)));
220 void compute_ddj2() {
221 const base_matrix &di2_ = grad_i2();
222 const base_matrix &di3_ = grad_i3();
223 scalar_type coeff1 = scalar_type(2) / (scalar_type(3)*i3());
224 scalar_type coeff2 = scalar_type(5) * coeff1 * coeff1 * i2()
226 ddj2 = sym_grad_grad_i2();
227 gmm::add(gmm::scaled(sym_grad_grad_i3().as_vector(), -i2() * coeff1),
235 (di3_(i, j) * di3_(k, l)) * coeff2
236 - (di2_(i, j) * di3_(k, l) + di2_(k, l) * di3_(i, j)) * coeff1;
238 gmm::scale(ddj2.as_vector(),
239 ::pow(gmm::abs(i3()), -scalar_type(2)/scalar_type(3)));
244 inline scalar_type j2()
245 {
if (!j2_c) compute_j2();
return j2_; }
247 inline const base_matrix &grad_j2()
248 {
if (!dj2_c) compute_dj2();
return dj2; }
250 inline const base_tensor &sym_grad_grad_j2()
251 {
if (!ddj2_c) compute_ddj2();
return ddj2; }
254 compute_invariants(
const base_matrix &EE)
255 : E(EE), i1_c(false), i2_c(false), i3_c(false),
256 j1_c(false), j2_c(false), di1_c(false), di2_c(false), di3_c(false),
257 dj1_c(false), dj2_c(false), ddi1_c(false), ddi2_c(false),
258 ddi3_c(false), ddj1_c(false), ddj2_c(false)
259 { N = gmm::mat_nrows(E); }
269 int check_symmetry(
const base_tensor &t) {
275 if (gmm::abs(t(n,m,l,k) - t(l,k,n,m))>1e-5) flags &= (~1);
276 if (gmm::abs(t(n,m,l,k) - t(m,n,l,k))>1e-5) flags &= (~2);
277 if (gmm::abs(t(n,m,l,k) - t(n,m,k,l))>1e-5) flags &= (~4);
284 void abstract_hyperelastic_law::random_E(base_matrix &E) {
286 base_matrix Phi(N,N);
290 d = bgeot::lu_det(&(*(Phi.begin())), N);
291 }
while (d < scalar_type(0.01));
292 gmm::mult(gmm::transposed(Phi),Phi,E);
293 gmm::scale(E,-1.); gmm::add(gmm::identity_matrix(),E);
297 void abstract_hyperelastic_law::test_derivatives
298 (
size_type N, scalar_type h,
const base_vector& param)
const {
299 base_matrix E(N,N), E2(N,N), DE(N,N);
302 for (
size_type count = 0; count < 100; ++count) {
303 random_E(E); random_E(DE);
307 base_matrix sigma1(N,N), sigma2(N,N);
308 getfem::base_tensor tdsigma(N,N,N,N);
309 base_matrix dsigma(N,N);
310 gmm::copy(E, E2); gmm::add(DE, E2);
311 sigma(E, sigma1, param, scalar_type(1));
312 sigma(E2, sigma2, param, scalar_type(1));
314 scalar_type d = strain_energy(E2, param, scalar_type(1))
315 - strain_energy(E, param, scalar_type(1));
318 for (
size_type j=0; j < N; ++j) d2 += sigma1(i,j)*DE(i,j);
319 if (gmm::abs(d-d2)/(gmm::abs(d)+1e-40) > 1e-4) {
320 cout <<
"Test " << count <<
" wrong derivative of strain_energy, d=" 321 << d/h <<
", d2=" << d2/h << endl;
325 grad_sigma(E,tdsigma,param, scalar_type(1));
331 dsigma(i,j) += tdsigma(i,j,k,m)*DE(k,m);
334 sigma2(i,j) -= sigma1(i,j);
335 if (gmm::abs(dsigma(i,j) - sigma2(i,j))
336 /(gmm::abs(dsigma(i,j)) + 1e-40) > 1.5e-4) {
337 cout <<
"Test " << count <<
" wrong derivative of sigma, i=" 338 << i <<
", j=" << j <<
", dsigma=" << dsigma(i,j)/h
339 <<
", var sigma = " << sigma2(i,j)/h << endl;
345 GMM_ASSERT1(ok,
"Derivative test has failed");
349 (
const base_matrix& F,
const base_matrix &E,
350 base_matrix &cauchy_stress,
const base_vector ¶ms,
351 scalar_type det_trans)
const 354 base_matrix PK2(N,N);
355 sigma(E,PK2,params,det_trans);
356 base_matrix aux(N,N);
357 gmm::mult(F,PK2,aux);
358 gmm::mult(aux,gmm::transposed(F),cauchy_stress);
359 gmm::scale(cauchy_stress,scalar_type(1.0/det_trans));
364 (
const base_matrix& F,
const base_matrix& E,
365 const base_vector ¶ms, scalar_type det_trans,
366 base_tensor &grad_sigma_ul)
const 369 base_tensor Cse(N,N,N,N);
370 grad_sigma(E,Cse,params,det_trans);
371 scalar_type mult = 1.0/det_trans;
379 grad_sigma_ul(i,j,k,l) = 0.0;
384 grad_sigma_ul(i,j,k,l)+=
385 F(i,m)*F(j,n)*F(k,p)*F(l,q)*Cse(m,n,p,q);
387 grad_sigma_ul(i,j,k,l) *= mult;
391 scalar_type SaintVenant_Kirchhoff_hyperelastic_law::strain_energy
392 (
const base_matrix &E,
const base_vector ¶ms, scalar_type det_trans)
const {
394 if (det_trans <= scalar_type(0))
397 return gmm::sqr(gmm::mat_trace(E)) * params[0] / scalar_type(2)
398 + gmm::mat_euclidean_norm_sqr(E) * params[1];
401 void SaintVenant_Kirchhoff_hyperelastic_law::sigma
402 (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
const {
403 gmm::copy(gmm::identity_matrix(), result);
404 gmm::scale(result, params[0] * gmm::mat_trace(E));
405 gmm::add(gmm::scaled(E, 2 * params[1]), result);
406 if (det_trans <= scalar_type(0)) {
407 gmm::add(gmm::scaled(E, 1e200), result);
410 void SaintVenant_Kirchhoff_hyperelastic_law::grad_sigma
411 (
const base_matrix &E, base_tensor &result,
const base_vector ¶ms, scalar_type)
const {
412 std::fill(result.begin(), result.end(), scalar_type(0));
416 result(i, i, l, l) += params[0];
417 result(i, l, i, l) += params[1]/scalar_type(2);
418 result(i, l, l, i) += params[1]/scalar_type(2);
419 result(l, i, i, l) += params[1]/scalar_type(2);
420 result(l, i, l, i) += params[1]/scalar_type(2);
425 const base_matrix& E,
426 const base_vector ¶ms,
427 scalar_type det_trans,
428 base_tensor &grad_sigma_ul)
const 431 base_tensor Cse(N,N,N,N);
432 grad_sigma(E,Cse,params,det_trans);
433 base_matrix Cinv(N,N);
434 gmm::mult(F,gmm::transposed(F),Cinv);
435 scalar_type mult=1.0/det_trans;
440 grad_sigma_ul(i, j, k, l)= (Cinv(i,j)*Cinv(k,l)*params[0] +
441 params[1]*(Cinv(i,k)*Cinv(j,l) + Cinv(i,l)*Cinv(j,k)))*mult;
444 SaintVenant_Kirchhoff_hyperelastic_law::SaintVenant_Kirchhoff_hyperelastic_law() {
448 scalar_type membrane_elastic_law::strain_energy
449 (
const base_matrix & ,
const base_vector & , scalar_type)
const {
451 GMM_ASSERT1(
false,
"To be done");
455 void membrane_elastic_law::sigma
456 (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
const {
458 base_tensor tt(2,2,2,2);
459 size_type N = (gmm::mat_nrows(E) > 2)? 2 : gmm::mat_nrows(E);
460 grad_sigma(E,tt,params, det_trans);
466 result(i,j)+=tt(i,j,k,l)*E(k,l);
469 if(params[4]!=0) result(0,0)+=params[4];
471 if(params[5]!=0) result(1,1)+=params[5];
475 void membrane_elastic_law::grad_sigma
476 (
const base_matrix & , base_tensor &result,
477 const base_vector ¶ms, scalar_type)
const {
479 std::fill(result.begin(), result.end(), scalar_type(0));
480 scalar_type poisonXY=params[0]*params[1]/params[2];
482 scalar_type G=( params[3] == 0) ? params[0]/(2*(1+params[1])) : params[3];
483 std::fill(result.begin(), result.end(), scalar_type(0));
484 result(0,0,0,0) = params[0]/(1-params[1]*poisonXY);
487 result(0,0,1,1) = params[1]*params[0]/(1-params[1]*poisonXY);
488 result(1,1,0,0) = params[1]*params[0]/(1-params[1]*poisonXY);
491 result(1,1,1,1) = params[2]/(1-params[1]*poisonXY);
502 scalar_type Mooney_Rivlin_hyperelastic_law::strain_energy
503 (
const base_matrix &E,
const base_vector ¶ms
504 ,scalar_type det_trans)
const {
506 if (compressible && det_trans <= scalar_type(0))
return 1e200;
508 GMM_ASSERT1(N == 3,
"Mooney Rivlin hyperelastic law only defined " 509 "on dimension 3, sorry");
511 gmm::scale(C, scalar_type(2));
512 gmm::add(gmm::identity_matrix(), C);
513 compute_invariants ci(C);
515 scalar_type C1 = params[i++];
516 scalar_type W = C1 * (ci.j1() - scalar_type(3));
518 scalar_type C2 = params[i++];
519 W += C2 * (ci.j2() - scalar_type(3));
522 scalar_type D1 = params[i++];
523 W += D1 * gmm::sqr(sqrt(gmm::abs(ci.i3())) - scalar_type(1));
528 void Mooney_Rivlin_hyperelastic_law::sigma
529 (
const base_matrix &E, base_matrix &result,
530 const base_vector ¶ms, scalar_type det_trans)
const {
532 GMM_ASSERT1(N == 3,
"Mooney Rivlin hyperelastic law only defined " 533 "on dimension 3, sorry");
535 gmm::scale(C, scalar_type(2));
536 gmm::add(gmm::identity_matrix(), C);
537 compute_invariants ci(C);
540 scalar_type C1 = params[i++];
541 gmm::copy(gmm::scaled(ci.grad_j1(), scalar_type(2) * C1), result);
543 scalar_type C2 = params[i++];
544 gmm::add(gmm::scaled(ci.grad_j2(), scalar_type(2) * C2), result);
547 scalar_type D1 = params[i++];
548 scalar_type di3 = D1 - D1 / sqrt(gmm::abs(ci.i3()));
549 gmm::add(gmm::scaled(ci.grad_i3(), scalar_type(2) * di3), result);
552 if (det_trans <= scalar_type(0))
553 gmm::add(gmm::scaled(C, 1e200), result);
557 void Mooney_Rivlin_hyperelastic_law::grad_sigma
558 (
const base_matrix &E, base_tensor &result,
559 const base_vector ¶ms, scalar_type)
const {
561 GMM_ASSERT1(N == 3,
"Mooney Rivlin hyperelastic law only defined " 562 "on dimension 3, sorry");
564 gmm::scale(C, scalar_type(2));
565 gmm::add(gmm::identity_matrix(), C);
566 compute_invariants ci(C);
569 scalar_type C1 = params[i++];
570 gmm::copy(gmm::scaled(ci.sym_grad_grad_j1().as_vector(),
571 scalar_type(4)*C1), result.as_vector());
573 scalar_type C2 = params[i++];
574 gmm::add(gmm::scaled(ci.sym_grad_grad_j2().as_vector(),
575 scalar_type(4)*C2), result.as_vector());
578 scalar_type D1 = params[i++];
579 scalar_type di3 = D1 - D1 / sqrt(gmm::abs(ci.i3()));
580 gmm::add(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
581 scalar_type(4)*di3), result.as_vector());
584 scalar_type A22 = D1 / (scalar_type(2) * pow(gmm::abs(ci.i3()), 1.5));
585 const base_matrix &di = ci.grad_i3();
590 result(l1, l2, l3, l4) +=
591 scalar_type(4) * A22 * di(l1, l2) * di(l3, l4);
598 Mooney_Rivlin_hyperelastic_law::Mooney_Rivlin_hyperelastic_law
599 (
bool compressible_,
bool neohookean_)
600 : compressible(compressible_), neohookean(neohookean_)
603 if (compressible) ++nb_params_;
604 if (neohookean) --nb_params_;
611 scalar_type Neo_Hookean_hyperelastic_law::strain_energy
612 (
const base_matrix &E,
const base_vector ¶ms, scalar_type det_trans)
const {
613 if (det_trans <= scalar_type(0))
return 1e200;
615 GMM_ASSERT1(N == 3,
"Neo Hookean hyperelastic law only defined " 616 "on dimension 3, sorry");
618 gmm::scale(C, scalar_type(2));
619 gmm::add(gmm::identity_matrix(), C);
620 compute_invariants ci(C);
622 scalar_type lambda = params[0];
623 scalar_type mu = params[1];
624 scalar_type logi3 = log(ci.i3());
625 scalar_type W = mu/2 * (ci.i1() - scalar_type(3) - logi3);
627 W += lambda/8 * gmm::sqr(logi3);
629 W += lambda/4 * (ci.i3() - scalar_type(1) - logi3);
634 void Neo_Hookean_hyperelastic_law::sigma
635 (
const base_matrix &E, base_matrix &result,
636 const base_vector ¶ms , scalar_type det_trans)
const {
638 GMM_ASSERT1(N == 3,
"Neo Hookean hyperelastic law only defined " 639 "on dimension 3, sorry");
641 gmm::scale(C, scalar_type(2));
642 gmm::add(gmm::identity_matrix(), C);
643 compute_invariants ci(C);
645 scalar_type lambda = params[0];
646 scalar_type mu = params[1];
647 gmm::copy(gmm::scaled(ci.grad_i1(), mu), result);
649 gmm::add(gmm::scaled(ci.grad_i3(),
650 (lambda/2 * log(ci.i3()) - mu) / ci.i3()), result);
652 gmm::add(gmm::scaled(ci.grad_i3(),
653 lambda/2 - lambda/(2*ci.i3()) - mu / ci.i3()), result);
654 if (det_trans <= scalar_type(0))
655 gmm::add(gmm::scaled(C, 1e200), result);
658 void Neo_Hookean_hyperelastic_law::grad_sigma
659 (
const base_matrix &E, base_tensor &result,
660 const base_vector ¶ms, scalar_type)
const {
662 GMM_ASSERT1(N == 3,
"Neo Hookean hyperelastic law only defined " 663 "on dimension 3, sorry");
665 gmm::scale(C, scalar_type(2));
666 gmm::add(gmm::identity_matrix(), C);
667 compute_invariants ci(C);
669 scalar_type lambda = params[0];
670 scalar_type mu = params[1];
674 scalar_type logi3 = log(ci.i3());
675 gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
676 (lambda * logi3 - 2*mu) / ci.i3()),
678 coeff = (lambda + 2 * mu - lambda * logi3) / gmm::sqr(ci.i3());
680 gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
681 lambda - (lambda + 2 * mu) / ci.i3()),
683 coeff = (lambda + 2 * mu) / gmm::sqr(ci.i3());
686 const base_matrix &di = ci.grad_i3();
691 result(l1, l2, l3, l4) += coeff * di(l1, l2) * di(l3, l4);
697 Neo_Hookean_hyperelastic_law::Neo_Hookean_hyperelastic_law(
bool bonet_)
705 scalar_type generalized_Blatz_Ko_hyperelastic_law::strain_energy
706 (
const base_matrix &E,
const base_vector ¶ms, scalar_type det_trans)
const {
707 if (det_trans <= scalar_type(0))
return 1e200;
708 scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
709 scalar_type n = params[4];
711 GMM_ASSERT1(N == 3,
"Generalized Blatz Ko hyperelastic law only defined " 712 "on dimension 3, sorry");
714 gmm::scale(C, scalar_type(2));
715 gmm::add(gmm::identity_matrix(), C);
716 compute_invariants ci(C);
718 return pow(a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
719 + c*ci.i2() / ci.i3() + d, n);
722 void generalized_Blatz_Ko_hyperelastic_law::sigma
723 (
const base_matrix &E, base_matrix &result,
724 const base_vector ¶ms, scalar_type det_trans)
const {
725 scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
726 scalar_type n = params[4];
728 GMM_ASSERT1(N == 3,
"Generalized Blatz Ko hyperelastic law only defined " 729 "on dimension 3, sorry");
731 gmm::scale(C, scalar_type(2));
732 gmm::add(gmm::identity_matrix(), C);
733 compute_invariants ci(C);
735 scalar_type z = a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
736 + c*ci.i2() / ci.i3() + d;
737 scalar_type nz = n * pow(z, n-1.);
738 scalar_type di1 = nz * a;
739 scalar_type di2 = nz * c / ci.i3();
740 scalar_type di3 = nz *
741 (b / (2. * sqrt(gmm::abs(ci.i3()))) - c * ci.i2() / gmm::sqr(ci.i3()));
743 gmm::copy(gmm::scaled(ci.grad_i1(), di1 * 2.0), result);
744 gmm::add(gmm::scaled(ci.grad_i2(), di2 * 2.0), result);
745 gmm::add(gmm::scaled(ci.grad_i3(), di3 * 2.0), result);
746 if (det_trans <= scalar_type(0))
747 gmm::add(gmm::scaled(C, 1e200), result);
751 void generalized_Blatz_Ko_hyperelastic_law::grad_sigma
752 (
const base_matrix &E, base_tensor &result,
753 const base_vector ¶ms, scalar_type)
const {
754 scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
755 scalar_type n = params[4];
757 GMM_ASSERT1(N == 3,
"Generalized Blatz Ko hyperelastic law only defined " 758 "on dimension 3, sorry");
760 gmm::scale(C, scalar_type(2));
761 gmm::add(gmm::identity_matrix(), C);
762 compute_invariants ci(C);
765 scalar_type z = a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
766 + c*ci.i2() / ci.i3() + d;
767 scalar_type nz = n * pow(z, n-1.);
768 scalar_type di1 = nz * a;
769 scalar_type di2 = nz * c / ci.i3();
770 scalar_type y = (b / (2. * sqrt(gmm::abs(ci.i3()))) - c * ci.i2() / gmm::sqr(ci.i3()));
771 scalar_type di3 = nz * y;
773 gmm::copy(gmm::scaled(ci.sym_grad_grad_i1().as_vector(),
774 scalar_type(4)*di1), result.as_vector());
775 gmm::add(gmm::scaled(ci.sym_grad_grad_i2().as_vector(),
776 scalar_type(4)*di2), result.as_vector());
777 gmm::add(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
778 scalar_type(4)*di3), result.as_vector());
780 scalar_type nnz = n * (n-1.) * pow(z, n-2.);
782 A(0, 0) = nnz * a * a;
783 A(1, 0) = A(0, 1) = nnz * a * c / ci.i3();
784 A(2, 0) = A(0, 2) = nnz * a * y;
785 A(1, 1) = nnz * c * c / gmm::sqr(ci.i3());
786 A(2, 1) = A(1, 2) = nnz * y * c / ci.i3() - nz * c / gmm::sqr(ci.i3());
787 A(2, 2) = nnz * y * y + nz * (2. * c * ci.i2() / pow(ci.i3(), 3.) - b / (4. * pow(ci.i3(), 1.5)));
789 typedef const base_matrix * pointer_base_matrix__;
790 pointer_base_matrix__ di[3];
791 di[0] = &(ci.grad_i1());
792 di[1] = &(ci.grad_i2());
793 di[2] = &(ci.grad_i3());
801 result(l1, l2, l3, l4)
802 += 4. * A(j, k) * (*di[j])(l1, l2) * (*di[k])(l3, l4);
809 generalized_Blatz_Ko_hyperelastic_law::generalized_Blatz_Ko_hyperelastic_law() {
812 V[0] = 1.0; V[1] = 1.0, V[2] = 1.5; V[3] = -0.5; V[4] = 1.5;
816 scalar_type Ciarlet_Geymonat_hyperelastic_law::strain_energy
817 (
const base_matrix &E,
const base_vector ¶ms, scalar_type det_trans)
const {
818 if (det_trans <= scalar_type(0))
return 1e200;
820 scalar_type a = params[2];
821 scalar_type b = params[1]/scalar_type(2) - params[2];
822 scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
824 scalar_type d = params[0]/scalar_type(2) + params[1];
825 scalar_type e = -(scalar_type(3)*(a+b) + c);
827 gmm::copy(gmm::scaled(E, scalar_type(2)), C);
828 gmm::add(gmm::identity_matrix(), C);
829 scalar_type det = bgeot::lu_det(&(*(C.begin())), N);
830 return a * gmm::mat_trace(C)
831 + b * (gmm::sqr(gmm::mat_trace(C)) -
832 gmm::mat_euclidean_norm_sqr(C))/scalar_type(2)
833 + c * det - d * log(det) / scalar_type(2) + e;
836 void Ciarlet_Geymonat_hyperelastic_law::sigma
837 (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
const {
839 scalar_type a = params[2];
840 scalar_type b = params[1]/scalar_type(2) - params[2];
841 scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
843 scalar_type d = params[0]/scalar_type(2) + params[1];
845 if (a > params[1]/scalar_type(2)
846 || a < params[1]/scalar_type(2) - params[0]/scalar_type(4) || a < 0)
847 GMM_WARNING1(
"Inconsistent third parameter for Ciarlet-Geymonat " 849 gmm::copy(gmm::scaled(E, scalar_type(2)), C);
850 gmm::add(gmm::identity_matrix(), C);
851 gmm::copy(gmm::identity_matrix(), result);
852 gmm::scale(result, scalar_type(2) * (a + b * gmm::mat_trace(C)));
853 gmm::add(gmm::scaled(C, -scalar_type(2) * b), result);
854 if (det_trans <= scalar_type(0))
855 gmm::add(gmm::scaled(C, 1e200), result);
857 scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
858 gmm::add(gmm::scaled(C, scalar_type(2) * c * det - d), result);
862 void Ciarlet_Geymonat_hyperelastic_law::grad_sigma
863 (
const base_matrix &E, base_tensor &result,
const base_vector ¶ms, scalar_type)
const {
866 scalar_type b2 = params[1] - params[2]*scalar_type(2);
867 scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
869 scalar_type d = params[0]/scalar_type(2) + params[1];
871 gmm::copy(gmm::scaled(E, scalar_type(2)), C);
872 gmm::add(gmm::identity_matrix(), C);
873 scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
874 std::fill(result.begin(), result.end(), scalar_type(0));
877 result(i, i, j, j) += 2*b2;
878 result(i, j, i, j) -= b2;
879 result(i, j, j, i) -= b2;
882 result(i, j, k, l) +=
883 (C(i, k)*C(l, j) + C(i, l)*C(k, j)) * (d-scalar_type(2)*det*c)
884 + (C(i, j) * C(k, l)) * det*c*scalar_type(4);
892 int levi_civita(
int i,
int j,
int k) {
896 return static_cast<int> 897 (int(- 1)*(
static_cast<int>(pow(
double(ii-jj),2.))%3)
898 * (
static_cast<int> (pow(
double(ii-kk),
double(2)))%3 )
899 * (
static_cast<int> (pow(
double(jj-kk),
double(2)))%3)
900 * (pow(
double(jj-(ii%3))-
double(0.5),
double(2))-double(1.25)));
905 scalar_type plane_strain_hyperelastic_law::strain_energy
906 (
const base_matrix &E,
const base_vector ¶ms, scalar_type det_trans)
const {
907 GMM_ASSERT1(gmm::mat_nrows(E) == 2,
"Plane strain law is for 2D only.");
908 base_matrix E3D(3,3);
909 E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
910 return pl->strain_energy(E3D, params, det_trans);
913 void plane_strain_hyperelastic_law::sigma
914 (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
const {
915 GMM_ASSERT1(gmm::mat_nrows(E) == 2,
"Plane strain law is for 2D only.");
916 base_matrix E3D(3,3), result3D(3,3);
917 E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
918 pl->sigma(E3D, result3D, params, det_trans);
919 result(0,0) = result3D(0,0); result(1,0) = result3D(1,0);
920 result(0,1) = result3D(0,1); result(1,1) = result3D(1,1);
923 void plane_strain_hyperelastic_law::grad_sigma
924 (
const base_matrix &E, base_tensor &result,
const base_vector ¶ms, scalar_type det_trans)
const {
925 GMM_ASSERT1(gmm::mat_nrows(E) == 2,
"Plane strain law is for 2D only.");
926 base_matrix E3D(3,3);
927 base_tensor result3D(3,3,3,3);
928 E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
929 pl->grad_sigma(E3D, result3D, params, det_trans);
930 result(0,0,0,0) = result3D(0,0,0,0); result(1,0,0,0) = result3D(1,0,0,0);
931 result(0,1,0,0) = result3D(0,1,0,0); result(1,1,0,0) = result3D(1,1,0,0);
932 result(0,0,1,0) = result3D(0,0,1,0); result(1,0,1,0) = result3D(1,0,1,0);
933 result(0,1,1,0) = result3D(0,1,1,0); result(1,1,1,0) = result3D(1,1,1,0);
934 result(0,0,0,1) = result3D(0,0,0,1); result(1,0,0,1) = result3D(1,0,0,1);
935 result(0,1,0,1) = result3D(0,1,0,1); result(1,1,0,1) = result3D(1,1,0,1);
936 result(0,0,1,1) = result3D(0,0,1,1); result(1,0,1,1) = result3D(1,0,1,1);
937 result(0,1,1,1) = result3D(0,1,1,1); result(1,1,1,1) = result3D(1,1,1,1);
954 phyperelastic_law AHL;
957 const model::varnamelist &vl,
958 const model::varnamelist &dl,
959 const model::mimlist &mims,
960 model::real_matlist &matl,
961 model::real_veclist &vecl,
962 model::real_veclist &,
964 build_version version)
const {
965 GMM_ASSERT1(mims.size() == 1,
966 "Nonlinear elasticity brick need a single mesh_im");
967 GMM_ASSERT1(vl.size() == 1,
968 "Nonlinear elasticity brick need a single variable");
969 GMM_ASSERT1(dl.size() == 1,
970 "Wrong number of data for nonlinear elasticity brick, " 971 << dl.size() <<
" should be 1 (vector).");
972 GMM_ASSERT1(matl.size() == 1,
"Wrong number of terms for nonlinear " 979 const model_real_plain_vector ¶ms = md.
real_variable(dl[0]);
983 if (mf_params) sl = sl * mf_params->
get_qdim() / mf_params->
nb_dof();
984 GMM_ASSERT1(sl == AHL->nb_params(),
"Wrong number of coefficients for the " 985 "nonlinear constitutive elastic law");
990 if (version & model::BUILD_MATRIX) {
992 GMM_TRACE2(
"Nonlinear elasticity stiffness matrix assembly");
994 (matl[0], mim, mf_u, u, mf_params, params, *AHL, rg);
998 if (version & model::BUILD_RHS) {
999 asm_nonlinear_elasticity_rhs(vecl[0], mim,
1000 mf_u, u, mf_params, params, *AHL, rg);
1001 gmm::scale(vecl[0], scalar_type(-1));
1006 nonlinear_elasticity_brick(
const phyperelastic_law &AHL_)
1008 set_flags(
"Nonlinear elasticity brick",
false ,
1022 const phyperelastic_law &AHL,
const std::string &dataname,
1024 pbrick pbr = std::make_shared<nonlinear_elasticity_brick>(AHL);
1027 tl.push_back(model::term_description(varname, varname,
true));
1028 model::varnamelist dl(1, dataname);
1029 model::varnamelist vl(1, varname);
1030 return md.
add_brick(pbr, vl, dl, tl, model::mimlist(1,&mim), region);
1037 void compute_Von_Mises_or_Tresca(
model &md,
1038 const std::string &varname,
1039 const phyperelastic_law &AHL,
1040 const std::string &dataname,
1042 model_real_plain_vector &VM,
1044 GMM_ASSERT1(gmm::vect_size(VM) == mf_vm.
nb_dof(),
1045 "The vector has not the good size");
1047 const model_real_plain_vector &u = md.
real_variable(varname);
1049 const model_real_plain_vector ¶ms = md.
real_variable(dataname);
1052 if (mf_params) sl = sl * mf_params->
get_qdim() / mf_params->
nb_dof();
1053 GMM_ASSERT1(sl == AHL->nb_params(),
"Wrong number of coefficients for " 1054 "the nonlinear constitutive elastic law");
1056 unsigned N = unsigned(mf_u.linked_mesh().dim());
1057 unsigned NP = unsigned(AHL->nb_params()), NFem = mf_u.get_qdim();
1058 model_real_plain_vector GRAD(mf_vm.
nb_dof()*NFem*N);
1059 model_real_plain_vector PARAMS(mf_vm.
nb_dof()*NP);
1060 if (mf_params)
interpolation(*mf_params, mf_vm, params, PARAMS);
1062 base_matrix E(N, N), gradphi(NFem,N),gradphit(N,NFem), Id(N, N),
1063 sigmahathat(N,N),aux(NFem,N), sigma(NFem,NFem),
1066 if (!mf_params) gmm::copy(params, p);
1067 base_vector eig(NFem);
1068 base_vector ez(NFem);
1070 gmm::copy(gmm::identity_matrix(), Id);
1071 gmm::copy(gmm::identity_matrix(), IdNFem);
1073 gmm::resize(gradphi,NFem,N);
1074 std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
1076 gmm::copy(gmm::transposed(gradphit),gradphi);
1077 for (
unsigned int alpha = 0; alpha <N; ++alpha)
1078 gradphi(alpha, alpha)+=1;
1079 gmm::mult(gmm::transposed(gradphi), gradphi, E);
1080 gmm::add(gmm::scaled(Id, -scalar_type(1)), E);
1081 gmm::scale(E, scalar_type(1)/scalar_type(2));
1083 gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*NP,NP)), p);
1084 AHL->sigma(E, sigmahathat, p, scalar_type(1));
1085 if (NFem == 3 && N == 2) {
1087 for (
unsigned int l = 0; l <NFem; ++l) {
1089 for (
unsigned int m = 0; m <NFem; ++m)
1090 for (
unsigned int n = 0; n <NFem; ++n){
1091 ez[l]+=levi_civita(l,m,n)*gradphi(m,0)*gradphi(n,1);
1093 normEz= gmm::vect_norm2(ez);
1097 gmm::mult(gradphi, sigmahathat, aux);
1098 gmm::mult(aux, gmm::transposed(gradphi), sigma);
1102 if (NFem == 3 && N == 2) {
1103 gmm::resize(gradphi,NFem,NFem);
1104 for (
unsigned int ll = 0; ll <NFem; ++ll)
1105 for (
unsigned int ii = 0; ii <NFem; ++ii)
1106 for (
unsigned int jj = 0; jj <NFem; ++jj)
1107 gradphi(ll,2)+=(levi_civita(ll,ii,jj)*gradphi(ii,0)
1108 *gradphi(jj,1))/normEz;
1112 gmm::scale(sigma, scalar_type(1) / bgeot::lu_det(&(*(gradphi.begin())), NFem));
1116 gmm::add(gmm::scaled(IdNFem, -gmm::mat_trace(sigma) / NFem), sigma);
1119 VM[i] = sqrt(3.0/2)*gmm::mat_euclidean_norm(sigma);
1123 gmm::symmetric_qr_algorithm(sigma, eig);
1124 std::sort(eig.begin(), eig.end());
1125 VM[i] = eig.back() - eig.front();
1131 void compute_sigmahathat(
model &md,
1132 const std::string &varname,
1133 const phyperelastic_law &AHL,
1134 const std::string &dataname,
1136 model_real_plain_vector &SIGMA) {
1138 const model_real_plain_vector &u = md.
real_variable(varname);
1140 const model_real_plain_vector ¶ms = md.
real_variable(dataname);
1143 if (mf_params) sl = sl * mf_params->
get_qdim() / mf_params->
nb_dof();
1144 GMM_ASSERT1(sl == AHL->nb_params(),
"Wrong number of coefficients for " 1145 "the nonlinear constitutive elastic law");
1148 unsigned NP = unsigned(AHL->nb_params()), NFem = mf_u.
get_qdim();
1149 GMM_ASSERT1(mf_sigma.
nb_dof() > 0,
"Bad mf_sigma");
1153 GMM_ASSERT1(((ratio == 1) || (ratio == N*N)) &&
1154 (gmm::vect_size(SIGMA) == mf_sigma.
nb_dof()*ratio),
1155 "The vector has not the good size");
1157 model_real_plain_vector GRAD(mf_sigma.
nb_dof()*ratio*NFem/N);
1158 model_real_plain_vector PARAMS(mf_sigma.
nb_dof()*NP);
1161 getfem::mesh_trans_inv mti(mf_sigma.
linked_mesh());
1169 base_matrix E(N, N), gradphi(NFem,N),gradphit(N,NFem), Id(N, N),
1170 sigmahathat(N,N),aux(NFem,N), sigma(NFem,NFem),
1175 if (!mf_params) gmm::copy(params, p);
1176 base_vector eig(NFem);
1177 base_vector ez(NFem);
1178 gmm::copy(gmm::identity_matrix(), Id);
1179 gmm::copy(gmm::identity_matrix(), IdNFem);
1189 gmm::resize(gradphi,NFem,N);
1190 std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
1193 gmm::copy(gmm::transposed(gradphit),gradphi);
1194 for (
unsigned int alpha = 0; alpha <N; ++alpha)
1195 gradphi(alpha, alpha) += scalar_type(1);
1196 gmm::mult(gmm::transposed(gradphi), gradphi, E);
1197 gmm::add(gmm::scaled(Id, -scalar_type(1)), E);
1198 gmm::scale(E, scalar_type(1)/scalar_type(2));
1200 gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*ratio*NP,NP)),p);
1202 AHL->sigma(E, sigmahathat, p, scalar_type(1));
1204 std::copy(sigmahathat.begin(), sigmahathat.end(), SIGMA.begin()+i*N*N);
1216 struct nonlinear_incompressibility_brick :
public virtual_brick {
1219 const model::varnamelist &vl,
1220 const model::varnamelist &dl,
1221 const model::mimlist &mims,
1222 model::real_matlist &matl,
1223 model::real_veclist &vecl,
1224 model::real_veclist &veclsym,
1226 build_version version)
const {
1228 GMM_ASSERT1(matl.size() == 2,
"Wrong number of terms for nonlinear " 1229 "incompressibility brick");
1230 GMM_ASSERT1(dl.size() == 0,
"Nonlinear incompressibility brick need no " 1232 GMM_ASSERT1(mims.size() == 1,
"Nonlinear incompressibility brick need a " 1234 GMM_ASSERT1(vl.size() == 2,
"Wrong number of variables for nonlinear " 1235 "incompressibility brick");
1241 const mesh_im &mim = *mims[0];
1245 if (version & model::BUILD_MATRIX) {
1246 gmm::clear(matl[0]);
1247 gmm::clear(matl[1]);
1248 asm_nonlinear_incomp_tangent_matrix(matl[0], matl[1],
1249 mim, mf_u, mf_p, u, p, rg);
1252 if (version & model::BUILD_RHS) {
1253 asm_nonlinear_incomp_rhs(vecl[0], veclsym[1], mim, mf_u, mf_p,u,p, rg);
1254 gmm::scale(vecl[0], scalar_type(-1));
1255 gmm::scale(veclsym[1], scalar_type(-1));
1260 nonlinear_incompressibility_brick() {
1261 set_flags(
"Nonlinear incompressibility brick",
1272 const std::string &multname,
size_type region) {
1273 pbrick pbr = std::make_shared<nonlinear_incompressibility_brick>();
1275 tl.push_back(model::term_description(varname, varname,
true));
1276 tl.push_back(model::term_description(varname, multname,
true));
1277 model::varnamelist vl(1, varname);
1278 vl.push_back(multname);
1279 model::varnamelist dl;
1280 return md.
add_brick(pbr, vl, dl, tl, model::mimlist(1, &mim), region);
1294 static void ga_init_scalar_(bgeot::multi_index &mi) { mi.resize(0); }
1295 static void ga_init_square_matrix_(bgeot::multi_index &mi,
size_type N)
1296 { mi.resize(2); mi[0] = mi[1] = N; }
1301 struct matrix_i2_operator :
public ga_nonlinear_operator {
1302 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1303 if (args.size() != 1 || args[0]->sizes().size() != 2
1304 || args[0]->sizes()[0] != args[0]->sizes()[1])
return false;
1305 ga_init_scalar_(sizes);
1310 void value(
const arg_list &args, base_tensor &result)
const {
1312 const base_tensor &t = *args[0];
1313 scalar_type tr = scalar_type(0);
1314 for (
size_type i = 0; i < N; ++i) tr += t[i*N+i];
1315 scalar_type tr2 = scalar_type(0);
1318 tr2 += t[i+ j*N] * t[j + i*N];
1319 result[0] = (tr*tr - tr2)/2;
1323 void derivative(
const arg_list &args,
size_type,
1324 base_tensor &result)
const {
1326 const base_tensor &t = *args[0];
1327 scalar_type tr = scalar_type(0);
1328 for (
size_type i = 0; i < N; ++i) tr += t[i*N+i];
1329 base_tensor::iterator it = result.begin();
1332 *it = ((i == j) ? tr : scalar_type(0)) - t[i*N+j];
1333 GMM_ASSERT1(it == result.end(),
"Internal error");
1338 base_tensor &result)
const {
1340 gmm::clear(result.as_vector());
1343 result[(N+1)*(i+N*N*j)] += scalar_type(1);
1344 result[(N+1)*N*j + i*(N*N*N + 1)] -= scalar_type(1);
1351 struct matrix_j1_operator :
public ga_nonlinear_operator {
1352 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1353 if (args.size() != 1 || args[0]->sizes().size() != 2
1354 || args[0]->sizes()[0] != args[0]->sizes()[1])
return false;
1355 ga_init_scalar_(sizes);
1360 void value(
const arg_list &args, base_tensor &result)
const {
1362 base_matrix M(N, N);
1363 gmm::copy(args[0]->as_vector(), M.as_vector());
1364 scalar_type det = bgeot::lu_det(&(*(M.begin())), N);
1365 scalar_type tr = scalar_type(0);
1366 for (
size_type i = 0; i < N; ++i) tr += M(i,i);
1368 result[0] = tr / pow(det, scalar_type(1)/scalar_type(3));
1374 void derivative(
const arg_list &args,
size_type,
1375 base_tensor &result)
const {
1377 base_matrix M(N, N);
1378 gmm::copy(args[0]->as_vector(), M.as_vector());
1379 scalar_type tr = scalar_type(0);
1380 for (
size_type i = 0; i < N; ++i) tr += M(i,i);
1381 scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1383 base_tensor::iterator it = result.begin();
1386 *it = (((i == j) ? scalar_type(1) : scalar_type(0))
1387 - tr*M(j,i)/scalar_type(3))
1388 / pow(det, scalar_type(1)/scalar_type(3));
1389 GMM_ASSERT1(it == result.end(),
"Internal error");
1391 std::fill(result.begin(), result.end(), 1.E200);
1397 base_tensor &result)
const {
1399 base_matrix M(N, N);
1400 gmm::copy(args[0]->as_vector(), M.as_vector());
1401 scalar_type tr = scalar_type(0);
1402 for (
size_type i = 0; i < N; ++i) tr += M(i,i);
1403 scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1405 base_tensor::iterator it = result.begin();
1410 *it = (- ((k == l) ? M(j, i) : scalar_type(0))
1412 - ((i == j) ? M(l, k) : scalar_type(0))
1413 + tr*M(j,i)*M(k,l)/ scalar_type(3))
1414 / (scalar_type(3)*pow(det, scalar_type(1)/scalar_type(3)));
1415 GMM_ASSERT1(it == result.end(),
"Internal error");
1417 std::fill(result.begin(), result.end(), 1.E200);
1423 struct matrix_j2_operator :
public ga_nonlinear_operator {
1424 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1425 if (args.size() != 1 || args[0]->sizes().size() != 2
1426 || args[0]->sizes()[0] != args[0]->sizes()[1])
return false;
1427 ga_init_scalar_(sizes);
1432 void value(
const arg_list &args, base_tensor &result)
const {
1434 base_matrix M(N, N);
1435 gmm::copy(args[0]->as_vector(), M.as_vector());
1436 scalar_type tr = scalar_type(0);
1437 for (
size_type i = 0; i < N; ++i) tr += M(i,i);
1438 scalar_type tr2 = scalar_type(0);
1441 tr2 += M(i,j)*M(j,i);
1442 scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
1443 scalar_type det = bgeot::lu_det(&(*(M.begin())), N);
1446 result[0] = i2 / pow(det, scalar_type(2)/scalar_type(3));
1452 void derivative(
const arg_list &args,
size_type,
1453 base_tensor &result)
const {
1455 const base_tensor &t = *args[0];
1456 base_matrix M(N, N);
1457 gmm::copy(t.as_vector(), M.as_vector());
1458 scalar_type tr = scalar_type(0);
1459 for (
size_type i = 0; i < N; ++i) tr += M(i,i);
1460 scalar_type tr2 = scalar_type(0);
1463 tr2 += M(i,j)*M(j,i);
1464 scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
1465 scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1466 base_tensor::iterator it = result.begin();
1469 *it = (((i == j) ? tr : scalar_type(0)) - t[j+N*i]
1470 - scalar_type(2)*i2*M(j,i)/scalar_type(3))
1471 / pow(det, scalar_type(2)/scalar_type(3));
1472 GMM_ASSERT1(it == result.end(),
"Internal error");
1477 base_tensor &result)
const {
1479 const base_tensor &t = *args[0];
1480 base_matrix M(N, N);
1481 gmm::copy(t.as_vector(), M.as_vector());
1482 scalar_type tr = scalar_type(0);
1483 for (
size_type i = 0; i < N; ++i) tr += M(i,i);
1484 scalar_type tr2 = scalar_type(0);
1487 tr2 += M(i,j)*M(j,i);
1488 scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
1489 scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1490 base_tensor::iterator it = result.begin();
1495 *it = ( ((i==j) ? 1. : 0.) * ((k==l) ? 1. : 0.)
1496 - ((i==l) ? 1. : 0.) * ((k==j) ? 1. : 0.)
1497 - 2.*tr*M(j,i)*((k==l) ? 1. : 0.)/3.
1498 + 2.*tr*M(j,i)*M(l,k)/3.
1499 - 2.*i2*M(i,k)*M(l,j)/3.
1500 - 2.*((tr*((i==j) ? 1. : 0.))-t[j+N*i]
1501 - 2.*i2*M(j,i)/3)*M(l,k)/3.)
1502 / pow(det, scalar_type(2)/scalar_type(3));
1507 struct Right_Cauchy_Green_operator :
public ga_nonlinear_operator {
1508 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1509 if (args.size() != 1 || args[0]->sizes().size() != 2)
return false;
1510 ga_init_square_matrix_(sizes, args[0]->sizes()[1]);
1515 void value(
const arg_list &args, base_tensor &result)
const {
1517 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1518 base_tensor::iterator it = result.begin();
1520 for (
size_type i = 0; i < n; ++i, ++it) {
1521 *it = scalar_type(0);
1523 *it += (*(args[0]))[i*m+k] * (*(args[0]))[j*m+k];
1529 void derivative(
const arg_list &args,
size_type,
1530 base_tensor &result)
const {
1531 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1532 base_tensor::iterator it = result.begin();
1536 for (
size_type i = 0; i < n; ++i, ++it) {
1537 *it = scalar_type(0);
1538 if (l == i) *it += (*(args[0]))[j*m+k];
1539 if (l == j) *it += (*(args[0]))[i*m+k];
1541 GMM_ASSERT1(it == result.end(),
"Internal error");
1548 base_tensor &result)
const {
1549 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1550 base_tensor::iterator it = result.begin();
1556 for (
size_type i = 0; i < n; ++i, ++it) {
1557 *it = scalar_type(0);
1559 if (l == i && p == j) *it += scalar_type(1);
1560 if (p == i && l == j) *it += scalar_type(1);
1563 GMM_ASSERT1(it == result.end(),
"Internal error");
1568 struct Left_Cauchy_Green_operator :
public ga_nonlinear_operator {
1569 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1570 if (args.size() != 1 || args[0]->sizes().size() != 2)
return false;
1571 ga_init_square_matrix_(sizes, args[0]->sizes()[0]);
1576 void value(
const arg_list &args, base_tensor &result)
const {
1578 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1579 base_tensor::iterator it = result.begin();
1581 for (
size_type i = 0; i < m; ++i, ++it) {
1582 *it = scalar_type(0);
1584 *it += (*(args[0]))[k*m+i] * (*(args[0]))[k*m+j];
1590 void derivative(
const arg_list &args,
size_type,
1591 base_tensor &result)
const {
1592 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1593 base_tensor::iterator it = result.begin();
1597 for (
size_type i = 0; i < m; ++i, ++it) {
1598 *it = scalar_type(0);
1599 if (k == i) *it += (*(args[0]))[l*m+j];
1600 if (k == j) *it += (*(args[0]))[l*m+i];
1602 GMM_ASSERT1(it == result.end(),
"Internal error");
1609 base_tensor &result)
const {
1610 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1611 base_tensor::iterator it = result.begin();
1617 for (
size_type i = 0; i < m; ++i, ++it) {
1618 *it = scalar_type(0);
1620 if (k == i && o == j) *it += scalar_type(1);
1621 if (o == i && k == j) *it += scalar_type(1);
1624 GMM_ASSERT1(it == result.end(),
"Internal error");
1630 struct Green_Lagrangian_operator :
public ga_nonlinear_operator {
1631 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1632 if (args.size() != 1 || args[0]->sizes().size() != 2)
return false;
1633 ga_init_square_matrix_(sizes, args[0]->sizes()[1]);
1638 void value(
const arg_list &args, base_tensor &result)
const {
1640 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1641 base_tensor::iterator it = result.begin();
1643 for (
size_type i = 0; i < n; ++i, ++it) {
1644 *it = scalar_type(0);
1646 *it += (*(args[0]))[i*m+k]*(*(args[0]))[j*m+k]*scalar_type(0.5);
1647 if (i == j) *it -= scalar_type(0.5);
1653 void derivative(
const arg_list &args,
size_type,
1654 base_tensor &result)
const {
1655 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1656 base_tensor::iterator it = result.begin();
1660 for (
size_type i = 0; i < n; ++i, ++it) {
1661 *it = scalar_type(0);
1662 if (l == i) *it += (*(args[0]))[j*m+k]*scalar_type(0.5);
1663 if (l == j) *it += (*(args[0]))[i*m+k]*scalar_type(0.5);
1665 GMM_ASSERT1(it == result.end(),
"Internal error");
1672 base_tensor &result)
const {
1673 size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1674 base_tensor::iterator it = result.begin();
1680 for (
size_type i = 0; i < n; ++i, ++it) {
1681 *it = scalar_type(0);
1683 if (l == i && p == j) *it += scalar_type(0.5);
1684 if (p == i && l == j) *it += scalar_type(0.5);
1687 GMM_ASSERT1(it == result.end(),
"Internal error");
1693 struct Cauchy_stress_from_PK2 :
public ga_nonlinear_operator {
1694 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1695 if (args.size() != 2 || args[0]->sizes().size() != 2
1696 || args[1]->sizes().size() != 2
1697 || args[0]->sizes()[0] != args[0]->sizes()[1]
1698 || args[1]->sizes()[0] != args[0]->sizes()[1]
1699 || args[1]->sizes()[1] != args[0]->sizes()[1])
return false;
1700 ga_init_square_matrix_(sizes, args[0]->sizes()[1]);
1705 void value(
const arg_list &args, base_tensor &result)
const {
1707 base_matrix F(N, N), sigma(N,N), aux(N, N);
1708 gmm::copy(args[0]->as_vector(), sigma.as_vector());
1709 gmm::copy(args[1]->as_vector(), F.as_vector());
1710 gmm::add(gmm::identity_matrix(), F);
1711 gmm::mult(F, sigma, aux);
1712 gmm::mult(aux, gmm::transposed(F), sigma);
1713 scalar_type det = bgeot::lu_det(&(*(F.begin())), N);
1714 gmm::scale(sigma, scalar_type(1)/det);
1715 gmm::copy(sigma.as_vector(), result.as_vector());
1719 void derivative(
const arg_list &args,
size_type nder,
1720 base_tensor &result)
const {
1722 base_matrix F(N, N);
1723 gmm::copy(args[1]->as_vector(), F.as_vector());
1724 gmm::add(gmm::identity_matrix(), F);
1725 scalar_type det = bgeot::lu_det(&(*(F.begin())), N);
1727 base_tensor::iterator it = result.begin();
1735 *it = F(i,k) * F(j,l) / det;
1743 base_matrix sigma(N,N), aux(N,N), aux2(N,N);
1744 gmm::copy(args[0]->as_vector(), sigma.as_vector());
1745 gmm::mult(sigma, gmm::transposed(F), aux);
1746 gmm::mult(F, aux, aux2);
1747 bgeot::lu_inverse(&(*(F.begin())), N);
1751 for (
size_type i = 0; i < N; ++i, ++it) {
1752 *it = scalar_type(0);
1753 if (i == k) *it += aux(l, j) / det;
1754 if (l == j) *it += aux(k, i) / det;
1755 *it -= aux2(i,j) * F(l,k) / det;
1760 GMM_ASSERT1(it == result.end(),
"Internal error");
1765 base_tensor &)
const {
1766 GMM_ASSERT1(
false,
"Sorry, not implemented");
1771 struct AHL_wrapper_sigma :
public ga_nonlinear_operator {
1772 phyperelastic_law AHL;
1773 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1774 if (args.size() != 2 || args[0]->sizes().size() != 2
1775 || args[1]->size() != AHL->nb_params()
1776 || args[0]->sizes()[0] != args[0]->sizes()[1])
return false;
1777 ga_init_square_matrix_(sizes, args[0]->sizes()[0]);
1782 void value(
const arg_list &args, base_tensor &result)
const {
1784 base_vector params(AHL->nb_params());
1785 gmm::copy(args[1]->as_vector(), params);
1786 base_matrix Gu(N, N), E(N,N), sigma(N,N);
1787 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1788 gmm::mult(gmm::transposed(Gu), Gu, E);
1789 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1790 gmm::scale(E, scalar_type(0.5));
1791 gmm::add(gmm::identity_matrix(), Gu);
1792 scalar_type det = bgeot::lu_det(&(*(Gu.begin())), N);
1794 AHL->sigma(E, sigma, params, det);
1795 gmm::copy(sigma.as_vector(), result.as_vector());
1799 void derivative(
const arg_list &args,
size_type nder,
1800 base_tensor &result)
const {
1802 base_vector params(AHL->nb_params());
1803 gmm::copy(args[1]->as_vector(), params);
1804 base_tensor grad_sigma(N, N, N, N);
1805 base_matrix Gu(N, N), E(N,N);
1806 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1807 gmm::mult(gmm::transposed(Gu), Gu, E);
1808 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1809 gmm::scale(E, scalar_type(0.5));
1810 gmm::add(gmm::identity_matrix(), Gu);
1811 scalar_type det = bgeot::lu_det(&(*(Gu.begin())), N);
1813 GMM_ASSERT1(nder == 1,
"Sorry, the derivative of this hyperelastic " 1814 "law with respect to its parameters is not available.");
1816 AHL->grad_sigma(E, grad_sigma, params, det);
1818 base_tensor::iterator it = result.begin();
1822 for (
size_type i = 0; i < N; ++i, ++it) {
1823 *it = scalar_type(0);
1825 *it += grad_sigma(i,j,m,l) * Gu(k, m);
1827 GMM_ASSERT1(it == result.end(),
"Internal error");
1833 base_tensor &)
const {
1834 GMM_ASSERT1(
false,
"Sorry, second derivative not implemented");
1837 AHL_wrapper_sigma(
const phyperelastic_law &A) : AHL(A) {}
1842 struct AHL_wrapper_potential :
public ga_nonlinear_operator {
1843 phyperelastic_law AHL;
1844 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1845 if (args.size() != 2 || args[0]->sizes().size() != 2
1846 || args[1]->size() != AHL->nb_params()
1847 || args[0]->sizes()[0] != args[0]->sizes()[1])
return false;
1848 ga_init_scalar_(sizes);
1853 void value(
const arg_list &args, base_tensor &result)
const {
1855 base_vector params(AHL->nb_params());
1856 gmm::copy(args[1]->as_vector(), params);
1857 base_matrix Gu(N, N), E(N,N);
1858 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1859 gmm::mult(gmm::transposed(Gu), Gu, E);
1860 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1861 gmm::scale(E, scalar_type(0.5));
1862 gmm::add(gmm::identity_matrix(), Gu);
1863 scalar_type det = bgeot::lu_det(&(*(Gu.begin())), N);
1865 result[0] = AHL->strain_energy(E, params, det);
1869 void derivative(
const arg_list &args,
size_type nder,
1870 base_tensor &result)
const {
1872 base_vector params(AHL->nb_params());
1873 gmm::copy(args[1]->as_vector(), params);
1874 base_matrix Gu(N, N), E(N,N), sigma(N,N);
1875 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1876 gmm::mult(gmm::transposed(Gu), Gu, E);
1877 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1878 gmm::scale(E, scalar_type(0.5));
1879 gmm::add(gmm::identity_matrix(), Gu);
1880 scalar_type det = bgeot::lu_det(&(*(Gu.begin())), N);
1882 GMM_ASSERT1(nder == 1,
"Sorry, Cannot derive the potential with " 1883 "respect to law parameters.");
1885 AHL->sigma(E, sigma, params, det);
1886 gmm::mult(Gu, sigma, E);
1887 gmm::copy(E.as_vector(), result.as_vector());
1892 void second_derivative(
const arg_list &args,
size_type nder1,
1893 size_type nder2, base_tensor &result)
const {
1896 base_vector params(AHL->nb_params());
1897 gmm::copy(args[1]->as_vector(), params);
1898 base_tensor grad_sigma(N, N, N, N);
1899 base_matrix Gu(N, N), E(N,N), sigma(N,N);
1900 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1901 gmm::mult(gmm::transposed(Gu), Gu, E);
1902 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1903 gmm::scale(E, scalar_type(0.5));
1904 gmm::add(gmm::identity_matrix(), Gu);
1905 scalar_type det = bgeot::lu_det(&(*(Gu.begin())), N);
1907 GMM_ASSERT1(nder1 == 1 && nder2 == 1,
"Sorry, Cannot derive the " 1908 "potential with respect to law parameters.");
1910 AHL->sigma(E, sigma, params, det);
1911 AHL->grad_sigma(E, grad_sigma, params, det);
1913 base_tensor::iterator it = result.begin();
1917 for (
size_type i = 0; i < N; ++i, ++it) {
1918 *it = scalar_type(0);
1919 if (i == k) *it += sigma(l,j);
1923 *it += grad_sigma(n,j,m,l) * Gu(k, m) * Gu(i, n);
1925 GMM_ASSERT1(it == result.end(),
"Internal error");
1929 AHL_wrapper_potential(
const phyperelastic_law &A) : AHL(A) {}
1936 struct Saint_Venant_Kirchhoff_sigma :
public ga_nonlinear_operator {
1937 bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
const {
1938 if (args.size() != 2 || args[0]->sizes().size() != 2
1939 || args[1]->size() != 2
1940 || args[0]->sizes()[0] != args[0]->sizes()[1])
return false;
1941 ga_init_square_matrix_(sizes, args[0]->sizes()[0]);
1946 void value(
const arg_list &args, base_tensor &result)
const {
1948 scalar_type lambda = (*(args[1]))[0], mu = (*(args[1]))[1];
1949 base_matrix Gu(N, N), E(N,N);
1950 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1951 gmm::mult(gmm::transposed(Gu), Gu, E);
1952 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1953 gmm::scale(E, scalar_type(0.5));
1954 scalar_type trE = gmm::mat_trace(E);
1956 base_tensor::iterator it = result.begin();
1958 for (
size_type i = 0; i < N; ++i, ++it) {
1960 if (i == j) *it += lambda*trE;
1969 void derivative(
const arg_list &args,
size_type nder,
1970 base_tensor &result)
const {
1972 scalar_type lambda = (*(args[1]))[0], mu = (*(args[1]))[1], trE;
1973 base_matrix Gu(N, N), E(N,N);
1974 gmm::copy(args[0]->as_vector(), Gu.as_vector());
1976 gmm::mult(gmm::transposed(Gu), Gu, E);
1977 gmm::add(Gu, E); gmm::add(gmm::transposed(Gu), E);
1978 gmm::scale(E, scalar_type(0.5));
1980 base_tensor::iterator it = result.begin();
1986 for (
size_type i = 0; i < N; ++i, ++it) {
1987 *it = scalar_type(0);
1988 if (i == j && k == l) *it += lambda;
1989 if (i == j) *it += lambda*Gu(k,l);
1990 if (i == k && j == l) *it += mu;
1991 if (i == l && j == k) *it += mu;
1992 if (i == l) *it += mu*Gu(k,j);
1993 if (l == j) *it += mu*Gu(k,i);
1998 trE = gmm::mat_trace(E);
2000 for (
size_type i = 0; i < N; ++i, ++it) {
2001 *it = scalar_type(0);
if (i == j) *it += trE;
2005 for (
size_type i = 0; i < N; ++i, ++it) {
2009 default: GMM_ASSERT1(
false,
"Internal error");
2011 GMM_ASSERT1(it == result.end(),
"Internal error");
2016 base_tensor &)
const {
2017 GMM_ASSERT1(
false,
"Sorry, second derivative not implemented");
2023 static bool init_predef_operators() {
2025 ga_predef_operator_tab &PREDEF_OPERATORS
2028 PREDEF_OPERATORS.add_method
2029 (
"Matrix_i2", std::make_shared<matrix_i2_operator>());
2030 PREDEF_OPERATORS.add_method
2031 (
"Matrix_j1", std::make_shared<matrix_j1_operator>());
2032 PREDEF_OPERATORS.add_method
2033 (
"Matrix_j2", std::make_shared<matrix_j2_operator>());
2034 PREDEF_OPERATORS.add_method
2035 (
"Right_Cauchy_Green", std::make_shared<Right_Cauchy_Green_operator>());
2036 PREDEF_OPERATORS.add_method
2037 (
"Left_Cauchy_Green", std::make_shared<Left_Cauchy_Green_operator>());
2038 PREDEF_OPERATORS.add_method
2039 (
"Green_Lagrangian", std::make_shared<Green_Lagrangian_operator>());
2041 PREDEF_OPERATORS.add_method
2042 (
"Cauchy_stress_from_PK2", std::make_shared<Cauchy_stress_from_PK2>());
2044 PREDEF_OPERATORS.add_method
2045 (
"Saint_Venant_Kirchhoff_sigma",
2046 std::make_shared<Saint_Venant_Kirchhoff_sigma>());
2047 PREDEF_OPERATORS.add_method
2048 (
"Saint_Venant_Kirchhoff_potential",
2049 std::make_shared<AHL_wrapper_potential>
2050 (std::make_shared<SaintVenant_Kirchhoff_hyperelastic_law>()));
2051 PREDEF_OPERATORS.add_method
2052 (
"Plane_Strain_Saint_Venant_Kirchhoff_sigma",
2053 std::make_shared<Saint_Venant_Kirchhoff_sigma>());
2054 PREDEF_OPERATORS.add_method
2055 (
"Plane_Strain_Saint_Venant_Kirchhoff_potential",
2056 std::make_shared<AHL_wrapper_potential>
2057 (std::make_shared<SaintVenant_Kirchhoff_hyperelastic_law>()));
2059 phyperelastic_law gbklaw
2060 = std::make_shared<generalized_Blatz_Ko_hyperelastic_law>();
2061 PREDEF_OPERATORS.add_method
2062 (
"Generalized_Blatz_Ko_sigma",
2063 std::make_shared<AHL_wrapper_sigma>(gbklaw));
2064 PREDEF_OPERATORS.add_method
2065 (
"Generalized_Blatz_Ko_potential",
2066 std::make_shared<AHL_wrapper_potential>
2067 (std::make_shared<generalized_Blatz_Ko_hyperelastic_law>()));
2068 PREDEF_OPERATORS.add_method
2069 (
"Plane_Strain_Generalized_Blatz_Ko_sigma",
2070 std::make_shared<AHL_wrapper_sigma>
2071 (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
2072 PREDEF_OPERATORS.add_method
2073 (
"Plane_Strain_Generalized_Blatz_Ko_potential",
2074 std::make_shared<AHL_wrapper_potential>
2075 (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
2077 phyperelastic_law cigelaw
2078 = std::make_shared<Ciarlet_Geymonat_hyperelastic_law>();
2079 PREDEF_OPERATORS.add_method
2080 (
"Ciarlet_Geymonat_sigma", std::make_shared<AHL_wrapper_sigma>(cigelaw));
2081 PREDEF_OPERATORS.add_method
2082 (
"Ciarlet_Geymonat_potential",
2083 std::make_shared<AHL_wrapper_potential>
2084 (std::make_shared<Ciarlet_Geymonat_hyperelastic_law>()));
2085 PREDEF_OPERATORS.add_method
2086 (
"Plane_Strain_Ciarlet_Geymonat_sigma",
2087 std::make_shared<AHL_wrapper_sigma>
2088 (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
2089 PREDEF_OPERATORS.add_method
2090 (
"Plane_Strain_Ciarlet_Geymonat_potential",
2091 std::make_shared<AHL_wrapper_potential>
2092 (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
2094 phyperelastic_law morilaw
2095 = std::make_shared<Mooney_Rivlin_hyperelastic_law>();
2096 PREDEF_OPERATORS.add_method
2097 (
"Incompressible_Mooney_Rivlin_sigma",
2098 std::make_shared<AHL_wrapper_sigma>(morilaw));
2099 PREDEF_OPERATORS.add_method
2100 (
"Incompressible_Mooney_Rivlin_potential",
2101 std::make_shared<AHL_wrapper_potential>
2102 (std::make_shared<Mooney_Rivlin_hyperelastic_law>()));
2103 PREDEF_OPERATORS.add_method
2104 (
"Plane_Strain_Incompressible_Mooney_Rivlin_sigma",
2105 std::make_shared<AHL_wrapper_sigma>
2106 (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
2107 PREDEF_OPERATORS.add_method
2108 (
"Plane_Strain_Incompressible_Mooney_Rivlin_potential",
2109 std::make_shared<AHL_wrapper_potential>
2110 (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
2112 phyperelastic_law cmorilaw
2113 = std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true);
2114 PREDEF_OPERATORS.add_method
2115 (
"Compressible_Mooney_Rivlin_sigma",
2116 std::make_shared<AHL_wrapper_sigma>(cmorilaw));
2117 PREDEF_OPERATORS.add_method
2118 (
"Compressible_Mooney_Rivlin_potential",
2119 std::make_shared<AHL_wrapper_potential>
2120 (std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true)));
2121 PREDEF_OPERATORS.add_method
2122 (
"Plane_Strain_Compressible_Mooney_Rivlin_sigma",
2123 std::make_shared<AHL_wrapper_sigma>
2124 (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
2125 PREDEF_OPERATORS.add_method
2126 (
"Plane_Strain_Compressible_Mooney_Rivlin_potential",
2127 std::make_shared<AHL_wrapper_potential>
2128 (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
2130 phyperelastic_law ineolaw
2131 = std::make_shared<Mooney_Rivlin_hyperelastic_law>(
false,
true);
2132 PREDEF_OPERATORS.add_method
2133 (
"Incompressible_Neo_Hookean_sigma",
2134 std::make_shared<AHL_wrapper_sigma>(ineolaw));
2135 PREDEF_OPERATORS.add_method
2136 (
"Incompressible_Neo_Hookean_potential",
2137 std::make_shared<AHL_wrapper_potential>
2138 (std::make_shared<Mooney_Rivlin_hyperelastic_law>(
false,
true)));
2139 PREDEF_OPERATORS.add_method
2140 (
"Plane_Strain_Incompressible_Neo_Hookean_sigma",
2141 std::make_shared<AHL_wrapper_sigma>
2142 (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
2143 PREDEF_OPERATORS.add_method
2144 (
"Plane_Strain_Incompressible_Neo_Hookean_potential",
2145 std::make_shared<AHL_wrapper_potential>
2146 (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
2148 phyperelastic_law cneolaw
2149 = std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true,
true);
2150 PREDEF_OPERATORS.add_method
2151 (
"Compressible_Neo_Hookean_sigma",
2152 std::make_shared<AHL_wrapper_sigma>(cneolaw));
2153 PREDEF_OPERATORS.add_method
2154 (
"Compressible_Neo_Hookean_potential",
2155 std::make_shared<AHL_wrapper_potential>
2156 (std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true,
true)));
2157 PREDEF_OPERATORS.add_method
2158 (
"Plane_Strain_Compressible_Neo_Hookean_sigma",
2159 std::make_shared<AHL_wrapper_sigma>
2160 (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
2161 PREDEF_OPERATORS.add_method
2162 (
"Plane_Strain_Compressible_Neo_Hookean_potential",
2163 std::make_shared<AHL_wrapper_potential>
2164 (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
2166 phyperelastic_law cneobolaw
2167 = std::make_shared<Neo_Hookean_hyperelastic_law>(
true);
2168 PREDEF_OPERATORS.add_method
2169 (
"Compressible_Neo_Hookean_Bonet_sigma",
2170 std::make_shared<AHL_wrapper_sigma>(cneobolaw));
2171 PREDEF_OPERATORS.add_method
2172 (
"Compressible_Neo_Hookean_Bonet_potential",
2173 std::make_shared<AHL_wrapper_potential>
2174 (std::make_shared<Neo_Hookean_hyperelastic_law>(
true)));
2175 PREDEF_OPERATORS.add_method
2176 (
"Plane_Strain_Compressible_Neo_Hookean_Bonet_sigma",
2177 std::make_shared<AHL_wrapper_sigma>
2178 (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
2179 PREDEF_OPERATORS.add_method
2180 (
"Plane_Strain_Compressible_Neo_Hookean_Bonet_potential",
2181 std::make_shared<AHL_wrapper_potential>
2182 (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
2184 phyperelastic_law cneocilaw
2185 = std::make_shared<Neo_Hookean_hyperelastic_law>(
false);
2186 PREDEF_OPERATORS.add_method
2187 (
"Compressible_Neo_Hookean_Ciarlet_sigma",
2188 std::make_shared<AHL_wrapper_sigma>(cneocilaw));
2189 PREDEF_OPERATORS.add_method
2190 (
"Compressible_Neo_Hookean_Ciarlet_potential",
2191 std::make_shared<AHL_wrapper_potential>
2192 (std::make_shared<Neo_Hookean_hyperelastic_law>(
false)));
2193 PREDEF_OPERATORS.add_method
2194 (
"Plane_Strain_Compressible_Neo_Hookean_Ciarlet_sigma",
2195 std::make_shared<AHL_wrapper_sigma>
2196 (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
2197 PREDEF_OPERATORS.add_method
2198 (
"Plane_Strain_Compressible_Neo_Hookean_Ciarlet_potential",
2199 std::make_shared<AHL_wrapper_potential>
2200 (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
2206 bool predef_operators_nonlinear_elasticity_initialized
2207 = init_predef_operators();
2210 std::string adapt_law_name(
const std::string &lawname,
size_type N) {
2211 std::string adapted_lawname = lawname;
2213 for (
size_type i = 0; i < lawname.size(); ++i)
2214 if (adapted_lawname[i] ==
' ') adapted_lawname[i] =
'_';
2216 if (adapted_lawname.compare(
"SaintVenant_Kirchhoff") == 0) {
2217 adapted_lawname =
"Saint_Venant_Kirchhoff";
2218 }
else if (adapted_lawname.compare(
"Saint_Venant_Kirchhoff") == 0) {
2220 }
else if (adapted_lawname.compare(
"Generalized_Blatz_Ko") == 0) {
2221 if (N == 2) adapted_lawname =
"Plane_Strain_" + adapted_lawname;
2222 }
else if (adapted_lawname.compare(
"Ciarlet_Geymonat") == 0) {
2223 if (N == 2) adapted_lawname =
"Plane_Strain_" + adapted_lawname;
2224 }
else if (adapted_lawname.compare(
"Incompressible_Mooney_Rivlin") == 0) {
2225 if (N == 2) adapted_lawname =
"Plane_Strain_" + adapted_lawname;
2226 }
else if (adapted_lawname.compare(
"Compressible_Mooney_Rivlin") == 0) {
2227 if (N == 2) adapted_lawname =
"Plane_Strain_" + adapted_lawname;
2228 }
else if (adapted_lawname.compare(
"Incompressible_Neo_Hookean") == 0) {
2229 if (N == 2) adapted_lawname =
"Plane_Strain_" + adapted_lawname;
2230 }
else if (adapted_lawname.compare(
"Compressible_Neo_Hookean") == 0 ||
2231 adapted_lawname.compare(
"Compressible_Neo_Hookean_Bonet") == 0 ||
2232 adapted_lawname.compare(
"Compressible_Neo_Hookean_Ciarlet") == 0 ) {
2233 if (N == 2) adapted_lawname =
"Plane_Strain_" + adapted_lawname;
2235 GMM_ASSERT1(
false, lawname <<
" is not a known hyperelastic law");
2237 return adapted_lawname;
2243 const std::string &varname,
const std::string ¶ms,
2245 std::string test_varname =
"Test_" + sup_previous_and_dot_to_varname(varname);
2247 GMM_ASSERT1(N >= 2 && N <= 3,
2248 "Finite strain elasticity brick works only in 2D or 3D");
2251 GMM_ASSERT1(mf,
"Finite strain elasticity brick can only be applied on " 2254 GMM_ASSERT1(Q == N,
"Finite strain elasticity brick can only be applied " 2255 "on a fem variable having the same dimension than the mesh");
2257 std::string adapted_lawname = adapt_law_name(lawname, N);
2259 std::string expr =
"((Id(meshdim)+Grad_"+varname+
")*(" + adapted_lawname
2260 +
"_sigma(Grad_"+varname+
","+params+
"))):Grad_" + test_varname;
2262 return add_nonlinear_generic_assembly_brick
2263 (md, mim, expr, region,
true,
false,
2264 "Finite strain elasticity brick for " + adapted_lawname +
" law");
2269 const std::string &multname,
size_type region) {
2270 std::string test_varname =
"Test_" + sup_previous_and_dot_to_varname(varname);
2271 std::string test_multname =
"Test_" + sup_previous_and_dot_to_varname(multname);
2274 =
"(" + test_multname+
")*(1-Det(Id(meshdim)+Grad_" + varname +
"))" 2275 +
"-(" + multname +
")*(Det(Id(meshdim)+Grad_" + varname +
")" 2276 +
"*((Inv(Id(meshdim)+Grad_" + varname +
"))':Grad_" 2277 + test_varname +
"))" ;
2278 return add_nonlinear_generic_assembly_brick
2279 (md, mim, expr, region,
true,
false,
2280 "Finite strain incompressibility brick");
2284 (
model &md,
const std::string &lawname,
const std::string &varname,
2285 const std::string ¶ms,
const mesh_fem &mf_vm,
2286 model_real_plain_vector &VM,
const mesh_region &rg) {
2289 std::string adapted_lawname = adapt_law_name(lawname, N);
2291 std::string expr =
"sqrt(3/2)*Norm(Deviator(Cauchy_stress_from_PK2(" 2292 + adapted_lawname +
"_sigma(Grad_" + varname +
"," + params +
"),Grad_" 2294 ga_interpolation_Lagrange_fem(md, expr, mf_vm, VM, rg);
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
structure used to hold a set of convexes and/or convex faces.
size_type add_nonlinear_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a nonlinear incompressibility term (for large strain elasticity) to the model with respect to the...
Non-linear elasticty and incompressibility bricks.
virtual void cauchy_updated_lagrangian(const base_matrix &F, const base_matrix &E, base_matrix &cauchy_stress, const base_vector ¶ms, scalar_type det_trans) const
True Cauchy stress (for Updated Lagrangian formulation)
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
void compute_gradient(const mesh_fem &mf, const mesh_fem &mf_target, const VECT1 &UU, VECT2 &VV)
Compute the gradient of a field on a getfem::mesh_fem.
Describe an integration method linked to a mesh.
void fill_random(L &l, double cfill)
*/
static T & instance()
Instance from the current thread.
size_type add_finite_strain_elasticity_brick(model &md, const mesh_im &mim, const std::string &lawname, const std::string &varname, const std::string ¶ms, size_type region=size_type(-1))
Add a finite strain elasticity brick to the model with respect to the variable varname (the displacem...
const mesh_fem & mesh_fem_of_variable(const std::string &name) const
Gives the access to the mesh_fem of a variable if any.
void compute_finite_strain_elasticity_Von_Mises(model &md, const std::string &lawname, const std::string &varname, const std::string ¶ms, const mesh_fem &mf_vm, model_real_plain_vector &VM, const mesh_region &rg=mesh_region::all_convexes())
Interpolate the Von-Mises stress of a field varname with respect to the nonlinear elasticity constitu...
``Model'' variables store the variables, the data and the description of a model. ...
size_t size_type
used as the common size type in the library
Model representation in Getfem.
size_type add_brick(pbrick pbr, const varnamelist &varnames, const varnamelist &datanames, const termlist &terms, const mimlist &mims, size_type region)
Add a brick to the model.
virtual dim_type get_qdim() const
Return the Q dimension.
virtual base_node point_of_basic_dof(size_type cv, size_type i) const
Return the geometrical location of a degree of freedom.
A langage for generic assembly of pde boundary value problems.
GEneric Tool for Finite Element Methods.
Describe a finite element method linked to a mesh.
The virtual brick has to be derived to describe real model bricks.
void asm_nonlinear_elasticity_tangent_matrix(const MAT &K_, const mesh_im &mim, const getfem::mesh_fem &mf, const VECT1 &U, const getfem::mesh_fem *mf_data, const VECT2 &PARAMS, const abstract_hyperelastic_law &AHL, const mesh_region &rg=mesh_region::all_convexes())
Tangent matrix for the non-linear elasticity.
linalg_traits< M >::value_type mat_trace(const M &m)
Trace of a matrix.
size_type add_nonlinear_elasticity_brick(model &md, const mesh_im &mim, const std::string &varname, const phyperelastic_law &AHL, const std::string &dataname, size_type region=size_type(-1))
Add a nonlinear (large strain) elasticity term to the model with respect to the variable varname (dep...
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
const mesh & linked_mesh() const
Give a reference to the linked mesh of type mesh.
const mesh_fem * pmesh_fem_of_variable(const std::string &name) const
Gives a pointer to the mesh_fem of a variable if any.
const model_real_plain_vector & real_variable(const std::string &name, size_type niter) const
Gives the access to the vector value of a variable.
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector ¶ms, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian
void resize(M &v, size_type m, size_type n)
*/
std::shared_ptr< const virtual_brick > pbrick
type of pointer on a brick
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector ¶ms, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian
size_type add_finite_strain_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a finite strain incompressibility term (for large strain elasticity) to the model with respect to...