GetFEM  5.4.3
getfem_nonlinear_elasticity.cc
1 /*===========================================================================
2 
3  Copyright (C) 2000-2020 Yves Renard
4 
5  This file is a part of GetFEM
6 
7  GetFEM is free software; you can redistribute it and/or modify it
8  under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation; either version 3 of the License, or
10  (at your option) any later version along with the GCC Runtime Library
11  Exception either version 3.1 or (at your option) any later version.
12  This program is distributed in the hope that it will be useful, but
13  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  License and GCC Runtime Library Exception for more details.
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program; if not, write to the Free Software Foundation,
18  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
19 
20 ===========================================================================*/
21 
22 
23 #include "getfem/getfem_models.h"
26 
27 namespace getfem {
28 
29 
30  /* Useful functions to compute the invariants and their derivatives
31  Note that the second derivative is symmetrized (see the user
32  documentation for more details). The matrix E is assumed to be symmetric.
33  */
34 
35 
36  static scalar_type frobenius_product_trans(const base_matrix &A,
37  const base_matrix &B) {
38  size_type N = gmm::mat_nrows(A);
39  scalar_type res = scalar_type(0);
40  for (size_type i = 0; i < N; ++i)
41  for (size_type j = 0; j < N; ++j)
42  res += A(i, j) * B(j, i);
43  return res;
44  }
45 
46  struct compute_invariants {
47 
48  const base_matrix &E;
49  base_matrix Einv;
50  size_type N;
51  scalar_type i1_, i2_, i3_, j1_, j2_;
52  bool i1_c, i2_c, i3_c, j1_c, j2_c;
53 
54  base_matrix di1, di2, di3, dj1, dj2;
55  bool di1_c, di2_c, di3_c, dj1_c, dj2_c;
56 
57  base_tensor ddi1, ddi2, ddi3, ddj1, ddj2;
58  bool ddi1_c, ddi2_c, ddi3_c, ddj1_c, ddj2_c;
59 
60 
61  /* First invariant tr(E) */
62  void compute_i1() {
63  i1_ = gmm::mat_trace(E);
64  i1_c = true;
65  }
66 
67  void compute_di1() {
68  gmm::resize(di1, N, N);
69  gmm::copy(gmm::identity_matrix(), di1);
70  di1_c = true;
71  }
72 
73  void compute_ddi1() { // not very useful, null tensor
74  ddi1 = base_tensor(N, N, N, N);
75  ddi1_c = true;
76  }
77 
78  inline scalar_type i1()
79  { if (!i1_c) compute_i1(); return i1_; }
80 
81  inline const base_matrix &grad_i1()
82  { if (!di1_c) compute_di1(); return di1; }
83 
84  inline const base_tensor &sym_grad_grad_i1()
85  { if (!ddi1_c) compute_ddi1(); return ddi1; }
86 
87 
88  /* Second invariant (tr(E)^2 - tr(E^2))/2 */
89  void compute_i2() {
90  i2_ = (gmm::sqr(gmm::mat_trace(E))
91  - frobenius_product_trans(E, E)) / scalar_type(2);
92  i2_c = true;
93  }
94 
95  void compute_di2() {
96  gmm::resize(di2, N, N);
97  gmm::copy(gmm::identity_matrix(), di2);
98  gmm::scale(di2, i1());
99  // gmm::add(gmm::scale(gmm::transposed(E), -scalar_type(1)), di2);
100  gmm::add(gmm::scaled(E, -scalar_type(1)), di2);
101  di2_c = true;
102  }
103 
104  void compute_ddi2() {
105  ddi2 = base_tensor(N, N, N, N);
106  for (size_type i = 0; i < N; ++i)
107  for (size_type k = 0; k < N; ++k)
108  ddi2(i,i,k,k) += scalar_type(1);
109  for (size_type i = 0; i < N; ++i)
110  for (size_type j = 0; j < N; ++j) {
111  ddi2(i,j,j,i) -= scalar_type(1)/scalar_type(2);
112  ddi2(j,i,j,i) -= scalar_type(1)/scalar_type(2);
113  }
114  ddi2_c = true;
115  }
116 
117  inline scalar_type i2()
118  { if (!i2_c) compute_i2(); return i2_; }
119 
120  inline const base_matrix &grad_i2()
121  { if (!di2_c) compute_di2(); return di2; }
122 
123  inline const base_tensor &sym_grad_grad_i2()
124  { if (!ddi2_c) compute_ddi2(); return ddi2; }
125 
126  /* Third invariant det(E) */
127  void compute_i3() {
128  Einv = E;
129  i3_ = bgeot::lu_inverse(&(*(Einv.begin())), gmm::mat_nrows(Einv));
130  i3_c = true;
131  }
132 
133  void compute_di3() {
134  scalar_type det = i3();
135  // gmm::resize(di3, N, N);
136  // gmm::copy(gmm::transposed(E), di3);
137  di3 = Einv;
138  // gmm::lu_inverse(di3);
139  gmm::scale(di3, det);
140  di3_c = true;
141  }
142 
143  void compute_ddi3() {
144  ddi3 = base_tensor(N, N, N, N);
145  scalar_type det = i3() / scalar_type(2); // computes also E inverse.
146  for (size_type i = 0; i < N; ++i)
147  for (size_type j = 0; j < N; ++j)
148  for (size_type k = 0; k < N; ++k)
149  for (size_type l = 0; l < N; ++l)
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));
152  ddi3_c = true;
153  }
154 
155  inline scalar_type i3()
156  { if (!i3_c) compute_i3(); return i3_; }
157 
158  inline const base_matrix &grad_i3()
159  { if (!di3_c) compute_di3(); return di3; }
160 
161  inline const base_tensor &sym_grad_grad_i3()
162  { if (!ddi3_c) compute_ddi3(); return ddi3; }
163 
164  /* Invariant j1(E) = i1(E)*i3(E)^(-1/3) */
165  void compute_j1() {
166  j1_ = i1() * ::pow(gmm::abs(i3()), -scalar_type(1) / scalar_type(3));
167  j1_c = true;
168  }
169 
170  void compute_dj1() {
171  dj1 = grad_i1();
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)));
174  dj1_c = true;
175  }
176 
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);
184 
185  for (size_type i = 0; i < N; ++i)
186  for (size_type j = 0; j < N; ++j)
187  for (size_type k = 0; k < N; ++k)
188  for (size_type l = 0; l < N; ++l)
189  ddj1(i,j,k,l) +=
190  (di3_(i, j) * di3_(k, l)) * coeff2
191  - (di1_(i, j) * di3_(k, l) + di1_(k, l) * di3_(i, j)) * coeff1;
192 
193  gmm::scale(ddj1.as_vector(),
194  ::pow(gmm::abs(i3()), -scalar_type(1)/scalar_type(3)));
195  ddj1_c = true;
196  }
197 
198  inline scalar_type j1()
199  { if (!j1_c) compute_j1(); return j1_; }
200 
201  inline const base_matrix &grad_j1()
202  { if (!dj1_c) compute_dj1(); return dj1; }
203 
204  inline const base_tensor &sym_grad_grad_j1()
205  { if (!ddj1_c) compute_ddj1(); return ddj1; }
206 
207  /* Invariant j2(E) = i2(E)*i3(E)^(-2/3) */
208  void compute_j2() {
209  j2_ = i2() * ::pow(gmm::abs(i3()), -scalar_type(2) / scalar_type(3));
210  j2_c = true;
211  }
212 
213  void compute_dj2() {
214  dj2 = grad_i2();
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)));
217  dj2_c = true;
218  }
219 
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),
227  ddj2.as_vector());
228 
229  for (size_type i = 0; i < N; ++i)
230  for (size_type j = 0; j < N; ++j)
231  for (size_type k = 0; k < N; ++k)
232  for (size_type l = 0; l < N; ++l)
233  ddj2(i,j,k,l) +=
234  (di3_(i, j) * di3_(k, l)) * coeff2
235  - (di2_(i, j) * di3_(k, l) + di2_(k, l) * di3_(i, j)) * coeff1;
236 
237  gmm::scale(ddj2.as_vector(),
238  ::pow(gmm::abs(i3()), -scalar_type(2)/scalar_type(3)));
239  ddj2_c = true;
240  }
241 
242 
243  inline scalar_type j2()
244  { if (!j2_c) compute_j2(); return j2_; }
245 
246  inline const base_matrix &grad_j2()
247  { if (!dj2_c) compute_dj2(); return dj2; }
248 
249  inline const base_tensor &sym_grad_grad_j2()
250  { if (!ddj2_c) compute_ddj2(); return ddj2; }
251 
252 
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)
258  {
259  N = gmm::mat_nrows(E);
260  i1_=i2_=i3_=j1_=j2_=0.;
261  }
262 
263  };
264 
265 
266 
267 
268 
269  /* Symmetry check */
270 
271  int check_symmetry(const base_tensor &t) {
272  int flags = 7; size_type N = 3;
273  for (size_type n = 0; n < N; ++n)
274  for (size_type m = 0; m < N; ++m)
275  for (size_type l = 0; l < N; ++l)
276  for (size_type k = 0; k < N; ++k) {
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);
280  }
281  return flags;
282  }
283 
284  /* Member functions of hyperelastic laws */
285 
286  void abstract_hyperelastic_law::random_E(base_matrix &E) {
287  size_type N = gmm::mat_nrows(E);
288  base_matrix Phi(N,N);
289  scalar_type d;
290  do {
291  gmm::fill_random(Phi);
292  d = bgeot::lu_det(&(*(Phi.begin())), N);
293  } while (d < scalar_type(0.01));
294  gmm::mult(gmm::transposed(Phi),Phi,E);
295  gmm::scale(E,-1.); gmm::add(gmm::identity_matrix(),E);
296  gmm::scale(E,-0.5);
297  }
298 
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);
302  bool ok = true;
303 
304  for (size_type count = 0; count < 100; ++count) {
305  random_E(E); random_E(DE);
306  gmm::scale(DE, h);
307  gmm::add(E, DE, E2);
308 
309  base_matrix sigma1(N,N), sigma2(N,N);
310  getfem::base_tensor tdsigma(N,N,N,N);
311  base_matrix dsigma(N,N);
312  gmm::copy(E, E2); gmm::add(DE, E2);
313  sigma(E, sigma1, param, scalar_type(1));
314  sigma(E2, sigma2, param, scalar_type(1));
315 
316  scalar_type d = strain_energy(E2, param, scalar_type(1))
317  - strain_energy(E, param, scalar_type(1));
318  scalar_type d2 = 0;
319  for (size_type i=0; i < N; ++i)
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;
324  ok = false;
325  }
326 
327  grad_sigma(E,tdsigma,param, scalar_type(1));
328  for (size_type i=0; i < N; ++i) {
329  for (size_type j=0; j < N; ++j) {
330  dsigma(i,j) = 0;
331  for (size_type k=0; k < N; ++k) {
332  for (size_type m=0; m < N; ++m) {
333  dsigma(i,j) += tdsigma(i,j,k,m)*DE(k,m);
334  }
335  }
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;
342  ok = false;
343  }
344  }
345  }
346  }
347  GMM_ASSERT1(ok, "Derivative test has failed");
348  }
349 
351  (const base_matrix& F, const base_matrix &E,
352  base_matrix &cauchy_stress, const base_vector &params,
353  scalar_type det_trans) const
354  {
355  size_type N = E.ncols();
356  base_matrix PK2(N,N);
357  sigma(E,PK2,params,det_trans);//second Piola-Kirchhoff stress
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)); //cauchy = 1/J*F*PK2*F^T
362  }
363 
364 
366  (const base_matrix& F, const base_matrix& E,
367  const base_vector &params, scalar_type det_trans,
368  base_tensor &grad_sigma_ul) const
369  {
370  size_type N = E.ncols();
371  base_tensor Cse(N,N,N,N);
372  grad_sigma(E,Cse,params,det_trans);
373  scalar_type mult = 1.0/det_trans;
374  // this is a general transformation for an anisotropic material, very non-efficient;
375  // more effiecient calculations can be overloaded for every specific material
376  for(size_type i = 0; i < N; ++i)
377  for(size_type j = 0; j < N; ++j)
378  for(size_type k = 0; k < N; ++k)
379  for(size_type l = 0; l < N; ++l)
380  {
381  grad_sigma_ul(i,j,k,l) = 0.0;
382  for(size_type m = 0; m < N; ++m)
383  { for(size_type n = 0; n < N; ++n)
384  for(size_type p = 0; p < N; ++p)
385  for(size_type q = 0; q < N; ++q)
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);
388  }
389  grad_sigma_ul(i,j,k,l) *= mult;
390  }
391  }
392 
393  scalar_type SaintVenant_Kirchhoff_hyperelastic_law::strain_energy
394  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
395  // should be optimized, maybe deriving sigma from strain energy
396  if (det_trans <= scalar_type(0))
397  { return 1e200; }
398 
399  return gmm::sqr(gmm::mat_trace(E)) * params[0] / scalar_type(2)
400  + gmm::mat_euclidean_norm_sqr(E) * params[1];
401  }
402 
403  void SaintVenant_Kirchhoff_hyperelastic_law::sigma
404  (const base_matrix &E, base_matrix &result,const base_vector &params, 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);
410  }
411  }
412  void SaintVenant_Kirchhoff_hyperelastic_law::grad_sigma
413  (const base_matrix &E, base_tensor &result,const base_vector &params, scalar_type) const {
414  std::fill(result.begin(), result.end(), scalar_type(0));
415  size_type N = gmm::mat_nrows(E);
416  for (size_type i = 0; i < N; ++i)
417  for (size_type l = 0; l < N; ++l) {
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);
423  }
424  }
425 
427  const base_matrix& E,
428  const base_vector &params,
429  scalar_type det_trans,
430  base_tensor &grad_sigma_ul)const
431  {
432  size_type N = E.ncols();
433  base_tensor Cse(N,N,N,N);
434  grad_sigma(E,Cse,params,det_trans);
435  base_matrix Cinv(N,N); // left Cauchy-Green deform. tens.
436  gmm::mult(F,gmm::transposed(F),Cinv);
437  scalar_type mult=1.0/det_trans;
438  for(size_type i = 0; i < N; ++i)
439  for(size_type j = 0; j < N; ++j)
440  for(size_type k = 0; k < N; ++k)
441  for(size_type l = 0; l < N; ++l)
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;
444  }
445 
446  SaintVenant_Kirchhoff_hyperelastic_law::SaintVenant_Kirchhoff_hyperelastic_law() {
447  nb_params_ = 2;
448  }
449 
450  scalar_type membrane_elastic_law::strain_energy
451  (const base_matrix & /* E */, const base_vector & /* params */, scalar_type) const {
452  // to be done if needed
453  GMM_ASSERT1(false, "To be done");
454  return 0;
455  }
456 
457  void membrane_elastic_law::sigma
458  (const base_matrix &E, base_matrix &result,const base_vector &params, scalar_type det_trans) const {
459  // should be optimized, maybe deriving sigma from strain energy
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);
463  for (size_type i = 0; i < N; ++i)
464  for (size_type j = 0; j < N; ++j) {
465  result(i,j)=0.0;
466  for (size_type k = 0; k < N; ++k)
467  for (size_type l = 0; l < N; ++l)
468  result(i,j)+=tt(i,j,k,l)*E(k,l);
469  }
470  // add pretension in X' direction
471  if(params[4]!=0) result(0,0)+=params[4];
472  // add pretension in Y' direction
473  if(params[5]!=0) result(1,1)+=params[5];
474  // cout<<"sigma="<<result<<endl;
475  }
476 
477  void membrane_elastic_law::grad_sigma
478  (const base_matrix & /* E */, base_tensor &result,
479  const base_vector &params, scalar_type) const {
480  // to be optimized!!
481  std::fill(result.begin(), result.end(), scalar_type(0));
482  scalar_type poisonXY=params[0]*params[1]/params[2]; //Ex*vYX=Ey*vXY
483  //if no G entered, compute G=E/(2*(1+v))
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);
487  // result(0,0,0,1) = 0;
488  // result(0,0,1,0) = 0;
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);
491  // result(1,1,0,1) = 0;
492  // result(1,1,1,0) = 0;
493  result(1,1,1,1) = params[2]/(1-params[1]*poisonXY);
494  // result(0,1,0,0) = 0;
495  result(0,1,0,1) = G;
496  result(0,1,1,0) = G;
497  // result(0,1,1,1) = 0;
498  // result(1,0,0,0) = 0;
499  result(1,0,0,1) = G;
500  result(1,0,1,0) = G;
501  // result(1,0,1,1) = 0;
502  }
503 
504  scalar_type Mooney_Rivlin_hyperelastic_law::strain_energy
505  (const base_matrix &E, const base_vector &params
506  ,scalar_type det_trans) const {
507  //shouldn't negative det_trans be handled here???
508  if (compressible && det_trans <= scalar_type(0)) return 1e200;
509  size_type N = gmm::mat_nrows(E);
510  GMM_ASSERT1(N == 3, "Mooney Rivlin hyperelastic law only defined "
511  "on dimension 3, sorry");
512  base_matrix C = E;
513  gmm::scale(C, scalar_type(2));
514  gmm::add(gmm::identity_matrix(), C);
515  compute_invariants ci(C);
516  size_type i=0;
517  scalar_type C1 = params[i++]; // C10
518  scalar_type W = C1 * (ci.j1() - scalar_type(3));
519  if (!neohookean) {
520  scalar_type C2 = params[i++]; // C01
521  W += C2 * (ci.j2() - scalar_type(3));
522  }
523  if (compressible) {
524  scalar_type D1 = params[i++];
525  W += D1 * gmm::sqr(sqrt(gmm::abs(ci.i3())) - scalar_type(1));
526  }
527  return W;
528  }
529 
530  void Mooney_Rivlin_hyperelastic_law::sigma
531  (const base_matrix &E, base_matrix &result,
532  const base_vector &params, scalar_type det_trans) const {
533  size_type N = gmm::mat_nrows(E);
534  GMM_ASSERT1(N == 3, "Mooney Rivlin hyperelastic law only defined "
535  "on dimension 3, sorry");
536  base_matrix C = E;
537  gmm::scale(C, scalar_type(2));
538  gmm::add(gmm::identity_matrix(), C);
539  compute_invariants ci(C);
540 
541  size_type i=0;
542  scalar_type C1 = params[i++]; // C10
543  gmm::copy(gmm::scaled(ci.grad_j1(), scalar_type(2) * C1), result);
544  if (!neohookean) {
545  scalar_type C2 = params[i++]; // C01
546  gmm::add(gmm::scaled(ci.grad_j2(), scalar_type(2) * C2), result);
547  }
548  if (compressible) {
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);
552 
553 // shouldn't negative det_trans be handled here???
554  if (det_trans <= scalar_type(0))
555  gmm::add(gmm::scaled(C, 1e200), result);
556  }
557  }
558 
559  void Mooney_Rivlin_hyperelastic_law::grad_sigma
560  (const base_matrix &E, base_tensor &result,
561  const base_vector &params, scalar_type) const {
562  size_type N = gmm::mat_nrows(E);
563  GMM_ASSERT1(N == 3, "Mooney Rivlin hyperelastic law only defined "
564  "on dimension 3, sorry");
565  base_matrix C = E;
566  gmm::scale(C, scalar_type(2));
567  gmm::add(gmm::identity_matrix(), C);
568  compute_invariants ci(C);
569 
570  size_type i=0;
571  scalar_type C1 = params[i++]; // C10
572  gmm::copy(gmm::scaled(ci.sym_grad_grad_j1().as_vector(),
573  scalar_type(4)*C1), result.as_vector());
574  if (!neohookean) {
575  scalar_type C2 = params[i++]; // C01
576  gmm::add(gmm::scaled(ci.sym_grad_grad_j2().as_vector(),
577  scalar_type(4)*C2), result.as_vector());
578  }
579  if (compressible) {
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());
584 
585  // second derivatives of W with respect to the third invariant
586  scalar_type A22 = D1 / (scalar_type(2) * pow(gmm::abs(ci.i3()), 1.5));
587  const base_matrix &di = ci.grad_i3();
588  for (size_type l1 = 0; l1 < N; ++l1)
589  for (size_type l2 = 0; l2 < N; ++l2)
590  for (size_type l3 = 0; l3 < N; ++l3)
591  for (size_type l4 = 0; l4 < N; ++l4)
592  result(l1, l2, l3, l4) +=
593  scalar_type(4) * A22 * di(l1, l2) * di(l3, l4);
594  }
595 
596 // GMM_ASSERT1(check_symmetry(result) == 7,
597 // "Fourth order tensor not symmetric : " << result);
598  }
599 
600  Mooney_Rivlin_hyperelastic_law::Mooney_Rivlin_hyperelastic_law
601  (bool compressible_, bool neohookean_)
602  : compressible(compressible_), neohookean(neohookean_)
603  {
604  nb_params_ = 2;
605  if (compressible) ++nb_params_; // D1 != 0
606  if (neohookean) --nb_params_; // C2 == 0
607 
608  }
609 
610 
611 
612 
613  scalar_type Neo_Hookean_hyperelastic_law::strain_energy
614  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
615  if (det_trans <= scalar_type(0)) return 1e200;
616  size_type N = gmm::mat_nrows(E);
617  GMM_ASSERT1(N == 3, "Neo Hookean hyperelastic law only defined "
618  "on dimension 3, sorry");
619  base_matrix C = E;
620  gmm::scale(C, scalar_type(2));
621  gmm::add(gmm::identity_matrix(), C);
622  compute_invariants ci(C);
623 
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);
628  if (bonet)
629  W += lambda/8 * gmm::sqr(logi3);
630  else // Wriggers
631  W += lambda/4 * (ci.i3() - scalar_type(1) - logi3);
632 
633  return W;
634  }
635 
636  void Neo_Hookean_hyperelastic_law::sigma
637  (const base_matrix &E, base_matrix &result,
638  const base_vector &params , scalar_type det_trans) const {
639  size_type N = gmm::mat_nrows(E);
640  GMM_ASSERT1(N == 3, "Neo Hookean hyperelastic law only defined "
641  "on dimension 3, sorry");
642  base_matrix C = E;
643  gmm::scale(C, scalar_type(2));
644  gmm::add(gmm::identity_matrix(), C);
645  compute_invariants ci(C);
646 
647  scalar_type lambda = params[0];
648  scalar_type mu = params[1];
649  gmm::copy(gmm::scaled(ci.grad_i1(), mu), result);
650  if (bonet)
651  gmm::add(gmm::scaled(ci.grad_i3(),
652  (lambda/2 * log(ci.i3()) - mu) / ci.i3()), result);
653  else
654  gmm::add(gmm::scaled(ci.grad_i3(),
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);
658  }
659 
660  void Neo_Hookean_hyperelastic_law::grad_sigma
661  (const base_matrix &E, base_tensor &result,
662  const base_vector &params, scalar_type) const {
663  size_type N = gmm::mat_nrows(E);
664  GMM_ASSERT1(N == 3, "Neo Hookean hyperelastic law only defined "
665  "on dimension 3, sorry");
666  base_matrix C = E;
667  gmm::scale(C, scalar_type(2));
668  gmm::add(gmm::identity_matrix(), C);
669  compute_invariants ci(C);
670 
671  scalar_type lambda = params[0];
672  scalar_type mu = params[1];
673 
674  scalar_type coeff;
675  if (bonet) {
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()),
679  result.as_vector());
680  coeff = (lambda + 2 * mu - lambda * logi3) / gmm::sqr(ci.i3());
681  } else {
682  gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
683  lambda - (lambda + 2 * mu) / ci.i3()),
684  result.as_vector());
685  coeff = (lambda + 2 * mu) / gmm::sqr(ci.i3());
686  }
687 
688  const base_matrix &di = ci.grad_i3();
689  for (size_type l1 = 0; l1 < N; ++l1)
690  for (size_type l2 = 0; l2 < N; ++l2)
691  for (size_type l3 = 0; l3 < N; ++l3)
692  for (size_type l4 = 0; l4 < N; ++l4)
693  result(l1, l2, l3, l4) += coeff * di(l1, l2) * di(l3, l4);
694 
695 // GMM_ASSERT1(check_symmetry(result) == 7,
696 // "Fourth order tensor not symmetric : " << result);
697  }
698 
699  Neo_Hookean_hyperelastic_law::Neo_Hookean_hyperelastic_law(bool bonet_)
700  : bonet(bonet_)
701  {
702  nb_params_ = 2;
703  }
704 
705 
706 
707  scalar_type generalized_Blatz_Ko_hyperelastic_law::strain_energy
708  (const base_matrix &E, const base_vector &params, 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];
712  size_type N = gmm::mat_nrows(E);
713  GMM_ASSERT1(N == 3, "Generalized Blatz Ko hyperelastic law only defined "
714  "on dimension 3, sorry");
715  base_matrix C = E;
716  gmm::scale(C, scalar_type(2));
717  gmm::add(gmm::identity_matrix(), C);
718  compute_invariants ci(C);
719 
720  return pow(a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
721  + c*ci.i2() / ci.i3() + d, n);
722  }
723 
724  void generalized_Blatz_Ko_hyperelastic_law::sigma
725  (const base_matrix &E, base_matrix &result,
726  const base_vector &params, 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];
729  size_type N = gmm::mat_nrows(E);
730  GMM_ASSERT1(N == 3, "Generalized Blatz Ko hyperelastic law only defined "
731  "on dimension 3, sorry");
732  base_matrix C = E;
733  gmm::scale(C, scalar_type(2));
734  gmm::add(gmm::identity_matrix(), C);
735  compute_invariants ci(C);
736 
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()));
744 
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);
750 
751  }
752 
753  void generalized_Blatz_Ko_hyperelastic_law::grad_sigma
754  (const base_matrix &E, base_tensor &result,
755  const base_vector &params, scalar_type) const {
756  scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
757  scalar_type n = params[4];
758  size_type N = gmm::mat_nrows(E);
759  GMM_ASSERT1(N == 3, "Generalized Blatz Ko hyperelastic law only defined "
760  "on dimension 3, sorry");
761  base_matrix C = E;
762  gmm::scale(C, scalar_type(2));
763  gmm::add(gmm::identity_matrix(), C);
764  compute_invariants ci(C);
765 
766 
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;
774 
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());
781 
782  scalar_type nnz = n * (n-1.) * pow(z, n-2.);
783  base_matrix A(3, 3); // second derivatives of W with respect to invariants
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)));
790 
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());
796 
797  for (size_type j = 0; j < N; ++j)
798  for (size_type k = 0; k < N; ++k) {
799  for (size_type l1 = 0; l1 < N; ++l1)
800  for (size_type l2 = 0; l2 < N; ++l2)
801  for (size_type l3 = 0; l3 < N; ++l3)
802  for (size_type l4 = 0; l4 < N; ++l4)
803  result(l1, l2, l3, l4)
804  += 4. * A(j, k) * (*di[j])(l1, l2) * (*di[k])(l3, l4);
805  }
806 
807 // GMM_ASSERT1(check_symmetry(result) == 7,
808 // "Fourth order tensor not symmetric : " << result);
809  }
810 
811  generalized_Blatz_Ko_hyperelastic_law::generalized_Blatz_Ko_hyperelastic_law() {
812  nb_params_ = 5;
813  base_vector V(5);
814  V[0] = 1.0; V[1] = 1.0, V[2] = 1.5; V[3] = -0.5; V[4] = 1.5;
815  }
816 
817 
818  scalar_type Ciarlet_Geymonat_hyperelastic_law::strain_energy
819  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
820  if (det_trans <= scalar_type(0)) return 1e200;
821  size_type N = gmm::mat_nrows(E);
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)
825  + params[2];
826  scalar_type d = params[0]/scalar_type(2) + params[1];
827  scalar_type e = -(scalar_type(3)*(a+b) + c);
828  base_matrix C(N, N);
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);
832  return a * gmm::mat_trace(C)
833  + b * (gmm::sqr(gmm::mat_trace(C)) -
834  gmm::mat_euclidean_norm_sqr(C))/scalar_type(2)
835  + c * det - d * log(det) / scalar_type(2) + e;
836  }
837 
838  void Ciarlet_Geymonat_hyperelastic_law::sigma
839  (const base_matrix &E, base_matrix &result, const base_vector &params, scalar_type det_trans) const {
840  size_type N = gmm::mat_nrows(E);
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)
844  + params[2];
845  scalar_type d = params[0]/scalar_type(2) + params[1];
846  base_matrix C(N, N);
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 "
850  "hyperelastic law");
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);
858  else {
859  scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
860  gmm::add(gmm::scaled(C, scalar_type(2) * c * det - d), result);
861  }
862  }
863 
864  void Ciarlet_Geymonat_hyperelastic_law::grad_sigma
865  (const base_matrix &E, base_tensor &result,const base_vector &params, scalar_type) const {
866  size_type N = gmm::mat_nrows(E);
867  // scalar_type a = params[2];
868  scalar_type b2 = params[1] - params[2]*scalar_type(2); // b*2
869  scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
870  + params[2];
871  scalar_type d = params[0]/scalar_type(2) + params[1];
872  base_matrix C(N, N);
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));
877  for (size_type i = 0; i < N; ++i)
878  for (size_type j = 0; j < N; ++j) {
879  result(i, i, j, j) += 2*b2;
880  result(i, j, i, j) -= b2;
881  result(i, j, j, i) -= b2;
882  for (size_type k = 0; k < N; ++k)
883  for (size_type l = 0; l < N; ++l)
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);
887  }
888 
889 // GMM_ASSERT1(check_symmetry(result) == 7,
890 // "Fourth order tensor not symmetric : " << result);
891  }
892 
893 
894  int levi_civita(int i, int j, int k) {
895  int ii=i+1;
896  int jj=j+1;
897  int kk=k+1; //i,j,k from 0 to 2 !
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)));
903  }
904 
905 
906 
907  scalar_type plane_strain_hyperelastic_law::strain_energy
908  (const base_matrix &E, const base_vector &params, 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);
913  }
914 
915  void plane_strain_hyperelastic_law::sigma
916  (const base_matrix &E, base_matrix &result,const base_vector &params, 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);
923  }
924 
925  void plane_strain_hyperelastic_law::grad_sigma
926  (const base_matrix &E, base_tensor &result,const base_vector &params, 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);
940  }
941 
942 
943 
944 
945 
946 
947 
948  //=========================================================================
949  //
950  // Nonlinear elasticity (old) Brick
951  //
952  //=========================================================================
953 
954  struct nonlinear_elasticity_brick : public virtual_brick {
955 
956  phyperelastic_law AHL;
957 
958  virtual void asm_real_tangent_terms(const model &md, size_type /* ib */,
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 &,
965  size_type region,
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 "
975  "elasticity brick");
976 
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]));
979 
980  const mesh_fem *mf_params = md.pmesh_fem_of_variable(dl[0]);
981  const model_real_plain_vector &params = md.real_variable(dl[0]);
982  const mesh_im &mim = *mims[0];
983 
984  size_type sl = gmm::vect_size(params);
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");
988 
989  mesh_region rg(region);
990  mf_u.linked_mesh().intersect_with_mpi_region(rg);
991 
992  if (version & model::BUILD_MATRIX) {
993  gmm::clear(matl[0]);
994  GMM_TRACE2("Nonlinear elasticity stiffness matrix assembly");
996  (matl[0], mim, mf_u, u, mf_params, params, *AHL, rg);
997  }
998 
999 
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));
1004  }
1005 
1006  }
1007 
1008  nonlinear_elasticity_brick(const phyperelastic_law &AHL_)
1009  : AHL(AHL_) {
1010  set_flags("Nonlinear elasticity brick", false /* is linear*/,
1011  true /* is symmetric */, true /* is coercive */,
1012  true /* is real */, false /* is complex */);
1013  }
1014 
1015  };
1016 
1017  //=========================================================================
1018  // Add a nonlinear elasticity brick.
1019  //=========================================================================
1020 
1021  // Deprecated brick
1023  (model &md, const mesh_im &mim, const std::string &varname,
1024  const phyperelastic_law &AHL, const std::string &dataname,
1025  size_type region) {
1026  pbrick pbr = std::make_shared<nonlinear_elasticity_brick>(AHL);
1027 
1028  model::termlist tl;
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);
1033  }
1034 
1035  //=========================================================================
1036  // Von Mises or Tresca stress computation.
1037  //=========================================================================
1038 
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,
1045  bool tresca) {
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 &params = md.real_variable(dataname);
1052 
1053  size_type sl = gmm::vect_size(params);
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");
1057 
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);
1063  compute_gradient(mf_u, mf_vm, u, GRAD);
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),
1066  IdNFem(NFem, NFem);
1067  base_vector p(NP);
1068  if (!mf_params) gmm::copy(params, p);
1069  base_vector eig(NFem);
1070  base_vector ez(NFem); // vector normal at deformed surface, (ex X ey)
1071  double normEz(0); //norm of ez
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,
1077  gradphit.begin());
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));
1084  if (mf_params)
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) {
1088  //jyh : compute ez, normal on deformed surface
1089  for (unsigned int l = 0; l <NFem; ++l) {
1090  ez[l]=0;
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);
1094  }
1095  normEz= gmm::vect_norm2(ez);
1096  }
1097  //jyh : end compute ez
1098  }
1099  gmm::mult(gradphi, sigmahathat, aux);
1100  gmm::mult(aux, gmm::transposed(gradphi), sigma);
1101 
1102  /* jyh : complete gradphi for virtual 3rd dim (perpendicular to
1103  deformed surface, same thickness) */
1104  if (NFem == 3 && N == 2) {
1105  gmm::resize(gradphi,NFem,NFem);
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;
1111  //jyh : end complete graphi
1112  }
1113 
1114  gmm::scale(sigma, scalar_type(1) / bgeot::lu_det(&(*(gradphi.begin())), NFem));
1115 
1116  if (!tresca) {
1117  /* von mises: norm(deviator(sigma)) */
1118  gmm::add(gmm::scaled(IdNFem, -gmm::mat_trace(sigma) / NFem), sigma);
1119 
1120  //jyh : von mises stress=sqrt(3/2)* norm(sigma) ?
1121  VM[i] = sqrt(3.0/2)*gmm::mat_euclidean_norm(sigma);
1122  } else {
1123  /* else compute the tresca criterion */
1124  //jyh : to be adapted for membrane if necessary
1125  gmm::symmetric_qr_algorithm(sigma, eig);
1126  std::sort(eig.begin(), eig.end());
1127  VM[i] = eig.back() - eig.front();
1128  }
1129  }
1130  }
1131 
1132 
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 &params = md.real_variable(dataname);
1143 
1144  size_type sl = gmm::vect_size(params);
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");
1148 
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");
1152  size_type qqdim = mf_sigma.get_qdim();
1153  size_type ratio = N*N / qqdim;
1154 
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");
1158 
1159  model_real_plain_vector GRAD(mf_sigma.nb_dof()*ratio*NFem/N);
1160  model_real_plain_vector PARAMS(mf_sigma.nb_dof()*NP);
1161 
1162 
1163  getfem::mesh_trans_inv mti(mf_sigma.linked_mesh());
1164  if (mf_params) {
1165  for (size_type i = 0; i < mf_sigma.nb_dof(); ++i)
1166  mti.add_point(mf_sigma.point_of_basic_dof(i));
1167  interpolation(*mf_params, mti, params, PARAMS);
1168  }
1169 
1170  compute_gradient(mf_u, mf_sigma, u, GRAD);
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),
1173  IdNFem(NFem, NFem);
1174 
1175 
1176  base_vector p(NP);
1177  if (!mf_params) gmm::copy(params, p);
1178  base_vector eig(NFem);
1179  base_vector ez(NFem); // vector normal at deformed surface, (ex X ey)
1180  gmm::copy(gmm::identity_matrix(), Id);
1181  gmm::copy(gmm::identity_matrix(), IdNFem);
1182 
1183  // cout << "GRAD = " << GRAD << endl;
1184  // GMM_ASSERT1(false, "oops");
1185 
1186  for (size_type i = 0; i < mf_sigma.nb_dof()/qqdim; ++i) {
1187 // cout << "GRAD size = " << gmm::vect_size(GRAD) << endl;
1188 // cout << "i = " << i << endl;
1189 // cout << "i*NFem*N = " << i*NFem*N << endl;
1190 
1191  gmm::resize(gradphi,NFem,N);
1192  std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
1193  gradphit.begin());
1194  // cout << "GRAD = " << gradphit << endl;
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));
1201  if (mf_params)
1202  gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*ratio*NP,NP)),p);
1203  // cout << "E = " << E << endl;
1204  AHL->sigma(E, sigmahathat, p, scalar_type(1));
1205  // cout << "ok" << endl;
1206  std::copy(sigmahathat.begin(), sigmahathat.end(), SIGMA.begin()+i*N*N);
1207  }
1208  }
1209 
1210 
1211 
1212  // ----------------------------------------------------------------------
1213  //
1214  // Nonlinear incompressibility brick
1215  //
1216  // ----------------------------------------------------------------------
1217 
1218  struct nonlinear_incompressibility_brick : public virtual_brick {
1219 
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,
1227  size_type region,
1228  build_version version) const {
1229 
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 "
1233  "data");
1234  GMM_ASSERT1(mims.size() == 1, "Nonlinear incompressibility brick need a "
1235  "single mesh_im");
1236  GMM_ASSERT1(vl.size() == 2, "Wrong number of variables for nonlinear "
1237  "incompressibility brick");
1238 
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);
1246 
1247  if (version & model::BUILD_MATRIX) {
1248  gmm::clear(matl[0]);
1249  gmm::clear(matl[1]);
1250  asm_nonlinear_incomp_tangent_matrix(matl[0], matl[1],
1251  mim, mf_u, mf_p, u, p, rg);
1252  }
1253 
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));
1258  }
1259  }
1260 
1261 
1262  nonlinear_incompressibility_brick() {
1263  set_flags("Nonlinear incompressibility brick",
1264  false /* is linear*/,
1265  true /* is symmetric */, false /* is coercive */,
1266  true /* is real */, false /* is complex */);
1267  }
1268 
1269 
1270  };
1271 
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>();
1276  model::termlist tl;
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);
1283  }
1284 
1285 
1286 
1287 
1288 
1289  // ----------------------------------------------------------------------
1290  //
1291  // High-level Generic assembly operators
1292  //
1293  // ----------------------------------------------------------------------
1294 
1295 
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; }
1299 
1300 
1301  // Matrix_i2 Operator (second invariant of square matrix of size >= 2)
1302  // (For 2x2 matrices, it is equivalent to det(M))
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);
1308  return true;
1309  }
1310 
1311  // Value : (Trace(M))^2 - Trace(M^2))/2
1312  void value(const arg_list &args, base_tensor &result) const {
1313  size_type N = args[0]->sizes()[0];
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);
1318  for (size_type i = 0; i < N; ++i)
1319  for (size_type j = 0; j < N; ++j)
1320  tr2 += t[i+ j*N] * t[j + i*N];
1321  result[0] = (tr*tr-tr2)/scalar_type(2);
1322  }
1323 
1324  // Derivative : Trace(M)I - M^T
1325  void derivative(const arg_list &args, size_type,
1326  base_tensor &result) const {
1327  size_type N = args[0]->sizes()[0];
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();
1332  for (size_type j = 0; j < N; ++j)
1333  for (size_type i = 0; i < N; ++i, ++it)
1334  *it = ((i == j) ? tr : scalar_type(0)) - t[i*N+j];
1335  GMM_ASSERT1(it == result.end(), "Internal error");
1336  }
1337 
1338  // Second derivative : I@I - \delta_{il}\delta_{jk}
1339  void second_derivative(const arg_list &args, size_type, size_type,
1340  base_tensor &result) const {
1341  size_type N = args[0]->sizes()[0];
1342  gmm::clear(result.as_vector());
1343  for (size_type i = 0; i < N; ++i)
1344  for (size_type j = 0; j < N; ++j) {
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);
1347  }
1348  }
1349  };
1350 
1351 
1352  // Matrix_j1 Operator
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);
1358  return true;
1359  }
1360 
1361  // Value : Trace(M)/(det(M)^1/3)
1362  void value(const arg_list &args, base_tensor &result) const {
1363  size_type N = args[0]->sizes()[0];
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);
1369  if (det > 0)
1370  result[0] = tr / pow(det, scalar_type(1)/scalar_type(3));
1371  else
1372  result[0] = 1.E200;
1373  }
1374 
1375  // Derivative : (I-Trace(M)*M^{-T}/3)/(det(M)^1/3)
1376  void derivative(const arg_list &args, size_type,
1377  base_tensor &result) const {
1378  size_type N = args[0]->sizes()[0];
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);
1384  if (det > 0) {
1385  base_tensor::iterator it = result.begin();
1386  for (size_type j = 0; j < N; ++j)
1387  for (size_type i = 0; i < N; ++i, ++it)
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");
1392  } else
1393  std::fill(result.begin(), result.end(), 1.E200);
1394  }
1395 
1396  // Second derivative : (-M^{-T}@I + Trace(M)*M^{-T}_{ik}M^{-T}_{lj}
1397  // -I@M^{-T} + Trace(M)*M^{-T}@M^{-T}/3)/(3det(M)^1/3)
1398  void second_derivative(const arg_list &args, size_type, size_type,
1399  base_tensor &result) const {
1400  size_type N = args[0]->sizes()[0];
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);
1406  if (det > 0) {
1407  base_tensor::iterator it = result.begin();
1408  for (size_type l = 0; l < N; ++l)
1409  for (size_type k = 0; k < N; ++k)
1410  for (size_type j = 0; j < N; ++j)
1411  for (size_type i = 0; i < N; ++i, ++it)
1412  *it = (- ((k == l) ? M(j, i) : scalar_type(0))
1413  + tr*M(i,k)*M(l,j)
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");
1418  } else
1419  std::fill(result.begin(), result.end(), 1.E200);
1420  }
1421  };
1422 
1423 
1424  // Matrix_j2 Operator
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);
1430  return true;
1431  }
1432 
1433  // Value : i2(M)/(det(M)^2/3)
1434  void value(const arg_list &args, base_tensor &result) const {
1435  size_type N = args[0]->sizes()[0];
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);
1441  for (size_type i = 0; i < N; ++i)
1442  for (size_type j = 0; j < N; ++j)
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);
1446  if (det > 0)
1447  result[0] = i2 / pow(det, scalar_type(2)/scalar_type(3));
1448  else
1449  result[0] = 1.E200;
1450  }
1451 
1452  // Derivative : (Trace(M)*I-M^T-2i2(M)M^{-T}/3)/(det(M)^2/3)
1453  void derivative(const arg_list &args, size_type,
1454  base_tensor &result) const {
1455  size_type N = args[0]->sizes()[0];
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);
1462  for (size_type i = 0; i < N; ++i)
1463  for (size_type j = 0; j < N; ++j)
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();
1468  for (size_type j = 0; j < N; ++j)
1469  for (size_type i = 0; i < N; ++i, ++it)
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");
1474  }
1475 
1476  // Second derivative
1477  void second_derivative(const arg_list &args, size_type, size_type,
1478  base_tensor &result) const {
1479  size_type N = args[0]->sizes()[0];
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);
1486  for (size_type i = 0; i < N; ++i)
1487  for (size_type j = 0; j < N; ++j)
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();
1492  for (size_type l = 0; l < N; ++l)
1493  for (size_type k = 0; k < N; ++k)
1494  for (size_type j = 0; j < N; ++j)
1495  for (size_type i = 0; i < N; ++i, ++it)
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));
1503  }
1504  };
1505 
1506  // Right-Cauchy-Green operator (F^{T}F)
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]);
1511  return true;
1512  }
1513 
1514  // Value : F^{T}F
1515  void value(const arg_list &args, base_tensor &result) const {
1516  // to be verified
1517  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1518  base_tensor::iterator it = result.begin();
1519  for (size_type j = 0; j < n; ++j)
1520  for (size_type i = 0; i < n; ++i, ++it) {
1521  *it = scalar_type(0);
1522  for (size_type k = 0; k < m; ++k)
1523  *it += (*(args[0]))[i*m+k] * (*(args[0]))[j*m+k];
1524  }
1525  }
1526 
1527  // Derivative : F{kj}delta{li}+F{ki}delta{lj}
1528  // (comes from H -> H^{T}F + F^{T}H)
1529  void derivative(const arg_list &args, size_type,
1530  base_tensor &result) const { // to be verified
1531  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1532  base_tensor::iterator it = result.begin();
1533  for (size_type l = 0; l < n; ++l)
1534  for (size_type k = 0; k < m; ++k)
1535  for (size_type j = 0; j < n; ++j)
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];
1540  }
1541  GMM_ASSERT1(it == result.end(), "Internal error");
1542  }
1543 
1544  // Second derivative :
1545  // A{ijklop}=delta{ok}delta{li}delta{pj} + delta{ok}delta{pi}delta{lj}
1546  // comes from (H,K) -> H^{T}K + K^{T}H
1547  void second_derivative(const arg_list &args, size_type, size_type,
1548  base_tensor &result) const { // to be verified
1549  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1550  base_tensor::iterator it = result.begin();
1551  for (size_type p = 0; p < n; ++p)
1552  for (size_type o = 0; o < m; ++o)
1553  for (size_type l = 0; l < n; ++l)
1554  for (size_type k = 0; k < m; ++k)
1555  for (size_type j = 0; j < n; ++j)
1556  for (size_type i = 0; i < n; ++i, ++it) {
1557  *it = scalar_type(0);
1558  if (o == k) {
1559  if (l == i && p == j) *it += scalar_type(1);
1560  if (p == i && l == j) *it += scalar_type(1);
1561  }
1562  }
1563  GMM_ASSERT1(it == result.end(), "Internal error");
1564  }
1565  };
1566 
1567  // Left-Cauchy-Green operator (FF^{T})
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]);
1572  return true;
1573  }
1574 
1575  // Value : FF^{T}
1576  void value(const arg_list &args, base_tensor &result) const {
1577  // to be verified
1578  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1579  base_tensor::iterator it = result.begin();
1580  for (size_type j = 0; j < m; ++j)
1581  for (size_type i = 0; i < m; ++i, ++it) {
1582  *it = scalar_type(0);
1583  for (size_type k = 0; k < n; ++k)
1584  *it += (*(args[0]))[k*m+i] * (*(args[0]))[k*m+j];
1585  }
1586  }
1587 
1588  // Derivative : F{jl}delta{ik}+F{il}delta{kj}
1589  // (comes from H -> HF^{T} + FH^{T})
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();
1594  for (size_type l = 0; l < n; ++l)
1595  for (size_type k = 0; k < m; ++k)
1596  for (size_type j = 0; j < m; ++j)
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];
1601  }
1602  GMM_ASSERT1(it == result.end(), "Internal error");
1603  }
1604 
1605  // Second derivative :
1606  // A{ijklop}=delta{ki}delta{lp}delta{oj} + delta{oi}delta{pl}delta{kj}
1607  // comes from (H,K) -> HK^{T} + KH^{T}
1608  void second_derivative(const arg_list &args, size_type, size_type,
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();
1612  for (size_type p = 0; p < n; ++p)
1613  for (size_type o = 0; o < m; ++o)
1614  for (size_type l = 0; l < n; ++l)
1615  for (size_type k = 0; k < m; ++k)
1616  for (size_type j = 0; j < m; ++j)
1617  for (size_type i = 0; i < m; ++i, ++it) {
1618  *it = scalar_type(0);
1619  if (p == l) {
1620  if (k == i && o == j) *it += scalar_type(1);
1621  if (o == i && k == j) *it += scalar_type(1);
1622  }
1623  }
1624  GMM_ASSERT1(it == result.end(), "Internal error");
1625  }
1626  };
1627 
1628 
1629  // Green-Lagrangian operator (F^{T}F-I)/2
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]);
1634  return true;
1635  }
1636 
1637  // Value : F^{T}F
1638  void value(const arg_list &args, base_tensor &result) const {
1639  // to be verified
1640  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1641  base_tensor::iterator it = result.begin();
1642  for (size_type j = 0; j < n; ++j)
1643  for (size_type i = 0; i < n; ++i, ++it) {
1644  *it = scalar_type(0);
1645  for (size_type k = 0; k < m; ++k)
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);
1648  }
1649  }
1650 
1651  // Derivative : (F{kj}delta{li}+F{ki}delta{lj})/2
1652  // (comes from H -> (H^{T}F + F^{T}H)/2)
1653  void derivative(const arg_list &args, size_type,
1654  base_tensor &result) const { // to be verified
1655  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1656  base_tensor::iterator it = result.begin();
1657  for (size_type l = 0; l < n; ++l)
1658  for (size_type k = 0; k < m; ++k)
1659  for (size_type j = 0; j < n; ++j)
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);
1664  }
1665  GMM_ASSERT1(it == result.end(), "Internal error");
1666  }
1667 
1668  // Second derivative :
1669  // A{ijklop}=(delta{ok}delta{li}delta{pj} + delta{ok}delta{pi}delta{lj})/2
1670  // comes from (H,K) -> (H^{T}K + K^{T}H)/2
1671  void second_derivative(const arg_list &args, size_type, size_type,
1672  base_tensor &result) const { // to be verified
1673  size_type m = args[0]->sizes()[0], n = args[0]->sizes()[1];
1674  base_tensor::iterator it = result.begin();
1675  for (size_type p = 0; p < n; ++p)
1676  for (size_type o = 0; o < m; ++o)
1677  for (size_type l = 0; l < n; ++l)
1678  for (size_type k = 0; k < m; ++k)
1679  for (size_type j = 0; j < n; ++j)
1680  for (size_type i = 0; i < n; ++i, ++it) {
1681  *it = scalar_type(0);
1682  if (o == k) {
1683  if (l == i && p == j) *it += scalar_type(0.5);
1684  if (p == i && l == j) *it += scalar_type(0.5);
1685  }
1686  }
1687  GMM_ASSERT1(it == result.end(), "Internal error");
1688  }
1689  };
1690 
1691  // Cauchy stress tensor from second Piola-Kirchhoff stress tensor
1692  // (I+Grad_u)\sigma(I+Grad_u')/det(I+Grad_u)
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]);
1701  return true;
1702  }
1703 
1704  // Value : (I+Grad_u)\sigma(I+Grad_u')/det(I+Grad_u)
1705  void value(const arg_list &args, base_tensor &result) const {
1706  size_type N = args[0]->sizes()[0];
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());
1716  }
1717 
1718  // Derivative :
1719  void derivative(const arg_list &args, size_type nder,
1720  base_tensor &result) const { // to be verified
1721  size_type N = args[0]->sizes()[0];
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);
1726 
1727  base_tensor::iterator it = result.begin();
1728 
1729  switch (nder) {
1730  case 1: // Derivative with respect to sigma : F_ik F_jl / det(F)
1731  for (size_type l = 0; l < N; ++l)
1732  for (size_type k = 0; k < N; ++k)
1733  for (size_type j = 0; j < N; ++j)
1734  for (size_type i = 0; i < N; ++i, ++it)
1735  *it = F(i,k) * F(j,l) / det;
1736  break;
1737 
1738  case 2:
1739  {
1740  // Derivative with respect to Grad_u :
1741  // (delta_ik sigma_lm F_jm + F_im sigma_mk delta_lj
1742  // - F_im sigma_mn F_jn F^{-1}_lk) / det(F)
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);
1748  for (size_type l = 0; l < N; ++l)
1749  for (size_type k = 0; k < N; ++k)
1750  for (size_type j = 0; j < N; ++j)
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;
1756  }
1757  }
1758  break;
1759  }
1760  GMM_ASSERT1(it == result.end(), "Internal error");
1761  }
1762 
1763  // Second derivative : to be implemented
1764  void second_derivative(const arg_list &, size_type, size_type,
1765  base_tensor &) const {
1766  GMM_ASSERT1(false, "Sorry, not implemented");
1767  }
1768  };
1769 
1770 
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]);
1778  return true;
1779  }
1780 
1781  // Value :
1782  void value(const arg_list &args, base_tensor &result) const {
1783  size_type N = args[0]->sizes()[0];
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);
1793 
1794  AHL->sigma(E, sigma, params, det);
1795  gmm::copy(sigma.as_vector(), result.as_vector());
1796  }
1797 
1798  // Derivative :
1799  void derivative(const arg_list &args, size_type nder,
1800  base_tensor &result) const {
1801  size_type N = args[0]->sizes()[0];
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);
1812 
1813  GMM_ASSERT1(nder == 1, "Sorry, the derivative of this hyperelastic "
1814  "law with respect to its parameters is not available.");
1815 
1816  AHL->grad_sigma(E, grad_sigma, params, det);
1817 
1818  base_tensor::iterator it = result.begin();
1819  for (size_type l = 0; l < N; ++l)
1820  for (size_type k = 0; k < N; ++k)
1821  for (size_type j = 0; j < N; ++j)
1822  for (size_type i = 0; i < N; ++i, ++it) {
1823  *it = scalar_type(0);
1824  for (size_type m = 0; m < N; ++m)
1825  *it += grad_sigma(i,j,m,l) * Gu(k, m); // to be optimized
1826  }
1827  GMM_ASSERT1(it == result.end(), "Internal error");
1828  }
1829 
1830 
1831  // Second derivative : not implemented (not necessary)
1832  void second_derivative(const arg_list &, size_type, size_type,
1833  base_tensor &) const {
1834  GMM_ASSERT1(false, "Sorry, second derivative not implemented");
1835  }
1836 
1837  AHL_wrapper_sigma(const phyperelastic_law &A) : AHL(A) {}
1838 
1839  };
1840 
1841 
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);
1849  return true;
1850  }
1851 
1852  // Value :
1853  void value(const arg_list &args, base_tensor &result) const {
1854  size_type N = args[0]->sizes()[0];
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); // not necessary here
1864 
1865  result[0] = AHL->strain_energy(E, params, det);
1866  }
1867 
1868  // Derivative :
1869  void derivative(const arg_list &args, size_type nder,
1870  base_tensor &result) const {
1871  size_type N = args[0]->sizes()[0];
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); // not necessary here
1881 
1882  GMM_ASSERT1(nder == 1, "Sorry, Cannot derive the potential with "
1883  "respect to law parameters.");
1884 
1885  AHL->sigma(E, sigma, params, det);
1886  gmm::mult(Gu, sigma, E);
1887  gmm::copy(E.as_vector(), result.as_vector());
1888  }
1889 
1890 
1891  // Second derivative :
1892  void second_derivative(const arg_list &args, size_type nder1,
1893  size_type nder2, base_tensor &result) const {
1894 
1895  size_type N = args[0]->sizes()[0];
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);
1906 
1907  GMM_ASSERT1(nder1 == 1 && nder2 == 1, "Sorry, Cannot derive the "
1908  "potential with respect to law parameters.");
1909 
1910  AHL->sigma(E, sigma, params, det);
1911  AHL->grad_sigma(E, grad_sigma, params, det);
1912 
1913  base_tensor::iterator it = result.begin();
1914  for (size_type l = 0; l < N; ++l)
1915  for (size_type k = 0; k < N; ++k)
1916  for (size_type j = 0; j < N; ++j)
1917  for (size_type i = 0; i < N; ++i, ++it) {
1918  *it = scalar_type(0);
1919  if (i == k) *it += sigma(l,j);
1920 
1921  for (size_type m = 0; m < N; ++m) // to be optimized
1922  for (size_type n = 0; n < N; ++n)
1923  *it += grad_sigma(n,j,m,l) * Gu(k, m) * Gu(i, n);
1924  }
1925  GMM_ASSERT1(it == result.end(), "Internal error");
1926 
1927  }
1928 
1929  AHL_wrapper_potential(const phyperelastic_law &A) : AHL(A) {}
1930 
1931  };
1932 
1933 
1934  // Saint-Venant Kirchhoff tensor:
1935  // lambda Tr(E)Id + 2*mu*E with E = (Grad_u+Grad_u'+Grad_u'Grad_u)/2
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]);
1942  return true;
1943  }
1944 
1945  // Value : Tr(E) + 2*mu*E
1946  void value(const arg_list &args, base_tensor &result) const {
1947  size_type N = args[0]->sizes()[0];
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);
1955 
1956  base_tensor::iterator it = result.begin();
1957  for (size_type j = 0; j < N; ++j)
1958  for (size_type i = 0; i < N; ++i, ++it) {
1959  *it = 2*mu*E(i,j);
1960  if (i == j) *it += lambda*trE;
1961  }
1962  }
1963 
1964  // Derivative / u: lambda Tr(H + H'*U) Id + mu(H + H' + H'*U + U'*H)
1965  // with U = Grad_u, H = Grad_Test_u
1966  // Implementation: A{ijkl} = lambda(delta{kl}delta{ij} + U{kl}delta{ij})
1967  // + mu(delta{ik}delta{jl} + delta{il}delta{jk} + U{kj}delta{il} +
1968  // U{ki}delta{lj})
1969  void derivative(const arg_list &args, size_type nder,
1970  base_tensor &result) const {
1971  size_type N = args[0]->sizes()[0];
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());
1975  if (nder > 1) {
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));
1979  }
1980  base_tensor::iterator it = result.begin();
1981  switch (nder) {
1982  case 1: // Derivative with respect to u
1983  for (size_type l = 0; l < N; ++l)
1984  for (size_type k = 0; k < N; ++k)
1985  for (size_type j = 0; j < N; ++j)
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);
1994  }
1995  break;
1996  case 2:
1997  // Derivative with respect to lambda: Tr(E)Id
1998  trE = gmm::mat_trace(E);
1999  for (size_type j = 0; j < N; ++j)
2000  for (size_type i = 0; i < N; ++i, ++it) {
2001  *it = scalar_type(0); if (i == j) *it += trE;
2002  }
2003  // Derivative with respect to mu: 2E
2004  for (size_type j = 0; j < N; ++j)
2005  for (size_type i = 0; i < N; ++i, ++it) {
2006  *it += 2*E(i,j);
2007  }
2008  break;
2009  default: GMM_ASSERT1(false, "Internal error");
2010  }
2011  GMM_ASSERT1(it == result.end(), "Internal error");
2012  }
2013 
2014  // Second derivative : not implemented
2015  void second_derivative(const arg_list &, size_type, size_type,
2016  base_tensor &) const {
2017  GMM_ASSERT1(false, "Sorry, second derivative not implemented");
2018  }
2019  };
2020 
2021 
2022 
2023  static bool init_predef_operators() {
2024 
2025  ga_predef_operator_tab &PREDEF_OPERATORS
2027 
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>());
2040 
2041  PREDEF_OPERATORS.add_method
2042  ("Cauchy_stress_from_PK2", std::make_shared<Cauchy_stress_from_PK2>());
2043 
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>()));
2064 
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)));
2089 
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)));
2112 
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)));
2137 
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)));
2162 
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)));
2187 
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)));
2212 
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)));
2237 
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)));
2262 
2263  return true;
2264  }
2265 
2266  // declared in getfem_generic_assembly.cc
2267  bool predef_operators_nonlinear_elasticity_initialized
2268  = init_predef_operators();
2269 
2270 
2271  std::string adapt_law_name(const std::string &lawname, size_type N) {
2272  std::string adapted_lawname = lawname;
2273 
2274  for (size_type i = 0; i < lawname.size(); ++i)
2275  if (adapted_lawname[i] == ' ') adapted_lawname[i] = '_';
2276 
2277  if (adapted_lawname.compare("SaintVenant_Kirchhoff") == 0) {
2278  adapted_lawname = "Saint_Venant_Kirchhoff";
2279  } else if (adapted_lawname.compare("Saint_Venant_Kirchhoff") == 0) {
2280 
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;
2295  } else
2296  GMM_ASSERT1(false, lawname << " is not a known hyperelastic law");
2297 
2298  return adapted_lawname;
2299  }
2300 
2301 
2303  (model &md, const mesh_im &mim, const std::string &lawname,
2304  const std::string &varname, const std::string &params,
2305  size_type region) {
2306  std::string test_varname = "Test_" + sup_previous_and_dot_to_varname(varname);
2307  size_type N = mim.linked_mesh().dim();
2308  GMM_ASSERT1(N >= 2 && N <= 3,
2309  "Finite strain elasticity brick works only in 2D or 3D");
2310 
2311  const mesh_fem *mf = md.pmesh_fem_of_variable(varname);
2312  GMM_ASSERT1(mf, "Finite strain elasticity brick can only be applied on "
2313  "fem variables");
2314  size_type Q = mf->get_qdim();
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");
2317 
2318  std::string adapted_lawname = adapt_law_name(lawname, N);
2319 
2320  std::string expr = "((Id(meshdim)+Grad_"+varname+")*(" + adapted_lawname
2321  + "_PK2(Grad_"+varname+","+params+"))):Grad_" + test_varname;
2322 
2323  return add_nonlinear_term
2324  (md, mim, expr, region, true, false,
2325  "Finite strain elasticity brick for " + adapted_lawname + " law");
2326  }
2327 
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);
2333 
2334  std::string expr
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 + "))" ;
2339  return add_nonlinear_term
2340  (md, mim, expr, region, true, false,
2341  "Finite strain incompressibility brick");
2342  }
2343 
2345  (model &md, const std::string &lawname, const std::string &varname,
2346  const std::string &params, const mesh_fem &mf_vm,
2347  model_real_plain_vector &VM, const mesh_region &rg) {
2348  size_type N = mf_vm.linked_mesh().dim();
2349 
2350  std::string adapted_lawname = adapt_law_name(lawname, N);
2351 
2352  std::string expr = "sqrt(3/2)*Norm(Deviator(Cauchy_stress_from_PK2("
2353  + adapted_lawname + "_PK2(Grad_" + varname + "," + params + "),Grad_"
2354  + varname + ")))";
2355  ga_interpolation_Lagrange_fem(md, expr, mf_vm, VM, rg);
2356  }
2357 
2358 
2359 
2360 } /* end of namespace getfem. */
2361 
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 &params, 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 &params, 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.
Definition: gmm_blas.h:69
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:978
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
Definition: gmm_blas.h:558
void fill_random(L &l)
fill a vector or matrix with random value (uniform [-1,1]).
Definition: gmm_blas.h:130
linalg_traits< M >::value_type mat_trace(const M &m)
Trace of a matrix.
Definition: gmm_blas.h:529
number_traits< typename linalg_traits< M >::value_type >::magnitude_type mat_euclidean_norm(const M &m)
Euclidean norm of a matrix.
Definition: gmm_blas.h:637
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:59
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:210
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*‍/
Definition: gmm_blas.h:1664
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1277
number_traits< typename linalg_traits< M >::value_type >::magnitude_type mat_euclidean_norm_sqr(const M &m)
*‍/
Definition: gmm_blas.h:627
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
Definition: bgeot_poly.h:49
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 .
Definition: bgeot_poly.cc:47
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 &params, 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 &params, 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
Definition: getfem_models.h:49
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector &params, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian