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;
 
   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);
 
   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() / scalar_type(2);
 
  225       ddj2 = sym_grad_grad_i2();
 
  226       gmm::add(gmm::scaled(sym_grad_grad_i3().as_vector(), -i2() * coeff1),
 
  234                 (di3_(i, j) * di3_(k, l)) * coeff2
 
  235                 - (di2_(i, j) * di3_(k, l) + di2_(k, l) * di3_(i, j)) * coeff1;
 
  237       gmm::scale(ddj2.as_vector(),
 
  238                  ::pow(gmm::abs(i3()), -scalar_type(2)/scalar_type(3)));
 
  243     inline scalar_type j2()
 
  244     { 
if (!j2_c) compute_j2(); 
return j2_; }
 
  246     inline const base_matrix &grad_j2()
 
  247     { 
if (!dj2_c) compute_dj2(); 
return dj2; }
 
  249     inline const base_tensor &sym_grad_grad_j2()
 
  250     { 
if (!ddj2_c) compute_ddj2(); 
return ddj2; }
 
  253     compute_invariants(
const base_matrix &EE)
 
  254       : E(EE), i1_c(false), i2_c(false), i3_c(false),
 
  255         j1_c(false), j2_c(false), di1_c(false), di2_c(false), di3_c(false),
 
  256         dj1_c(false), dj2_c(false), ddi1_c(false), ddi2_c(false),
 
  257         ddi3_c(false), ddj1_c(false), ddj2_c(false)
 
  259       N = gmm::mat_nrows(E);
 
  260       i1_=i2_=i3_=j1_=j2_=0.;
 
  271   int check_symmetry(
const base_tensor &t) {
 
  277             if (gmm::abs(t(n,m,l,k) - t(l,k,n,m))>1e-5) flags &= (~1);
 
  278             if (gmm::abs(t(n,m,l,k) - t(m,n,l,k))>1e-5) flags &= (~2);
 
  279             if (gmm::abs(t(n,m,l,k) - t(n,m,k,l))>1e-5) flags &= (~4);
 
  286   void abstract_hyperelastic_law::random_E(base_matrix &E) {
 
  288     base_matrix Phi(N,N);
 
  292       d = bgeot::lu_det(&(*(Phi.begin())), N);
 
  293     } 
while (d < scalar_type(0.01));
 
  295     gmm::scale(E,-1.); 
gmm::add(gmm::identity_matrix(),E);
 
  299   void abstract_hyperelastic_law::test_derivatives
 
  300   (
size_type N, scalar_type h, 
const base_vector& param)
 const {
 
  301     base_matrix E(N,N), E2(N,N), DE(N,N);
 
  304     for (
size_type count = 0; count < 100; ++count) {
 
  305       random_E(E); random_E(DE);
 
  309       base_matrix sigma1(N,N), sigma2(N,N);
 
  310       getfem::base_tensor tdsigma(N,N,N,N);
 
  311       base_matrix dsigma(N,N);
 
  313       sigma(E, sigma1, param, scalar_type(1));
 
  314       sigma(E2, sigma2, param, scalar_type(1));
 
  316       scalar_type d = strain_energy(E2, param, scalar_type(1))
 
  317         - strain_energy(E, param, scalar_type(1));
 
  320         for (
size_type j=0; j < N; ++j) d2 += sigma1(i,j)*DE(i,j);
 
  321       if (gmm::abs(d-d2)/(gmm::abs(d)+1e-40) > 1e-4) {
 
  322         cout << 
"Test " << count << 
" wrong derivative of strain_energy, d=" 
  323              << d/h << 
", d2=" << d2/h << endl;
 
  327       grad_sigma(E,tdsigma,param, scalar_type(1));
 
  333               dsigma(i,j) += tdsigma(i,j,k,m)*DE(k,m);
 
  336           sigma2(i,j) -= sigma1(i,j);
 
  337           if (gmm::abs(dsigma(i,j) - sigma2(i,j))
 
  338               /(gmm::abs(dsigma(i,j)) + 1e-40) > 1.5e-4) {
 
  339             cout << 
"Test " << count << 
" wrong derivative of sigma, i=" 
  340                  << i << 
", j=" << j << 
", dsigma=" << dsigma(i,j)/h
 
  341                  << 
", var sigma = " << sigma2(i,j)/h << endl;
 
  347     GMM_ASSERT1(ok, 
"Derivative test has failed");
 
  351   (
const base_matrix& F, 
const base_matrix &E,
 
  352    base_matrix &cauchy_stress, 
const base_vector ¶ms,
 
  353    scalar_type det_trans)
 const 
  356     base_matrix PK2(N,N);
 
  357     sigma(E,PK2,params,det_trans);
 
  358     base_matrix aux(N,N);
 
  359     gmm::mult(F,PK2,aux);
 
  360     gmm::mult(aux,gmm::transposed(F),cauchy_stress);
 
  361     gmm::scale(cauchy_stress,scalar_type(1.0/det_trans)); 
 
  366   (
const base_matrix& F, 
const base_matrix& E,
 
  367    const base_vector ¶ms, scalar_type det_trans,
 
  368    base_tensor &grad_sigma_ul)
 const 
  371     base_tensor Cse(N,N,N,N);
 
  372     grad_sigma(E,Cse,params,det_trans);
 
  373     scalar_type mult = 1.0/det_trans;
 
  381               grad_sigma_ul(i,j,k,l) = 0.0;
 
  386                         grad_sigma_ul(i,j,k,l)+=
 
  387                           F(i,m)*F(j,n)*F(k,p)*F(l,q)*Cse(m,n,p,q);
 
  389               grad_sigma_ul(i,j,k,l) *= mult;
 
  393   scalar_type SaintVenant_Kirchhoff_hyperelastic_law::strain_energy
 
  394   (
const base_matrix &E, 
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  396     if (det_trans <= scalar_type(0))
 
  399     return gmm::sqr(gmm::mat_trace(E)) * params[0] / scalar_type(2)
 
  400     + gmm::mat_euclidean_norm_sqr(E) * params[1];
 
  403   void SaintVenant_Kirchhoff_hyperelastic_law::sigma
 
  404   (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  405     gmm::copy(gmm::identity_matrix(), result);
 
  406     gmm::scale(result, params[0] * gmm::mat_trace(E));
 
  407     gmm::add(gmm::scaled(E, 2 * params[1]), result);
 
  408     if (det_trans <= scalar_type(0)) {
 
  409       gmm::add(gmm::scaled(E, 1e200), result);
 
  412   void SaintVenant_Kirchhoff_hyperelastic_law::grad_sigma
 
  413   (
const base_matrix &E, base_tensor &result,
const base_vector ¶ms, scalar_type)
 const {
 
  414     std::fill(result.begin(), result.end(), scalar_type(0));
 
  418         result(i, i, l, l) += params[0];
 
  419         result(i, l, i, l) += params[1]/scalar_type(2);
 
  420         result(i, l, l, i) += params[1]/scalar_type(2);
 
  421         result(l, i, i, l) += params[1]/scalar_type(2);
 
  422         result(l, i, l, i) += params[1]/scalar_type(2);
 
  427     const base_matrix& E,
 
  428     const base_vector ¶ms,
 
  429     scalar_type det_trans,
 
  430     base_tensor &grad_sigma_ul)
const 
  433     base_tensor Cse(N,N,N,N);
 
  434     grad_sigma(E,Cse,params,det_trans);
 
  435     base_matrix Cinv(N,N); 
 
  436     gmm::mult(F,gmm::transposed(F),Cinv);
 
  437     scalar_type mult=1.0/det_trans;
 
  442             grad_sigma_ul(i, j, k, l)= (Cinv(i,j)*Cinv(k,l)*params[0] +
 
  443             params[1]*(Cinv(i,k)*Cinv(j,l) + Cinv(i,l)*Cinv(j,k)))*mult;
 
  446   SaintVenant_Kirchhoff_hyperelastic_law::SaintVenant_Kirchhoff_hyperelastic_law() {
 
  450   scalar_type membrane_elastic_law::strain_energy
 
  451   (
const base_matrix & , 
const base_vector & , scalar_type)
 const {
 
  453     GMM_ASSERT1(
false, 
"To be done");
 
  457   void membrane_elastic_law::sigma
 
  458   (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  460     base_tensor tt(2,2,2,2);
 
  461     size_type N = (gmm::mat_nrows(E) > 2)? 2 : gmm::mat_nrows(E);
 
  462     grad_sigma(E,tt,params, det_trans);
 
  468             result(i,j)+=tt(i,j,k,l)*E(k,l);
 
  471     if(params[4]!=0) result(0,0)+=params[4];
 
  473     if(params[5]!=0) result(1,1)+=params[5];
 
  477   void membrane_elastic_law::grad_sigma
 
  478   (
const base_matrix & , base_tensor &result,
 
  479    const base_vector ¶ms, scalar_type)
 const {
 
  481     std::fill(result.begin(), result.end(), scalar_type(0));
 
  482     scalar_type poisonXY=params[0]*params[1]/params[2]; 
 
  484     scalar_type G=( params[3] == 0) ? params[0]/(2*(1+params[1])) : params[3];
 
  485     std::fill(result.begin(), result.end(), scalar_type(0));
 
  486     result(0,0,0,0) = params[0]/(1-params[1]*poisonXY);
 
  489     result(0,0,1,1) = params[1]*params[0]/(1-params[1]*poisonXY);
 
  490     result(1,1,0,0) = params[1]*params[0]/(1-params[1]*poisonXY);
 
  493     result(1,1,1,1) = params[2]/(1-params[1]*poisonXY);
 
  504   scalar_type Mooney_Rivlin_hyperelastic_law::strain_energy
 
  505   (
const base_matrix &E, 
const base_vector ¶ms
 
  506    ,scalar_type  det_trans)
 const {
 
  508     if (compressible && det_trans <= scalar_type(0)) 
return 1e200;
 
  510     GMM_ASSERT1(N == 3, 
"Mooney Rivlin hyperelastic law only defined " 
  511                 "on dimension 3, sorry");
 
  513     gmm::scale(C, scalar_type(2));
 
  514     gmm::add(gmm::identity_matrix(), C);
 
  515     compute_invariants ci(C);
 
  517     scalar_type C1 = params[i++]; 
 
  518     scalar_type W = C1 * (ci.j1() - scalar_type(3));
 
  520       scalar_type C2 = params[i++]; 
 
  521       W += C2 * (ci.j2() - scalar_type(3));
 
  524       scalar_type D1 = params[i++];
 
  525       W += D1 * gmm::sqr(sqrt(gmm::abs(ci.i3())) - scalar_type(1));
 
  530   void Mooney_Rivlin_hyperelastic_law::sigma
 
  531   (
const base_matrix &E, base_matrix &result,
 
  532    const base_vector ¶ms, scalar_type det_trans)
 const {
 
  534     GMM_ASSERT1(N == 3, 
"Mooney Rivlin hyperelastic law only defined " 
  535                 "on dimension 3, sorry");
 
  537     gmm::scale(C, scalar_type(2));
 
  538     gmm::add(gmm::identity_matrix(), C);
 
  539     compute_invariants ci(C);
 
  542     scalar_type C1 = params[i++]; 
 
  543     gmm::copy(gmm::scaled(ci.grad_j1(), scalar_type(2) * C1), result);
 
  545       scalar_type C2 = params[i++]; 
 
  546       gmm::add(gmm::scaled(ci.grad_j2(), scalar_type(2) * C2), result);
 
  549       scalar_type D1 = params[i++];
 
  550       scalar_type di3 = D1 - D1 / sqrt(gmm::abs(ci.i3()));
 
  551       gmm::add(gmm::scaled(ci.grad_i3(), scalar_type(2) * di3), result);
 
  554       if (det_trans <= scalar_type(0))
 
  555         gmm::add(gmm::scaled(C, 1e200), result);
 
  559   void Mooney_Rivlin_hyperelastic_law::grad_sigma
 
  560   (
const base_matrix &E, base_tensor &result,
 
  561    const base_vector ¶ms, scalar_type)
 const {
 
  563     GMM_ASSERT1(N == 3, 
"Mooney Rivlin hyperelastic law only defined " 
  564                 "on dimension 3, sorry");
 
  566     gmm::scale(C, scalar_type(2));
 
  567     gmm::add(gmm::identity_matrix(), C);
 
  568     compute_invariants ci(C);
 
  571     scalar_type C1 = params[i++]; 
 
  572     gmm::copy(gmm::scaled(ci.sym_grad_grad_j1().as_vector(),
 
  573                           scalar_type(4)*C1), result.as_vector());
 
  575       scalar_type C2 = params[i++]; 
 
  576       gmm::add(gmm::scaled(ci.sym_grad_grad_j2().as_vector(),
 
  577                scalar_type(4)*C2), result.as_vector());
 
  580       scalar_type D1 = params[i++];
 
  581       scalar_type di3 = D1 - D1 / sqrt(gmm::abs(ci.i3()));
 
  582       gmm::add(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
 
  583                                scalar_type(4)*di3), result.as_vector());
 
  586       scalar_type A22 = D1 / (scalar_type(2) * pow(gmm::abs(ci.i3()), 1.5));
 
  587       const base_matrix &di = ci.grad_i3();
 
  592               result(l1, l2, l3, l4) +=
 
  593                 scalar_type(4) * A22 * di(l1, l2) * di(l3, l4);
 
  600   Mooney_Rivlin_hyperelastic_law::Mooney_Rivlin_hyperelastic_law
 
  601   (
bool compressible_, 
bool neohookean_)
 
  602   : compressible(compressible_), neohookean(neohookean_)
 
  605     if (compressible) ++nb_params_; 
 
  606     if (neohookean) --nb_params_;   
 
  613   scalar_type Neo_Hookean_hyperelastic_law::strain_energy
 
  614   (
const base_matrix &E, 
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  615    if (det_trans <= scalar_type(0)) 
return 1e200;
 
  617     GMM_ASSERT1(N == 3, 
"Neo Hookean hyperelastic law only defined " 
  618                 "on dimension 3, sorry");
 
  620     gmm::scale(C, scalar_type(2));
 
  621     gmm::add(gmm::identity_matrix(), C);
 
  622     compute_invariants ci(C);
 
  624     scalar_type lambda = params[0];
 
  625     scalar_type mu = params[1];
 
  626     scalar_type logi3 = log(ci.i3());
 
  627     scalar_type W = mu/2 * (ci.i1() - scalar_type(3) - logi3);
 
  629       W += lambda/8 * gmm::sqr(logi3);
 
  631       W += lambda/4 * (ci.i3() - scalar_type(1) - logi3);
 
  636   void Neo_Hookean_hyperelastic_law::sigma
 
  637   (
const base_matrix &E, base_matrix &result,
 
  638    const base_vector ¶ms , scalar_type det_trans)
 const {
 
  640     GMM_ASSERT1(N == 3, 
"Neo Hookean hyperelastic law only defined " 
  641                 "on dimension 3, sorry");
 
  643     gmm::scale(C, scalar_type(2));
 
  644     gmm::add(gmm::identity_matrix(), C);
 
  645     compute_invariants ci(C);
 
  647     scalar_type lambda = params[0];
 
  648     scalar_type mu = params[1];
 
  649     gmm::copy(gmm::scaled(ci.grad_i1(), mu), result);
 
  652                             (lambda/2 * log(ci.i3()) - mu) / ci.i3()), result);
 
  655                             lambda/2 - lambda/(2*ci.i3()) - mu / ci.i3()), result);
 
  656     if (det_trans <= scalar_type(0))
 
  657       gmm::add(gmm::scaled(C, 1e200), result);
 
  660   void Neo_Hookean_hyperelastic_law::grad_sigma
 
  661   (
const base_matrix &E, base_tensor &result,
 
  662    const base_vector ¶ms, scalar_type)
 const {
 
  664     GMM_ASSERT1(N == 3, 
"Neo Hookean hyperelastic law only defined " 
  665                 "on dimension 3, sorry");
 
  667     gmm::scale(C, scalar_type(2));
 
  668     gmm::add(gmm::identity_matrix(), C);
 
  669     compute_invariants ci(C);
 
  671     scalar_type lambda = params[0];
 
  672     scalar_type mu = params[1];
 
  676       scalar_type logi3 = log(ci.i3());
 
  677       gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
 
  678                             (lambda * logi3 - 2*mu) / ci.i3()),
 
  680       coeff = (lambda + 2 * mu - lambda * logi3) / gmm::sqr(ci.i3());
 
  682       gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
 
  683                             lambda  - (lambda + 2 * mu) / ci.i3()),
 
  685       coeff = (lambda + 2 * mu) / gmm::sqr(ci.i3());
 
  688     const base_matrix &di = ci.grad_i3();
 
  693             result(l1, l2, l3, l4) += coeff * di(l1, l2) * di(l3, l4);
 
  699   Neo_Hookean_hyperelastic_law::Neo_Hookean_hyperelastic_law(
bool bonet_)
 
  707   scalar_type generalized_Blatz_Ko_hyperelastic_law::strain_energy
 
  708   (
const base_matrix &E, 
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  709     if (det_trans <= scalar_type(0)) 
return 1e200;
 
  710     scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
 
  711     scalar_type n = params[4];
 
  713     GMM_ASSERT1(N == 3, 
"Generalized Blatz Ko hyperelastic law only defined " 
  714                 "on dimension 3, sorry");
 
  716     gmm::scale(C, scalar_type(2));
 
  717     gmm::add(gmm::identity_matrix(), C);
 
  718     compute_invariants ci(C);
 
  720     return pow(a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
 
  721                + c*ci.i2() / ci.i3() + d, n);
 
  724   void generalized_Blatz_Ko_hyperelastic_law::sigma
 
  725   (
const base_matrix &E, base_matrix &result,
 
  726    const base_vector ¶ms, scalar_type det_trans)
 const {
 
  727     scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
 
  728     scalar_type n = params[4];
 
  730     GMM_ASSERT1(N == 3, 
"Generalized Blatz Ko hyperelastic law only defined " 
  731                 "on dimension 3, sorry");
 
  733     gmm::scale(C, scalar_type(2));
 
  734     gmm::add(gmm::identity_matrix(), C);
 
  735     compute_invariants ci(C);
 
  737     scalar_type z = a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
 
  738       + c*ci.i2() / ci.i3() + d;
 
  739     scalar_type nz = n * pow(z, n-1.);
 
  740     scalar_type di1 = nz * a;
 
  741     scalar_type di2 = nz * c / ci.i3();
 
  742     scalar_type di3 = nz *
 
  743       (b / (2. * sqrt(gmm::abs(ci.i3()))) - c * ci.i2() / gmm::sqr(ci.i3()));
 
  745     gmm::copy(gmm::scaled(ci.grad_i1(), di1 * 2.0), result);
 
  746     gmm::add(gmm::scaled(ci.grad_i2(), di2 * 2.0), result);
 
  747     gmm::add(gmm::scaled(ci.grad_i3(), di3 * 2.0), result);
 
  748     if (det_trans <= scalar_type(0))
 
  749       gmm::add(gmm::scaled(C, 1e200), result);
 
  753   void generalized_Blatz_Ko_hyperelastic_law::grad_sigma
 
  754   (
const base_matrix &E, base_tensor &result,
 
  755    const base_vector ¶ms, scalar_type)
 const {
 
  756     scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
 
  757     scalar_type n = params[4];
 
  759     GMM_ASSERT1(N == 3, 
"Generalized Blatz Ko hyperelastic law only defined " 
  760                 "on dimension 3, sorry");
 
  762     gmm::scale(C, scalar_type(2));
 
  763     gmm::add(gmm::identity_matrix(), C);
 
  764     compute_invariants ci(C);
 
  767     scalar_type z = a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
 
  768       + c*ci.i2() / ci.i3() + d;
 
  769     scalar_type nz = n * pow(z, n-1.);
 
  770     scalar_type di1 = nz * a;
 
  771     scalar_type di2 = nz * c / ci.i3();
 
  772     scalar_type y = (b / (2. * sqrt(gmm::abs(ci.i3()))) - c * ci.i2() / gmm::sqr(ci.i3()));
 
  773     scalar_type di3 = nz * y;
 
  775     gmm::copy(gmm::scaled(ci.sym_grad_grad_i1().as_vector(),
 
  776                           scalar_type(4)*di1), result.as_vector());
 
  777     gmm::add(gmm::scaled(ci.sym_grad_grad_i2().as_vector(),
 
  778                          scalar_type(4)*di2), result.as_vector());
 
  779     gmm::add(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
 
  780                          scalar_type(4)*di3), result.as_vector());
 
  782     scalar_type 
nnz = n * (n-1.) * pow(z, n-2.);
 
  784     A(0, 0) = 
nnz * a * a;
 
  785     A(1, 0) = A(0, 1) = 
nnz * a * c / ci.i3();
 
  786     A(2, 0) = A(0, 2) = 
nnz * a * y;
 
  787     A(1, 1) = 
nnz * c * c / gmm::sqr(ci.i3());
 
  788     A(2, 1) = A(1, 2) = 
nnz * y * c / ci.i3() - nz * c / gmm::sqr(ci.i3());
 
  789     A(2, 2) = 
nnz * y * y + nz * (2. * c * ci.i2() / pow(ci.i3(), 3.) - b / (4. * pow(ci.i3(), 1.5)));
 
  791     typedef const base_matrix * pointer_base_matrix__;
 
  792     pointer_base_matrix__ di[3];
 
  793     di[0] = &(ci.grad_i1());
 
  794     di[1] = &(ci.grad_i2());
 
  795     di[2] = &(ci.grad_i3());
 
  803                 result(l1, l2, l3, l4)
 
  804                   += 4. * A(j, k) * (*di[j])(l1, l2) * (*di[k])(l3, l4);
 
  811   generalized_Blatz_Ko_hyperelastic_law::generalized_Blatz_Ko_hyperelastic_law() {
 
  814     V[0] = 1.0;  V[1] = 1.0, V[2] = 1.5; V[3] = -0.5; V[4] = 1.5;
 
  818   scalar_type Ciarlet_Geymonat_hyperelastic_law::strain_energy
 
  819   (
const base_matrix &E, 
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  820     if (det_trans <= scalar_type(0)) 
return 1e200;
 
  822     scalar_type a = params[2];
 
  823     scalar_type b = params[1]/scalar_type(2) - params[2];
 
  824     scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
 
  826     scalar_type d = params[0]/scalar_type(2) + params[1];
 
  827     scalar_type e = -(scalar_type(3)*(a+b) + c);
 
  829     gmm::copy(gmm::scaled(E, scalar_type(2)), C);
 
  830     gmm::add(gmm::identity_matrix(), C);
 
  831     scalar_type det = bgeot::lu_det(&(*(C.begin())), N);
 
  833       + b * (gmm::sqr(gmm::mat_trace(C)) -
 
  835       + c * det - d * log(det) / scalar_type(2) + e;
 
  838   void Ciarlet_Geymonat_hyperelastic_law::sigma
 
  839   (
const base_matrix &E, base_matrix &result, 
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  841     scalar_type a = params[2];
 
  842     scalar_type b = params[1]/scalar_type(2) - params[2];
 
  843     scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
 
  845     scalar_type d = params[0]/scalar_type(2) + params[1];
 
  847     if (a > params[1]/scalar_type(2)
 
  848         || a < params[1]/scalar_type(2) - params[0]/scalar_type(4) || a < 0)
 
  849       GMM_WARNING1(
"Inconsistent third parameter for Ciarlet-Geymonat " 
  851     gmm::copy(gmm::scaled(E, scalar_type(2)), C);
 
  852     gmm::add(gmm::identity_matrix(), C);
 
  853     gmm::copy(gmm::identity_matrix(), result);
 
  854     gmm::scale(result, scalar_type(2) * (a + b * gmm::mat_trace(C)));
 
  855     gmm::add(gmm::scaled(C, -scalar_type(2) * b), result);
 
  856     if (det_trans <= scalar_type(0))
 
  857       gmm::add(gmm::scaled(C, 1e200), result);
 
  859       scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
 
  860       gmm::add(gmm::scaled(C, scalar_type(2) * c * det - d), result);
 
  864   void Ciarlet_Geymonat_hyperelastic_law::grad_sigma
 
  865   (
const base_matrix &E, base_tensor &result,
const base_vector ¶ms, scalar_type)
 const {
 
  868     scalar_type b2 = params[1] - params[2]*scalar_type(2); 
 
  869     scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
 
  871     scalar_type d = params[0]/scalar_type(2) + params[1];
 
  873     gmm::copy(gmm::scaled(E, scalar_type(2)), C);
 
  874     gmm::add(gmm::identity_matrix(), C);
 
  875     scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
 
  876     std::fill(result.begin(), result.end(), scalar_type(0));
 
  879         result(i, i, j, j) += 2*b2;
 
  880         result(i, j, i, j) -= b2;
 
  881         result(i, j, j, i) -= b2;
 
  884             result(i, j, k, l) +=
 
  885               (C(i, k)*C(l, j) + C(i, l)*C(k, j)) * (d-scalar_type(2)*det*c)
 
  886               + (C(i, j) * C(k, l)) * det*c*scalar_type(4);
 
  894   int levi_civita(
int i, 
int j, 
int k) {
 
  898     return static_cast<int> 
  899       (int(- 1)*(
static_cast<int>(pow(
double(ii-jj),2.))%3)
 
  900        * (
static_cast<int> (pow(
double(ii-kk),
double(2)))%3 )
 
  901        * (
static_cast<int> (pow(
double(jj-kk),
double(2)))%3)
 
  902        * (pow(
double(jj-(ii%3))-
double(0.5),
double(2))-double(1.25)));
 
  907   scalar_type plane_strain_hyperelastic_law::strain_energy
 
  908   (
const base_matrix &E, 
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  909     GMM_ASSERT1(gmm::mat_nrows(E) == 2, 
"Plane strain law is for 2D only.");
 
  910     base_matrix E3D(3,3);
 
  911     E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
 
  912     return pl->strain_energy(E3D, params, det_trans);
 
  915   void plane_strain_hyperelastic_law::sigma
 
  916   (
const base_matrix &E, base_matrix &result,
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  917     GMM_ASSERT1(gmm::mat_nrows(E) == 2, 
"Plane strain law is for 2D only.");
 
  918     base_matrix E3D(3,3), result3D(3,3);
 
  919     E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
 
  920     pl->sigma(E3D, result3D, params, det_trans);
 
  921     result(0,0) = result3D(0,0); result(1,0) = result3D(1,0);
 
  922     result(0,1) = result3D(0,1); result(1,1) = result3D(1,1);
 
  925   void plane_strain_hyperelastic_law::grad_sigma
 
  926   (
const base_matrix &E, base_tensor &result,
const base_vector ¶ms, scalar_type det_trans)
 const {
 
  927     GMM_ASSERT1(gmm::mat_nrows(E) == 2, 
"Plane strain law is for 2D only.");
 
  928     base_matrix E3D(3,3);
 
  929     base_tensor result3D(3,3,3,3);
 
  930     E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
 
  931     pl->grad_sigma(E3D, result3D, params, det_trans);
 
  932     result(0,0,0,0) = result3D(0,0,0,0); result(1,0,0,0) = result3D(1,0,0,0);
 
  933     result(0,1,0,0) = result3D(0,1,0,0); result(1,1,0,0) = result3D(1,1,0,0);
 
  934     result(0,0,1,0) = result3D(0,0,1,0); result(1,0,1,0) = result3D(1,0,1,0);
 
  935     result(0,1,1,0) = result3D(0,1,1,0); result(1,1,1,0) = result3D(1,1,1,0);
 
  936     result(0,0,0,1) = result3D(0,0,0,1); result(1,0,0,1) = result3D(1,0,0,1);
 
  937     result(0,1,0,1) = result3D(0,1,0,1); result(1,1,0,1) = result3D(1,1,0,1);
 
  938     result(0,0,1,1) = result3D(0,0,1,1); result(1,0,1,1) = result3D(1,0,1,1);
 
  939     result(0,1,1,1) = result3D(0,1,1,1); result(1,1,1,1) = result3D(1,1,1,1);
 
  954   struct nonlinear_elasticity_brick : 
public virtual_brick {
 
  956     phyperelastic_law AHL;
 
  958     virtual void asm_real_tangent_terms(
const model &md, 
size_type ,
 
  959                                         const model::varnamelist &vl,
 
  960                                         const model::varnamelist &dl,
 
  961                                         const model::mimlist &mims,
 
  962                                         model::real_matlist &matl,
 
  963                                         model::real_veclist &vecl,
 
  964                                         model::real_veclist &,
 
  966                                         build_version version)
 const {
 
  967       GMM_ASSERT1(mims.size() == 1,
 
  968                   "Nonlinear elasticity brick need a single mesh_im");
 
  969       GMM_ASSERT1(vl.size() == 1,
 
  970                   "Nonlinear elasticity brick need a single variable");
 
  971       GMM_ASSERT1(dl.size() == 1,
 
  972                   "Wrong number of data for nonlinear elasticity brick, " 
  973                   << dl.size() << 
" should be 1 (vector).");
 
  974       GMM_ASSERT1(matl.size() == 1,  
"Wrong number of terms for nonlinear " 
  977       const model_real_plain_vector &u = md.real_variable(vl[0]);
 
  978       const mesh_fem &mf_u = *(md.pmesh_fem_of_variable(vl[0]));
 
  980       const mesh_fem *mf_params = md.pmesh_fem_of_variable(dl[0]);
 
  981       const model_real_plain_vector ¶ms = md.real_variable(dl[0]);
 
  982       const mesh_im &mim = *mims[0];
 
  985       if (mf_params) sl = sl * mf_params->get_qdim() / mf_params->nb_dof();
 
  986       GMM_ASSERT1(sl == AHL->nb_params(), 
"Wrong number of coefficients for the " 
  987                   "nonlinear constitutive elastic law");
 
  989       mesh_region rg(region);
 
  990       mf_u.linked_mesh().intersect_with_mpi_region(rg);
 
  992       if (version & model::BUILD_MATRIX) {
 
  994         GMM_TRACE2(
"Nonlinear elasticity stiffness matrix assembly");
 
  996           (matl[0], mim, mf_u, u, mf_params, params, *AHL, rg);
 
 1000       if (version & model::BUILD_RHS) {
 
 1001         asm_nonlinear_elasticity_rhs(vecl[0], mim,
 
 1002                                      mf_u, u, mf_params, params, *AHL, rg);
 
 1003         gmm::scale(vecl[0], scalar_type(-1));
 
 1008     nonlinear_elasticity_brick(
const phyperelastic_law &AHL_)
 
 1010       set_flags(
"Nonlinear elasticity brick", 
false ,
 
 1023   (
model &md, 
const mesh_im &mim, 
const std::string &varname,
 
 1024    const phyperelastic_law &AHL, 
const std::string &dataname,
 
 1026     pbrick pbr = std::make_shared<nonlinear_elasticity_brick>(AHL);
 
 1029     tl.push_back(model::term_description(varname, varname, 
true));
 
 1030     model::varnamelist dl(1, dataname);
 
 1031     model::varnamelist vl(1, varname);
 
 1032     return md.
add_brick(pbr, vl, dl, tl, model::mimlist(1,&mim), region);
 
 1039   void compute_Von_Mises_or_Tresca(model &md,
 
 1040                                    const std::string &varname,
 
 1041                                    const phyperelastic_law &AHL,
 
 1042                                    const std::string &dataname,
 
 1043                                    const mesh_fem &mf_vm,
 
 1044                                    model_real_plain_vector &VM,
 
 1046     GMM_ASSERT1(gmm::vect_size(VM) == mf_vm.nb_dof(),
 
 1047                 "The vector has not the good size");
 
 1048     const mesh_fem &mf_u = md.mesh_fem_of_variable(varname);
 
 1049     const model_real_plain_vector &u = md.real_variable(varname);
 
 1050     const mesh_fem *mf_params = md.pmesh_fem_of_variable(dataname);
 
 1051     const model_real_plain_vector ¶ms = md.real_variable(dataname);
 
 1054     if (mf_params) sl = sl * mf_params->get_qdim() / mf_params->nb_dof();
 
 1055     GMM_ASSERT1(sl == AHL->nb_params(), 
"Wrong number of coefficients for " 
 1056                 "the nonlinear constitutive elastic law");
 
 1058     unsigned N = unsigned(mf_u.linked_mesh().dim());
 
 1059     unsigned NP = unsigned(AHL->nb_params()), NFem = mf_u.get_qdim();
 
 1060     model_real_plain_vector GRAD(mf_vm.nb_dof()*NFem*N);
 
 1061     model_real_plain_vector PARAMS(mf_vm.nb_dof()*NP);
 
 1062     if (mf_params) 
interpolation(*mf_params, mf_vm, params, PARAMS);
 
 1064     base_matrix E(N, N), gradphi(NFem,N),gradphit(N,NFem), Id(N, N),
 
 1065       sigmahathat(N,N),aux(NFem,N), sigma(NFem,NFem),
 
 1068     if (!mf_params) gmm::copy(params, p);
 
 1069     base_vector eig(NFem);
 
 1070     base_vector ez(NFem);    
 
 1072     gmm::copy(gmm::identity_matrix(), Id);
 
 1073     gmm::copy(gmm::identity_matrix(), IdNFem);
 
 1074     for (
size_type i = 0; i < mf_vm.nb_dof(); ++i) {
 
 1075       gmm::resize(gradphi,NFem,N);
 
 1076       std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
 
 1078       gmm::copy(gmm::transposed(gradphit),gradphi);
 
 1079       for (
unsigned int alpha = 0; alpha <N; ++alpha)
 
 1080         gradphi(alpha, alpha)+=1;
 
 1081       gmm::mult(gmm::transposed(gradphi), gradphi, E);
 
 1082       gmm::add(gmm::scaled(Id, -scalar_type(1)), E);
 
 1083       gmm::scale(E, scalar_type(1)/scalar_type(2));
 
 1085         gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*NP,NP)), p);
 
 1086       AHL->sigma(E, sigmahathat, p, scalar_type(1));
 
 1087       if (NFem == 3 && N == 2) {
 
 1089         for (
unsigned int l = 0; l <NFem; ++l)  {
 
 1091           for (
unsigned int m = 0; m <NFem; ++m)
 
 1092             for (
unsigned int n = 0; n <NFem; ++n){
 
 1093               ez[l]+=levi_civita(l,m,n)*gradphi(m,0)*gradphi(n,1);
 
 1100       gmm::mult(aux, gmm::transposed(gradphi), sigma);
 
 1104       if (NFem == 3 && N == 2) {
 
 1106         for (
unsigned int ll = 0; ll <NFem; ++ll)
 
 1107           for (
unsigned int ii = 0; ii <NFem; ++ii)
 
 1108             for (
unsigned int jj = 0; jj <NFem; ++jj)
 
 1109               gradphi(ll,2)+=(levi_civita(ll,ii,jj)*gradphi(ii,0)
 
 1110                               *gradphi(jj,1))/normEz;
 
 1114       gmm::scale(sigma, scalar_type(1) / bgeot::lu_det(&(*(gradphi.begin())), NFem));
 
 1118         gmm::add(gmm::scaled(IdNFem, -gmm::mat_trace(sigma) / NFem), sigma);
 
 1125         gmm::symmetric_qr_algorithm(sigma, eig);
 
 1126         std::sort(eig.begin(), eig.end());
 
 1127         VM[i] = eig.back() - eig.front();
 
 1133   void compute_sigmahathat(model &md,
 
 1134                            const std::string &varname,
 
 1135                            const phyperelastic_law &AHL,
 
 1136                            const std::string &dataname,
 
 1137                            const mesh_fem &mf_sigma,
 
 1138                            model_real_plain_vector &SIGMA) {
 
 1139     const mesh_fem &mf_u = md.mesh_fem_of_variable(varname);
 
 1140     const model_real_plain_vector &u = md.real_variable(varname);
 
 1141     const mesh_fem *mf_params = md.pmesh_fem_of_variable(dataname);
 
 1142     const model_real_plain_vector ¶ms = md.real_variable(dataname);
 
 1145     if (mf_params) sl = sl * mf_params->get_qdim() / mf_params->nb_dof();
 
 1146     GMM_ASSERT1(sl == AHL->nb_params(), 
"Wrong number of coefficients for " 
 1147                 "the nonlinear constitutive elastic law");
 
 1149     unsigned N = unsigned(mf_u.linked_mesh().dim());
 
 1150     unsigned NP = unsigned(AHL->nb_params()), NFem = mf_u.get_qdim();
 
 1151     GMM_ASSERT1(mf_sigma.nb_dof() > 0, 
"Bad mf_sigma");
 
 1155     GMM_ASSERT1(((ratio == 1) || (ratio == N*N)) &&
 
 1156                 (gmm::vect_size(SIGMA) == mf_sigma.nb_dof()*ratio),
 
 1157                 "The vector has not the good size");
 
 1159     model_real_plain_vector GRAD(mf_sigma.nb_dof()*ratio*NFem/N);
 
 1160     model_real_plain_vector PARAMS(mf_sigma.nb_dof()*NP);
 
 1163     getfem::mesh_trans_inv mti(mf_sigma.linked_mesh());
 
 1165       for (
size_type i = 0; i < mf_sigma.nb_dof(); ++i)
 
 1166         mti.add_point(mf_sigma.point_of_basic_dof(i));
 
 1171     base_matrix E(N, N), gradphi(NFem,N),gradphit(N,NFem), Id(N, N),
 
 1172       sigmahathat(N,N),aux(NFem,N), sigma(NFem,NFem),
 
 1178     base_vector eig(NFem);
 
 1179     base_vector ez(NFem); 
 
 1181     gmm::copy(gmm::identity_matrix(), IdNFem);
 
 1186     for (
size_type i = 0; i < mf_sigma.nb_dof()/qqdim; ++i) {
 
 1192       std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
 
 1195       gmm::copy(gmm::transposed(gradphit),gradphi);
 
 1196       for (
unsigned int alpha = 0; 
alpha <N; ++
alpha)
 
 1197         gradphi(alpha, alpha) += scalar_type(1);
 
 1198       gmm::mult(gmm::transposed(gradphi), gradphi, E);
 
 1199       gmm::add(gmm::scaled(Id, -scalar_type(1)), E);
 
 1200       gmm::scale(E, scalar_type(1)/scalar_type(2));
 
 1202         gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*ratio*NP,NP)),p);
 
 1204       AHL->sigma(E, sigmahathat, p, scalar_type(1));
 
 1206       std::copy(sigmahathat.begin(), sigmahathat.end(), SIGMA.begin()+i*N*N);
 
 1218   struct nonlinear_incompressibility_brick : 
public virtual_brick {
 
 1220     virtual void asm_real_tangent_terms(
const model &md, 
size_type,
 
 1221                                         const model::varnamelist &vl,
 
 1222                                         const model::varnamelist &dl,
 
 1223                                         const model::mimlist &mims,
 
 1224                                         model::real_matlist &matl,
 
 1225                                         model::real_veclist &vecl,
 
 1226                                         model::real_veclist &veclsym,
 
 1228                                         build_version version)
 const {
 
 1230       GMM_ASSERT1(matl.size() == 2,  
"Wrong number of terms for nonlinear " 
 1231                   "incompressibility brick");
 
 1232       GMM_ASSERT1(dl.size() == 0, 
"Nonlinear incompressibility brick need no " 
 1234       GMM_ASSERT1(mims.size() == 1, 
"Nonlinear incompressibility brick need a " 
 1236       GMM_ASSERT1(vl.size() == 2, 
"Wrong number of variables for nonlinear " 
 1237                   "incompressibility brick");
 
 1239       const mesh_fem &mf_u = md.mesh_fem_of_variable(vl[0]);
 
 1240       const mesh_fem &mf_p = md.mesh_fem_of_variable(vl[1]);
 
 1241       const model_real_plain_vector &u = md.real_variable(vl[0]);
 
 1242       const model_real_plain_vector &p = md.real_variable(vl[1]);
 
 1243       const mesh_im &mim = *mims[0];
 
 1244       mesh_region rg(region);
 
 1245       mim.linked_mesh().intersect_with_mpi_region(rg);
 
 1247       if (version & model::BUILD_MATRIX) {
 
 1250         asm_nonlinear_incomp_tangent_matrix(matl[0], matl[1],
 
 1251                                             mim, mf_u, mf_p, u, p, rg);
 
 1254       if (version & model::BUILD_RHS) {
 
 1255         asm_nonlinear_incomp_rhs(vecl[0], veclsym[1], mim, mf_u, mf_p,u,p, rg);
 
 1256         gmm::scale(vecl[0], scalar_type(-1));
 
 1257         gmm::scale(veclsym[1], scalar_type(-1));
 
 1262     nonlinear_incompressibility_brick() {
 
 1263       set_flags(
"Nonlinear incompressibility brick",
 
 1273   (
model &md, 
const mesh_im &mim, 
const std::string &varname,
 
 1274    const std::string &multname, 
size_type region) {
 
 1275     pbrick pbr = std::make_shared<nonlinear_incompressibility_brick>();
 
 1277     tl.push_back(model::term_description(varname, varname, 
true));
 
 1278     tl.push_back(model::term_description(varname, multname, 
true));
 
 1279     model::varnamelist vl(1, varname);
 
 1280     vl.push_back(multname);
 
 1281     model::varnamelist dl;
 
 1282     return md.
add_brick(pbr, vl, dl, tl, model::mimlist(1, &mim), region);
 
 1296   static void ga_init_scalar_(bgeot::multi_index &mi) { mi.resize(0); }
 
 1297   static void ga_init_square_matrix_(bgeot::multi_index &mi, 
size_type N)
 
 1298   { mi.resize(2); mi[0] = mi[1] = N; }
 
 1303   struct matrix_i2_operator : 
public ga_nonlinear_operator {
 
 1304     bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
 const {
 
 1305       if (args.size() != 1 || args[0]->sizes().size() != 2
 
 1306           || args[0]->sizes()[0] != args[0]->sizes()[1]) 
return false;
 
 1307       ga_init_scalar_(sizes);
 
 1312     void value(
const arg_list &args, base_tensor &result)
 const {
 
 1314       const base_tensor &t = *args[0];
 
 1315       scalar_type tr = scalar_type(0);
 
 1316       for (
size_type i = 0; i < N; ++i) tr += t[i*N+i];
 
 1317       scalar_type tr2 = scalar_type(0);
 
 1320           tr2 += t[i+ j*N] * t[j + i*N];
 
 1321       result[0] = (tr*tr-tr2)/scalar_type(2);
 
 1325     void derivative(
const arg_list &args, 
size_type,
 
 1326                     base_tensor &result)
 const {
 
 1328       const base_tensor &t = *args[0];
 
 1329       scalar_type tr = scalar_type(0);
 
 1330       for (
size_type i = 0; i < N; ++i) tr += t[i*N+i];
 
 1331       base_tensor::iterator it = result.begin();
 
 1334           *it = ((i == j) ? tr : scalar_type(0)) - t[i*N+j];
 
 1335       GMM_ASSERT1(it == result.end(), 
"Internal error");
 
 1340                            base_tensor &result)
 const {
 
 1345            result[(N+1)*(i+N*N*j)] += scalar_type(1);
 
 1346            result[(N+1)*N*j + i*(N*N*N + 1)] -= scalar_type(1);
 
 1353   struct matrix_j1_operator : 
public ga_nonlinear_operator {
 
 1354     bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
 const {
 
 1355       if (args.size() != 1 || args[0]->sizes().size() != 2
 
 1356           || args[0]->sizes()[0] != args[0]->sizes()[1]) 
return false;
 
 1357       ga_init_scalar_(sizes);
 
 1362     void value(
const arg_list &args, base_tensor &result)
 const {
 
 1364       base_matrix M(N, N);
 
 1365       gmm::copy(args[0]->as_vector(), M.as_vector());
 
 1366       scalar_type det = bgeot::lu_det(&(*(M.begin())), N);
 
 1367       scalar_type tr = scalar_type(0);
 
 1368       for (
size_type i = 0; i < N; ++i) tr += M(i,i);
 
 1370         result[0] = tr / pow(det, scalar_type(1)/scalar_type(3));
 
 1376     void derivative(
const arg_list &args, 
size_type,
 
 1377                     base_tensor &result)
 const {
 
 1379       base_matrix M(N, N);
 
 1380       gmm::copy(args[0]->as_vector(), M.as_vector());
 
 1381       scalar_type tr = scalar_type(0);
 
 1382       for (
size_type i = 0; i < N; ++i) tr += M(i,i);
 
 1383       scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
 
 1385         base_tensor::iterator it = result.begin();
 
 1388             *it = (((i == j) ? scalar_type(1) : scalar_type(0))
 
 1389                    - tr*M(j,i)/scalar_type(3))
 
 1390               / pow(det, scalar_type(1)/scalar_type(3));
 
 1391         GMM_ASSERT1(it == result.end(), 
"Internal error");
 
 1393         std::fill(result.begin(), result.end(), 1.E200);
 
 1399                            base_tensor &result)
 const {
 
 1401       base_matrix M(N, N);
 
 1402       gmm::copy(args[0]->as_vector(), M.as_vector());
 
 1403       scalar_type tr = scalar_type(0);
 
 1404       for (
size_type i = 0; i < N; ++i) tr += M(i,i);
 
 1405       scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
 
 1407         base_tensor::iterator it = result.begin();
 
 1412                 *it = (- ((k == l) ? M(j, i) : scalar_type(0))
 
 1414                        - ((i == j) ? M(l, k) : scalar_type(0))
 
 1415                        + tr*M(j,i)*M(k,l)/ scalar_type(3))
 
 1416                   / (scalar_type(3)*pow(det, scalar_type(1)/scalar_type(3)));
 
 1417         GMM_ASSERT1(it == result.end(), 
"Internal error");
 
 1419         std::fill(result.begin(), result.end(), 1.E200);
 
 1425   struct matrix_j2_operator : 
public ga_nonlinear_operator {
 
 1426     bool result_size(
const arg_list &args, bgeot::multi_index &sizes)
 const {
 
 1427       if (args.size() != 1 || args[0]->sizes().size() != 2
 
 1428           || args[0]->sizes()[0] != args[0]->sizes()[1]) 
return false;
 
 1429       ga_init_scalar_(sizes);
 
 1434     void value(
const arg_list &args, base_tensor &result)
 const {
 
 1436       base_matrix M(N, N);
 
 1437       gmm::copy(args[0]->as_vector(), M.as_vector());
 
 1438       scalar_type tr = scalar_type(0);
 
 1439       for (
size_type i = 0; i < N; ++i) tr += M(i,i);
 
 1440       scalar_type tr2 = scalar_type(0);
 
 1443           tr2 += M(i,j)*M(j,i);
 
 1444       scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
 
 1445       scalar_type det = bgeot::lu_det(&(*(M.begin())), N);
 
 1447         result[0] = i2 / pow(det, scalar_type(2)/scalar_type(3));
 
 1453     void derivative(
const arg_list &args, 
size_type,
 
 1454                     base_tensor &result)
 const {
 
 1456       const base_tensor &t = *args[0];
 
 1457       base_matrix M(N, N);
 
 1458       gmm::copy(t.as_vector(), M.as_vector());
 
 1459       scalar_type tr = scalar_type(0);
 
 1460       for (
size_type i = 0; i < N; ++i) tr += M(i,i);
 
 1461       scalar_type tr2 = scalar_type(0);
 
 1464           tr2 += M(i,j)*M(j,i);
 
 1465       scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
 
 1466       scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
 
 1467       base_tensor::iterator it = result.begin();
 
 1470           *it = (((i == j) ? tr : scalar_type(0)) - t[j+N*i]
 
 1471                  - scalar_type(2)*i2*M(j,i)/(scalar_type(3)))
 
 1472             / pow(det, scalar_type(2)/scalar_type(3));
 
 1473       GMM_ASSERT1(it == result.end(), 
"Internal error");
 
 1478                            base_tensor &result)
 const {
 
 1480       const base_tensor &t = *args[0];
 
 1481       base_matrix M(N, N);
 
 1482       gmm::copy(t.as_vector(), M.as_vector());
 
 1483       scalar_type tr = scalar_type(0);
 
 1484       for (
size_type i = 0; i < N; ++i) tr += M(i,i);
 
 1485       scalar_type tr2 = scalar_type(0);
 
 1488           tr2 += M(i,j)*M(j,i);
 
 1489       scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
 
 1490       scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
 
 1491       base_tensor::iterator it = result.begin();
 
 1496               *it = ( + (((i==j) && (k==l)) ? 1. : 0.)
 
 1497                       - (((i==l) && (k==j)) ? 1. : 0.)
 
 1498                       + 10.*i2*M(j,i)*M(l,k)/(9.)
 
 1499                       - 2.*(M(j,i)*(tr*((k==l) ? 1.:0.)-t[l+N*k]))/(3.)
 
 1500                       - 2.*(M(l,k)*(tr*((i==j) ? 1.:0.)-t[j+N*i]))/(3.)
 
 1501                       - 2.*i2*(M(j,i)*M(l,k)-M(j,k)*M(l,i))/(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);
 
 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);
 
 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());
 
 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());
 
 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());
 
 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());
 
 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);
 
 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());
 
 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());
 
 1953       gmm::scale(E, scalar_type(0.5));
 
 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());
 
 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);
 
 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_PK2",
 
 2049        std::make_shared<Saint_Venant_Kirchhoff_sigma>());
 
 2050     PREDEF_OPERATORS.add_method
 
 2051       (
"Saint_Venant_Kirchhoff_potential",
 
 2052        std::make_shared<AHL_wrapper_potential>
 
 2053        (std::make_shared<SaintVenant_Kirchhoff_hyperelastic_law>()));
 
 2054     PREDEF_OPERATORS.add_method
 
 2055       (
"Plane_Strain_Saint_Venant_Kirchhoff_sigma",
 
 2056        std::make_shared<Saint_Venant_Kirchhoff_sigma>());
 
 2057     PREDEF_OPERATORS.add_method
 
 2058       (
"Plane_Strain_Saint_Venant_Kirchhoff_PK2",
 
 2059        std::make_shared<Saint_Venant_Kirchhoff_sigma>());
 
 2060     PREDEF_OPERATORS.add_method
 
 2061       (
"Plane_Strain_Saint_Venant_Kirchhoff_potential",
 
 2062        std::make_shared<AHL_wrapper_potential>
 
 2063        (std::make_shared<SaintVenant_Kirchhoff_hyperelastic_law>()));
 
 2065     phyperelastic_law gbklaw
 
 2066       = std::make_shared<generalized_Blatz_Ko_hyperelastic_law>();
 
 2067     PREDEF_OPERATORS.add_method
 
 2068       (
"Generalized_Blatz_Ko_sigma",
 
 2069        std::make_shared<AHL_wrapper_sigma>(gbklaw));
 
 2070      PREDEF_OPERATORS.add_method
 
 2071       (
"Generalized_Blatz_Ko_PK2",
 
 2072        std::make_shared<AHL_wrapper_sigma>(gbklaw));
 
 2073    PREDEF_OPERATORS.add_method
 
 2074       (
"Generalized_Blatz_Ko_potential",
 
 2075        std::make_shared<AHL_wrapper_potential>
 
 2076        (std::make_shared<generalized_Blatz_Ko_hyperelastic_law>()));
 
 2077     PREDEF_OPERATORS.add_method
 
 2078       (
"Plane_Strain_Generalized_Blatz_Ko_sigma",
 
 2079        std::make_shared<AHL_wrapper_sigma>
 
 2080        (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
 
 2081     PREDEF_OPERATORS.add_method
 
 2082       (
"Plane_Strain_Generalized_Blatz_Ko_PK2",
 
 2083        std::make_shared<AHL_wrapper_sigma>
 
 2084        (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
 
 2085     PREDEF_OPERATORS.add_method
 
 2086       (
"Plane_Strain_Generalized_Blatz_Ko_potential",
 
 2087        std::make_shared<AHL_wrapper_potential>
 
 2088        (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
 
 2090     phyperelastic_law cigelaw
 
 2091       = std::make_shared<Ciarlet_Geymonat_hyperelastic_law>();
 
 2092     PREDEF_OPERATORS.add_method
 
 2093       (
"Ciarlet_Geymonat_PK2", std::make_shared<AHL_wrapper_sigma>(cigelaw));
 
 2094     PREDEF_OPERATORS.add_method
 
 2095       (
"Ciarlet_Geymonat_sigma", std::make_shared<AHL_wrapper_sigma>(cigelaw));
 
 2096     PREDEF_OPERATORS.add_method
 
 2097       (
"Ciarlet_Geymonat_potential",
 
 2098        std::make_shared<AHL_wrapper_potential>
 
 2099        (std::make_shared<Ciarlet_Geymonat_hyperelastic_law>()));
 
 2100     PREDEF_OPERATORS.add_method
 
 2101       (
"Plane_Strain_Ciarlet_Geymonat_sigma",
 
 2102        std::make_shared<AHL_wrapper_sigma>
 
 2103        (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
 
 2104     PREDEF_OPERATORS.add_method
 
 2105       (
"Plane_Strain_Ciarlet_Geymonat_PK2",
 
 2106        std::make_shared<AHL_wrapper_sigma>
 
 2107        (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
 
 2108     PREDEF_OPERATORS.add_method
 
 2109       (
"Plane_Strain_Ciarlet_Geymonat_potential",
 
 2110        std::make_shared<AHL_wrapper_potential>
 
 2111        (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
 
 2113     phyperelastic_law morilaw
 
 2114       = std::make_shared<Mooney_Rivlin_hyperelastic_law>();
 
 2115     PREDEF_OPERATORS.add_method
 
 2116       (
"Incompressible_Mooney_Rivlin_sigma",
 
 2117        std::make_shared<AHL_wrapper_sigma>(morilaw));
 
 2118     PREDEF_OPERATORS.add_method
 
 2119       (
"Incompressible_Mooney_Rivlin_PK2",
 
 2120        std::make_shared<AHL_wrapper_sigma>(morilaw));
 
 2121     PREDEF_OPERATORS.add_method
 
 2122       (
"Incompressible_Mooney_Rivlin_potential",
 
 2123        std::make_shared<AHL_wrapper_potential>
 
 2124        (std::make_shared<Mooney_Rivlin_hyperelastic_law>()));
 
 2125     PREDEF_OPERATORS.add_method
 
 2126       (
"Plane_Strain_Incompressible_Mooney_Rivlin_PK2",
 
 2127        std::make_shared<AHL_wrapper_sigma>
 
 2128        (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
 
 2129     PREDEF_OPERATORS.add_method
 
 2130       (
"Plane_Strain_Incompressible_Mooney_Rivlin_sigma",
 
 2131        std::make_shared<AHL_wrapper_sigma>
 
 2132        (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
 
 2133     PREDEF_OPERATORS.add_method
 
 2134       (
"Plane_Strain_Incompressible_Mooney_Rivlin_potential",
 
 2135        std::make_shared<AHL_wrapper_potential>
 
 2136        (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
 
 2138     phyperelastic_law cmorilaw
 
 2139       = std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true);
 
 2140     PREDEF_OPERATORS.add_method
 
 2141       (
"Compressible_Mooney_Rivlin_sigma",
 
 2142        std::make_shared<AHL_wrapper_sigma>(cmorilaw));
 
 2143     PREDEF_OPERATORS.add_method
 
 2144       (
"Compressible_Mooney_Rivlin_PK2",
 
 2145        std::make_shared<AHL_wrapper_sigma>(cmorilaw));
 
 2146     PREDEF_OPERATORS.add_method
 
 2147       (
"Compressible_Mooney_Rivlin_potential",
 
 2148        std::make_shared<AHL_wrapper_potential>
 
 2149        (std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true)));
 
 2150     PREDEF_OPERATORS.add_method
 
 2151       (
"Plane_Strain_Compressible_Mooney_Rivlin_PK2",
 
 2152        std::make_shared<AHL_wrapper_sigma>
 
 2153        (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
 
 2154     PREDEF_OPERATORS.add_method
 
 2155       (
"Plane_Strain_Compressible_Mooney_Rivlin_sigma",
 
 2156        std::make_shared<AHL_wrapper_sigma>
 
 2157        (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
 
 2158     PREDEF_OPERATORS.add_method
 
 2159       (
"Plane_Strain_Compressible_Mooney_Rivlin_potential",
 
 2160        std::make_shared<AHL_wrapper_potential>
 
 2161        (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
 
 2163     phyperelastic_law ineolaw
 
 2164       = std::make_shared<Mooney_Rivlin_hyperelastic_law>(
false, 
true);
 
 2165     PREDEF_OPERATORS.add_method
 
 2166       (
"Incompressible_Neo_Hookean_sigma",
 
 2167        std::make_shared<AHL_wrapper_sigma>(ineolaw));
 
 2168     PREDEF_OPERATORS.add_method
 
 2169       (
"Incompressible_Neo_Hookean_PK2",
 
 2170        std::make_shared<AHL_wrapper_sigma>(ineolaw));
 
 2171     PREDEF_OPERATORS.add_method
 
 2172       (
"Incompressible_Neo_Hookean_potential",
 
 2173        std::make_shared<AHL_wrapper_potential>
 
 2174        (std::make_shared<Mooney_Rivlin_hyperelastic_law>(
false,
true)));
 
 2175     PREDEF_OPERATORS.add_method
 
 2176       (
"Plane_Strain_Incompressible_Neo_Hookean_sigma",
 
 2177        std::make_shared<AHL_wrapper_sigma>
 
 2178        (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
 
 2179     PREDEF_OPERATORS.add_method
 
 2180       (
"Plane_Strain_Incompressible_Neo_Hookean_PK2",
 
 2181        std::make_shared<AHL_wrapper_sigma>
 
 2182        (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
 
 2183     PREDEF_OPERATORS.add_method
 
 2184       (
"Plane_Strain_Incompressible_Neo_Hookean_potential",
 
 2185        std::make_shared<AHL_wrapper_potential>
 
 2186        (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
 
 2188     phyperelastic_law cneolaw
 
 2189       = std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true, 
true);
 
 2190     PREDEF_OPERATORS.add_method
 
 2191       (
"Compressible_Neo_Hookean_sigma",
 
 2192        std::make_shared<AHL_wrapper_sigma>(cneolaw));
 
 2193     PREDEF_OPERATORS.add_method
 
 2194       (
"Compressible_Neo_Hookean_PK2",
 
 2195        std::make_shared<AHL_wrapper_sigma>(cneolaw));
 
 2196     PREDEF_OPERATORS.add_method
 
 2197       (
"Compressible_Neo_Hookean_potential",
 
 2198        std::make_shared<AHL_wrapper_potential>
 
 2199        (std::make_shared<Mooney_Rivlin_hyperelastic_law>(
true,
true)));
 
 2200     PREDEF_OPERATORS.add_method
 
 2201       (
"Plane_Strain_Compressible_Neo_Hookean_sigma",
 
 2202        std::make_shared<AHL_wrapper_sigma>
 
 2203        (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
 
 2204     PREDEF_OPERATORS.add_method
 
 2205       (
"Plane_Strain_Compressible_Neo_Hookean_PK2",
 
 2206        std::make_shared<AHL_wrapper_sigma>
 
 2207        (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
 
 2208     PREDEF_OPERATORS.add_method
 
 2209       (
"Plane_Strain_Compressible_Neo_Hookean_potential",
 
 2210        std::make_shared<AHL_wrapper_potential>
 
 2211        (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
 
 2213     phyperelastic_law cneobolaw
 
 2214       = std::make_shared<Neo_Hookean_hyperelastic_law>(
true);
 
 2215     PREDEF_OPERATORS.add_method
 
 2216       (
"Compressible_Neo_Hookean_Bonet_sigma",
 
 2217        std::make_shared<AHL_wrapper_sigma>(cneobolaw));
 
 2218     PREDEF_OPERATORS.add_method
 
 2219       (
"Compressible_Neo_Hookean_Bonet_PK2",
 
 2220        std::make_shared<AHL_wrapper_sigma>(cneobolaw));
 
 2221     PREDEF_OPERATORS.add_method
 
 2222       (
"Compressible_Neo_Hookean_Bonet_potential",
 
 2223        std::make_shared<AHL_wrapper_potential>
 
 2224        (std::make_shared<Neo_Hookean_hyperelastic_law>(
true)));
 
 2225     PREDEF_OPERATORS.add_method
 
 2226       (
"Plane_Strain_Compressible_Neo_Hookean_Bonet_sigma",
 
 2227        std::make_shared<AHL_wrapper_sigma>
 
 2228        (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
 
 2229     PREDEF_OPERATORS.add_method
 
 2230       (
"Plane_Strain_Compressible_Neo_Hookean_Bonet_PK2",
 
 2231        std::make_shared<AHL_wrapper_sigma>
 
 2232        (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
 
 2233     PREDEF_OPERATORS.add_method
 
 2234       (
"Plane_Strain_Compressible_Neo_Hookean_Bonet_potential",
 
 2235        std::make_shared<AHL_wrapper_potential>
 
 2236        (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
 
 2238     phyperelastic_law cneocilaw
 
 2239       = std::make_shared<Neo_Hookean_hyperelastic_law>(
false);
 
 2240     PREDEF_OPERATORS.add_method
 
 2241       (
"Compressible_Neo_Hookean_Ciarlet_sigma",
 
 2242        std::make_shared<AHL_wrapper_sigma>(cneocilaw));
 
 2243     PREDEF_OPERATORS.add_method
 
 2244       (
"Compressible_Neo_Hookean_Ciarlet_PK2",
 
 2245        std::make_shared<AHL_wrapper_sigma>(cneocilaw));
 
 2246     PREDEF_OPERATORS.add_method
 
 2247       (
"Compressible_Neo_Hookean_Ciarlet_potential",
 
 2248        std::make_shared<AHL_wrapper_potential>
 
 2249        (std::make_shared<Neo_Hookean_hyperelastic_law>(
false)));
 
 2250     PREDEF_OPERATORS.add_method
 
 2251       (
"Plane_Strain_Compressible_Neo_Hookean_Ciarlet_sigma",
 
 2252        std::make_shared<AHL_wrapper_sigma>
 
 2253        (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
 
 2254     PREDEF_OPERATORS.add_method
 
 2255       (
"Plane_Strain_Compressible_Neo_Hookean_Ciarlet_PK2",
 
 2256        std::make_shared<AHL_wrapper_sigma>
 
 2257        (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
 
 2258     PREDEF_OPERATORS.add_method
 
 2259       (
"Plane_Strain_Compressible_Neo_Hookean_Ciarlet_potential",
 
 2260        std::make_shared<AHL_wrapper_potential>
 
 2261        (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
 
 2267   bool predef_operators_nonlinear_elasticity_initialized
 
 2268     = init_predef_operators();
 
 2271   std::string adapt_law_name(
const std::string &lawname, 
size_type N) {
 
 2272     std::string adapted_lawname = lawname;
 
 2274     for (
size_type i = 0; i < lawname.size(); ++i)
 
 2275       if (adapted_lawname[i] == 
' ') adapted_lawname[i] = 
'_';
 
 2277     if (adapted_lawname.compare(
"SaintVenant_Kirchhoff") == 0) {
 
 2278       adapted_lawname = 
"Saint_Venant_Kirchhoff";
 
 2279     } 
else if (adapted_lawname.compare(
"Saint_Venant_Kirchhoff") == 0) {
 
 2281     } 
else if (adapted_lawname.compare(
"Generalized_Blatz_Ko") == 0) {
 
 2282       if (N == 2) adapted_lawname = 
"Plane_Strain_" + adapted_lawname;
 
 2283     } 
else if (adapted_lawname.compare(
"Ciarlet_Geymonat") == 0) {
 
 2284       if (N == 2) adapted_lawname = 
"Plane_Strain_" + adapted_lawname;
 
 2285     } 
else if (adapted_lawname.compare(
"Incompressible_Mooney_Rivlin") == 0) {
 
 2286       if (N == 2) adapted_lawname = 
"Plane_Strain_" + adapted_lawname;
 
 2287     } 
else if (adapted_lawname.compare(
"Compressible_Mooney_Rivlin") == 0) {
 
 2288       if (N == 2) adapted_lawname = 
"Plane_Strain_" + adapted_lawname;
 
 2289     } 
else if (adapted_lawname.compare(
"Incompressible_Neo_Hookean") == 0) {
 
 2290       if (N == 2) adapted_lawname = 
"Plane_Strain_" + adapted_lawname;
 
 2291     } 
else if (adapted_lawname.compare(
"Compressible_Neo_Hookean") == 0 ||
 
 2292                adapted_lawname.compare(
"Compressible_Neo_Hookean_Bonet") == 0 ||
 
 2293                adapted_lawname.compare(
"Compressible_Neo_Hookean_Ciarlet") == 0 ) {
 
 2294       if (N == 2) adapted_lawname = 
"Plane_Strain_" + adapted_lawname;
 
 2296       GMM_ASSERT1(
false, lawname << 
" is not a known hyperelastic law");
 
 2298     return adapted_lawname;
 
 2303   (
model &md, 
const mesh_im &mim, 
const std::string &lawname,
 
 2304    const std::string &varname, 
const std::string ¶ms,
 
 2306     std::string test_varname = 
"Test_" + sup_previous_and_dot_to_varname(varname);
 
 2308     GMM_ASSERT1(N >= 2 && N <= 3,
 
 2309                 "Finite strain elasticity brick works only in 2D or 3D");
 
 2312     GMM_ASSERT1(mf, 
"Finite strain elasticity brick can only be applied on " 
 2315     GMM_ASSERT1(Q == N, 
"Finite strain elasticity brick can only be applied " 
 2316                 "on a fem variable having the same dimension as the mesh");
 
 2318     std::string adapted_lawname = adapt_law_name(lawname, N);
 
 2320     std::string expr = 
"((Id(meshdim)+Grad_"+varname+
")*(" + adapted_lawname
 
 2321       + 
"_PK2(Grad_"+varname+
","+params+
"))):Grad_" + test_varname;
 
 2324       (md, mim, expr, region, 
true, 
false,
 
 2325        "Finite strain elasticity brick for " + adapted_lawname + 
" law");
 
 2329   (
model &md, 
const mesh_im &mim, 
const std::string &varname,
 
 2330    const std::string &multname, 
size_type region) {
 
 2331     std::string test_varname = 
"Test_" + sup_previous_and_dot_to_varname(varname);
 
 2332     std::string test_multname = 
"Test_" + sup_previous_and_dot_to_varname(multname);
 
 2335       = 
"(" + test_multname+ 
")*(1-Det(Id(meshdim)+Grad_" + varname + 
"))" 
 2336       + 
"-(" + multname + 
")*(Det(Id(meshdim)+Grad_" + varname + 
")" 
 2337       + 
"*((Inv(Id(meshdim)+Grad_" + varname + 
"))':Grad_" 
 2338       + test_varname + 
"))" ;
 
 2340       (md, mim, expr, region, 
true, 
false,
 
 2341        "Finite strain incompressibility brick");
 
 2345     (
model &md, 
const std::string &lawname, 
const std::string &varname,
 
 2346      const std::string ¶ms, 
const mesh_fem &mf_vm,
 
 2347      model_real_plain_vector &VM, 
const mesh_region &rg) {
 
 2350     std::string adapted_lawname = adapt_law_name(lawname, N);
 
 2352     std::string expr = 
"sqrt(3/2)*Norm(Deviator(Cauchy_stress_from_PK2(" 
 2353       + adapted_lawname + 
"_PK2(Grad_" + varname + 
"," + params + 
"),Grad_" 
 2355     ga_interpolation_Lagrange_fem(md, expr, mf_vm, VM, rg);
 
static T & instance()
Instance from the current thread.
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
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)
Describe a finite element method linked to a mesh.
virtual dim_type get_qdim() const
Return the Q dimension.
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
Describe an integration method linked to a mesh.
const mesh & linked_mesh() const
Give a reference to the linked mesh of type mesh.
structure used to hold a set of convexes and/or convex faces.
`‘Model’' variables store the variables, the data and the description of a model.
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.
const mesh_fem * pmesh_fem_of_variable(const std::string &name) const
Gives a pointer to the mesh_fem of a variable if any.
A language for generic assembly of pde boundary value problems.
Model representation in Getfem.
Non-linear elasticty and incompressibility bricks.
size_type nnz(const L &l)
count the number of non-zero entries of a vector or matrix.
void copy(const L1 &l1, L2 &l2)
*/
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
void fill_random(L &l)
fill a vector or matrix with random value (uniform [-1,1]).
linalg_traits< M >::value_type mat_trace(const M &m)
Trace of a matrix.
number_traits< typename linalg_traits< M >::value_type >::magnitude_type mat_euclidean_norm(const M &m)
Euclidean norm of a matrix.
void clear(L &l)
clear (fill with zeros) a vector or matrix.
void resize(V &v, size_type n)
*/
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*/
void add(const L1 &l1, L2 &l2)
*/
number_traits< typename linalg_traits< M >::value_type >::magnitude_type mat_euclidean_norm_sqr(const M &m)
*/
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.
size_t size_type
used as the common size type in the library
size_type alpha(short_type n, short_type d)
Return the value of  which is the number of monomials of a polynomial of  variables and degree .
GEneric Tool for Finite Element Methods.
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.
size_type APIDECL add_nonlinear_term(model &md, const mesh_im &mim, const std::string &expr, size_type region=size_type(-1), bool is_sym=false, bool is_coercive=false, const std::string &brickname="")
Add a nonlinear term given by the weak form language expression expr which will be assembled in regio...
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...
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...
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...
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...
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...
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.
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