GetFEM++  5.3
getfem_nonlinear_elasticity.cc
1 /*===========================================================================
2 
3  Copyright (C) 2000-2017 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()
225  / scalar_type(2);
226  ddj2 = sym_grad_grad_i2();
227  gmm::add(gmm::scaled(sym_grad_grad_i3().as_vector(), -i2() * coeff1),
228  ddj2.as_vector());
229 
230  for (size_type i = 0; i < N; ++i)
231  for (size_type j = 0; j < N; ++j)
232  for (size_type k = 0; k < N; ++k)
233  for (size_type l = 0; l < N; ++l)
234  ddj2(i,j,k,l) +=
235  (di3_(i, j) * di3_(k, l)) * coeff2
236  - (di2_(i, j) * di3_(k, l) + di2_(k, l) * di3_(i, j)) * coeff1;
237 
238  gmm::scale(ddj2.as_vector(),
239  ::pow(gmm::abs(i3()), -scalar_type(2)/scalar_type(3)));
240  ddj2_c = true;
241  }
242 
243 
244  inline scalar_type j2()
245  { if (!j2_c) compute_j2(); return j2_; }
246 
247  inline const base_matrix &grad_j2()
248  { if (!dj2_c) compute_dj2(); return dj2; }
249 
250  inline const base_tensor &sym_grad_grad_j2()
251  { if (!ddj2_c) compute_ddj2(); return ddj2; }
252 
253 
254  compute_invariants(const base_matrix &EE)
255  : E(EE), i1_c(false), i2_c(false), i3_c(false),
256  j1_c(false), j2_c(false), di1_c(false), di2_c(false), di3_c(false),
257  dj1_c(false), dj2_c(false), ddi1_c(false), ddi2_c(false),
258  ddi3_c(false), ddj1_c(false), ddj2_c(false)
259  { N = gmm::mat_nrows(E); }
260 
261  };
262 
263 
264 
265 
266 
267  /* Symmetry check */
268 
269  int check_symmetry(const base_tensor &t) {
270  int flags = 7; size_type N = 3;
271  for (size_type n = 0; n < N; ++n)
272  for (size_type m = 0; m < N; ++m)
273  for (size_type l = 0; l < N; ++l)
274  for (size_type k = 0; k < N; ++k) {
275  if (gmm::abs(t(n,m,l,k) - t(l,k,n,m))>1e-5) flags &= (~1);
276  if (gmm::abs(t(n,m,l,k) - t(m,n,l,k))>1e-5) flags &= (~2);
277  if (gmm::abs(t(n,m,l,k) - t(n,m,k,l))>1e-5) flags &= (~4);
278  }
279  return flags;
280  }
281 
282  /* Member functions of hyperelastic laws */
283 
284  void abstract_hyperelastic_law::random_E(base_matrix &E) {
285  size_type N = gmm::mat_nrows(E);
286  base_matrix Phi(N,N);
287  scalar_type d;
288  do {
289  gmm::fill_random(Phi);
290  d = bgeot::lu_det(&(*(Phi.begin())), N);
291  } while (d < scalar_type(0.01));
292  gmm::mult(gmm::transposed(Phi),Phi,E);
293  gmm::scale(E,-1.); gmm::add(gmm::identity_matrix(),E);
294  gmm::scale(E,-0.5);
295  }
296 
297  void abstract_hyperelastic_law::test_derivatives
298  (size_type N, scalar_type h, const base_vector& param) const {
299  base_matrix E(N,N), E2(N,N), DE(N,N);
300  bool ok = true;
301 
302  for (size_type count = 0; count < 100; ++count) {
303  random_E(E); random_E(DE);
304  gmm::scale(DE, h);
305  gmm::add(E, DE, E2);
306 
307  base_matrix sigma1(N,N), sigma2(N,N);
308  getfem::base_tensor tdsigma(N,N,N,N);
309  base_matrix dsigma(N,N);
310  gmm::copy(E, E2); gmm::add(DE, E2);
311  sigma(E, sigma1, param, scalar_type(1));
312  sigma(E2, sigma2, param, scalar_type(1));
313 
314  scalar_type d = strain_energy(E2, param, scalar_type(1))
315  - strain_energy(E, param, scalar_type(1));
316  scalar_type d2 = 0;
317  for (size_type i=0; i < N; ++i)
318  for (size_type j=0; j < N; ++j) d2 += sigma1(i,j)*DE(i,j);
319  if (gmm::abs(d-d2)/(gmm::abs(d)+1e-40) > 1e-4) {
320  cout << "Test " << count << " wrong derivative of strain_energy, d="
321  << d/h << ", d2=" << d2/h << endl;
322  ok = false;
323  }
324 
325  grad_sigma(E,tdsigma,param, scalar_type(1));
326  for (size_type i=0; i < N; ++i) {
327  for (size_type j=0; j < N; ++j) {
328  dsigma(i,j) = 0;
329  for (size_type k=0; k < N; ++k) {
330  for (size_type m=0; m < N; ++m) {
331  dsigma(i,j) += tdsigma(i,j,k,m)*DE(k,m);
332  }
333  }
334  sigma2(i,j) -= sigma1(i,j);
335  if (gmm::abs(dsigma(i,j) - sigma2(i,j))
336  /(gmm::abs(dsigma(i,j)) + 1e-40) > 1.5e-4) {
337  cout << "Test " << count << " wrong derivative of sigma, i="
338  << i << ", j=" << j << ", dsigma=" << dsigma(i,j)/h
339  << ", var sigma = " << sigma2(i,j)/h << endl;
340  ok = false;
341  }
342  }
343  }
344  }
345  GMM_ASSERT1(ok, "Derivative test has failed");
346  }
347 
349  (const base_matrix& F, const base_matrix &E,
350  base_matrix &cauchy_stress, const base_vector &params,
351  scalar_type det_trans) const
352  {
353  size_type N = E.ncols();
354  base_matrix PK2(N,N);
355  sigma(E,PK2,params,det_trans);//second Piola-Kirchhoff stress
356  base_matrix aux(N,N);
357  gmm::mult(F,PK2,aux);
358  gmm::mult(aux,gmm::transposed(F),cauchy_stress);
359  gmm::scale(cauchy_stress,scalar_type(1.0/det_trans)); //cauchy = 1/J*F*PK2*F^T
360  }
361 
362 
364  (const base_matrix& F, const base_matrix& E,
365  const base_vector &params, scalar_type det_trans,
366  base_tensor &grad_sigma_ul) const
367  {
368  size_type N = E.ncols();
369  base_tensor Cse(N,N,N,N);
370  grad_sigma(E,Cse,params,det_trans);
371  scalar_type mult = 1.0/det_trans;
372  // this is a general transformation for an anisotropic material, very non-efficient;
373  // more effiecient calculations can be overloaded for every specific material
374  for(size_type i = 0; i < N; ++i)
375  for(size_type j = 0; j < N; ++j)
376  for(size_type k = 0; k < N; ++k)
377  for(size_type l = 0; l < N; ++l)
378  {
379  grad_sigma_ul(i,j,k,l) = 0.0;
380  for(size_type m = 0; m < N; ++m)
381  { for(size_type n = 0; n < N; ++n)
382  for(size_type p = 0; p < N; ++p)
383  for(size_type q = 0; q < N; ++q)
384  grad_sigma_ul(i,j,k,l)+=
385  F(i,m)*F(j,n)*F(k,p)*F(l,q)*Cse(m,n,p,q);
386  }
387  grad_sigma_ul(i,j,k,l) *= mult;
388  }
389  }
390 
391  scalar_type SaintVenant_Kirchhoff_hyperelastic_law::strain_energy
392  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
393  // should be optimized, maybe deriving sigma from strain energy
394  if (det_trans <= scalar_type(0))
395  { return 1e200; }
396 
397  return gmm::sqr(gmm::mat_trace(E)) * params[0] / scalar_type(2)
398  + gmm::mat_euclidean_norm_sqr(E) * params[1];
399  }
400 
401  void SaintVenant_Kirchhoff_hyperelastic_law::sigma
402  (const base_matrix &E, base_matrix &result,const base_vector &params, scalar_type det_trans) const {
403  gmm::copy(gmm::identity_matrix(), result);
404  gmm::scale(result, params[0] * gmm::mat_trace(E));
405  gmm::add(gmm::scaled(E, 2 * params[1]), result);
406  if (det_trans <= scalar_type(0)) {
407  gmm::add(gmm::scaled(E, 1e200), result);
408  }
409  }
410  void SaintVenant_Kirchhoff_hyperelastic_law::grad_sigma
411  (const base_matrix &E, base_tensor &result,const base_vector &params, scalar_type) const {
412  std::fill(result.begin(), result.end(), scalar_type(0));
413  size_type N = gmm::mat_nrows(E);
414  for (size_type i = 0; i < N; ++i)
415  for (size_type l = 0; l < N; ++l) {
416  result(i, i, l, l) += params[0];
417  result(i, l, i, l) += params[1]/scalar_type(2);
418  result(i, l, l, i) += params[1]/scalar_type(2);
419  result(l, i, i, l) += params[1]/scalar_type(2);
420  result(l, i, l, i) += params[1]/scalar_type(2);
421  }
422  }
423 
425  const base_matrix& E,
426  const base_vector &params,
427  scalar_type det_trans,
428  base_tensor &grad_sigma_ul)const
429  {
430  size_type N = E.ncols();
431  base_tensor Cse(N,N,N,N);
432  grad_sigma(E,Cse,params,det_trans);
433  base_matrix Cinv(N,N); // left Cauchy-Green deform. tens.
434  gmm::mult(F,gmm::transposed(F),Cinv);
435  scalar_type mult=1.0/det_trans;
436  for(size_type i = 0; i < N; ++i)
437  for(size_type j = 0; j < N; ++j)
438  for(size_type k = 0; k < N; ++k)
439  for(size_type l = 0; l < N; ++l)
440  grad_sigma_ul(i, j, k, l)= (Cinv(i,j)*Cinv(k,l)*params[0] +
441  params[1]*(Cinv(i,k)*Cinv(j,l) + Cinv(i,l)*Cinv(j,k)))*mult;
442  }
443 
444  SaintVenant_Kirchhoff_hyperelastic_law::SaintVenant_Kirchhoff_hyperelastic_law() {
445  nb_params_ = 2;
446  }
447 
448  scalar_type membrane_elastic_law::strain_energy
449  (const base_matrix & /* E */, const base_vector & /* params */, scalar_type) const {
450  // to be done if needed
451  GMM_ASSERT1(false, "To be done");
452  return 0;
453  }
454 
455  void membrane_elastic_law::sigma
456  (const base_matrix &E, base_matrix &result,const base_vector &params, scalar_type det_trans) const {
457  // should be optimized, maybe deriving sigma from strain energy
458  base_tensor tt(2,2,2,2);
459  size_type N = (gmm::mat_nrows(E) > 2)? 2 : gmm::mat_nrows(E);
460  grad_sigma(E,tt,params, det_trans);
461  for (size_type i = 0; i < N; ++i)
462  for (size_type j = 0; j < N; ++j) {
463  result(i,j)=0.0;
464  for (size_type k = 0; k < N; ++k)
465  for (size_type l = 0; l < N; ++l)
466  result(i,j)+=tt(i,j,k,l)*E(k,l);
467  }
468  // add pretension in X' direction
469  if(params[4]!=0) result(0,0)+=params[4];
470  // add pretension in Y' direction
471  if(params[5]!=0) result(1,1)+=params[5];
472  // cout<<"sigma="<<result<<endl;
473  }
474 
475  void membrane_elastic_law::grad_sigma
476  (const base_matrix & /* E */, base_tensor &result,
477  const base_vector &params, scalar_type) const {
478  // to be optimized!!
479  std::fill(result.begin(), result.end(), scalar_type(0));
480  scalar_type poisonXY=params[0]*params[1]/params[2]; //Ex*vYX=Ey*vXY
481  //if no G entered, compute G=E/(2*(1+v))
482  scalar_type G=( params[3] == 0) ? params[0]/(2*(1+params[1])) : params[3];
483  std::fill(result.begin(), result.end(), scalar_type(0));
484  result(0,0,0,0) = params[0]/(1-params[1]*poisonXY);
485  // result(0,0,0,1) = 0;
486  // result(0,0,1,0) = 0;
487  result(0,0,1,1) = params[1]*params[0]/(1-params[1]*poisonXY);
488  result(1,1,0,0) = params[1]*params[0]/(1-params[1]*poisonXY);
489  // result(1,1,0,1) = 0;
490  // result(1,1,1,0) = 0;
491  result(1,1,1,1) = params[2]/(1-params[1]*poisonXY);
492  // result(0,1,0,0) = 0;
493  result(0,1,0,1) = G;
494  result(0,1,1,0) = G;
495  // result(0,1,1,1) = 0;
496  // result(1,0,0,0) = 0;
497  result(1,0,0,1) = G;
498  result(1,0,1,0) = G;
499  // result(1,0,1,1) = 0;
500  }
501 
502  scalar_type Mooney_Rivlin_hyperelastic_law::strain_energy
503  (const base_matrix &E, const base_vector &params
504  ,scalar_type det_trans) const {
505  //shouldn't negative det_trans be handled here???
506  if (compressible && det_trans <= scalar_type(0)) return 1e200;
507  size_type N = gmm::mat_nrows(E);
508  GMM_ASSERT1(N == 3, "Mooney Rivlin hyperelastic law only defined "
509  "on dimension 3, sorry");
510  base_matrix C = E;
511  gmm::scale(C, scalar_type(2));
512  gmm::add(gmm::identity_matrix(), C);
513  compute_invariants ci(C);
514  size_type i=0;
515  scalar_type C1 = params[i++]; // C10
516  scalar_type W = C1 * (ci.j1() - scalar_type(3));
517  if (!neohookean) {
518  scalar_type C2 = params[i++]; // C01
519  W += C2 * (ci.j2() - scalar_type(3));
520  }
521  if (compressible) {
522  scalar_type D1 = params[i++];
523  W += D1 * gmm::sqr(sqrt(gmm::abs(ci.i3())) - scalar_type(1));
524  }
525  return W;
526  }
527 
528  void Mooney_Rivlin_hyperelastic_law::sigma
529  (const base_matrix &E, base_matrix &result,
530  const base_vector &params, scalar_type det_trans) const {
531  size_type N = gmm::mat_nrows(E);
532  GMM_ASSERT1(N == 3, "Mooney Rivlin hyperelastic law only defined "
533  "on dimension 3, sorry");
534  base_matrix C = E;
535  gmm::scale(C, scalar_type(2));
536  gmm::add(gmm::identity_matrix(), C);
537  compute_invariants ci(C);
538 
539  size_type i=0;
540  scalar_type C1 = params[i++]; // C10
541  gmm::copy(gmm::scaled(ci.grad_j1(), scalar_type(2) * C1), result);
542  if (!neohookean) {
543  scalar_type C2 = params[i++]; // C01
544  gmm::add(gmm::scaled(ci.grad_j2(), scalar_type(2) * C2), result);
545  }
546  if (compressible) {
547  scalar_type D1 = params[i++];
548  scalar_type di3 = D1 - D1 / sqrt(gmm::abs(ci.i3()));
549  gmm::add(gmm::scaled(ci.grad_i3(), scalar_type(2) * di3), result);
550 
551 // shouldn't negative det_trans be handled here???
552  if (det_trans <= scalar_type(0))
553  gmm::add(gmm::scaled(C, 1e200), result);
554  }
555  }
556 
557  void Mooney_Rivlin_hyperelastic_law::grad_sigma
558  (const base_matrix &E, base_tensor &result,
559  const base_vector &params, scalar_type) const {
560  size_type N = gmm::mat_nrows(E);
561  GMM_ASSERT1(N == 3, "Mooney Rivlin hyperelastic law only defined "
562  "on dimension 3, sorry");
563  base_matrix C = E;
564  gmm::scale(C, scalar_type(2));
565  gmm::add(gmm::identity_matrix(), C);
566  compute_invariants ci(C);
567 
568  size_type i=0;
569  scalar_type C1 = params[i++]; // C10
570  gmm::copy(gmm::scaled(ci.sym_grad_grad_j1().as_vector(),
571  scalar_type(4)*C1), result.as_vector());
572  if (!neohookean) {
573  scalar_type C2 = params[i++]; // C01
574  gmm::add(gmm::scaled(ci.sym_grad_grad_j2().as_vector(),
575  scalar_type(4)*C2), result.as_vector());
576  }
577  if (compressible) {
578  scalar_type D1 = params[i++];
579  scalar_type di3 = D1 - D1 / sqrt(gmm::abs(ci.i3()));
580  gmm::add(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
581  scalar_type(4)*di3), result.as_vector());
582 
583  // second derivatives of W with respect to the third invariant
584  scalar_type A22 = D1 / (scalar_type(2) * pow(gmm::abs(ci.i3()), 1.5));
585  const base_matrix &di = ci.grad_i3();
586  for (size_type l1 = 0; l1 < N; ++l1)
587  for (size_type l2 = 0; l2 < N; ++l2)
588  for (size_type l3 = 0; l3 < N; ++l3)
589  for (size_type l4 = 0; l4 < N; ++l4)
590  result(l1, l2, l3, l4) +=
591  scalar_type(4) * A22 * di(l1, l2) * di(l3, l4);
592  }
593 
594 // GMM_ASSERT1(check_symmetry(result) == 7,
595 // "Fourth order tensor not symmetric : " << result);
596  }
597 
598  Mooney_Rivlin_hyperelastic_law::Mooney_Rivlin_hyperelastic_law
599  (bool compressible_, bool neohookean_)
600  : compressible(compressible_), neohookean(neohookean_)
601  {
602  nb_params_ = 2;
603  if (compressible) ++nb_params_; // D1 != 0
604  if (neohookean) --nb_params_; // C2 == 0
605 
606  }
607 
608 
609 
610 
611  scalar_type Neo_Hookean_hyperelastic_law::strain_energy
612  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
613  if (det_trans <= scalar_type(0)) return 1e200;
614  size_type N = gmm::mat_nrows(E);
615  GMM_ASSERT1(N == 3, "Neo Hookean hyperelastic law only defined "
616  "on dimension 3, sorry");
617  base_matrix C = E;
618  gmm::scale(C, scalar_type(2));
619  gmm::add(gmm::identity_matrix(), C);
620  compute_invariants ci(C);
621 
622  scalar_type lambda = params[0];
623  scalar_type mu = params[1];
624  scalar_type logi3 = log(ci.i3());
625  scalar_type W = mu/2 * (ci.i1() - scalar_type(3) - logi3);
626  if (bonet)
627  W += lambda/8 * gmm::sqr(logi3);
628  else // Wriggers
629  W += lambda/4 * (ci.i3() - scalar_type(1) - logi3);
630 
631  return W;
632  }
633 
634  void Neo_Hookean_hyperelastic_law::sigma
635  (const base_matrix &E, base_matrix &result,
636  const base_vector &params , scalar_type det_trans) const {
637  size_type N = gmm::mat_nrows(E);
638  GMM_ASSERT1(N == 3, "Neo Hookean hyperelastic law only defined "
639  "on dimension 3, sorry");
640  base_matrix C = E;
641  gmm::scale(C, scalar_type(2));
642  gmm::add(gmm::identity_matrix(), C);
643  compute_invariants ci(C);
644 
645  scalar_type lambda = params[0];
646  scalar_type mu = params[1];
647  gmm::copy(gmm::scaled(ci.grad_i1(), mu), result);
648  if (bonet)
649  gmm::add(gmm::scaled(ci.grad_i3(),
650  (lambda/2 * log(ci.i3()) - mu) / ci.i3()), result);
651  else
652  gmm::add(gmm::scaled(ci.grad_i3(),
653  lambda/2 - lambda/(2*ci.i3()) - mu / ci.i3()), result);
654  if (det_trans <= scalar_type(0))
655  gmm::add(gmm::scaled(C, 1e200), result);
656  }
657 
658  void Neo_Hookean_hyperelastic_law::grad_sigma
659  (const base_matrix &E, base_tensor &result,
660  const base_vector &params, scalar_type) const {
661  size_type N = gmm::mat_nrows(E);
662  GMM_ASSERT1(N == 3, "Neo Hookean hyperelastic law only defined "
663  "on dimension 3, sorry");
664  base_matrix C = E;
665  gmm::scale(C, scalar_type(2));
666  gmm::add(gmm::identity_matrix(), C);
667  compute_invariants ci(C);
668 
669  scalar_type lambda = params[0];
670  scalar_type mu = params[1];
671 
672  scalar_type coeff;
673  if (bonet) {
674  scalar_type logi3 = log(ci.i3());
675  gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
676  (lambda * logi3 - 2*mu) / ci.i3()),
677  result.as_vector());
678  coeff = (lambda + 2 * mu - lambda * logi3) / gmm::sqr(ci.i3());
679  } else {
680  gmm::copy(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
681  lambda - (lambda + 2 * mu) / ci.i3()),
682  result.as_vector());
683  coeff = (lambda + 2 * mu) / gmm::sqr(ci.i3());
684  }
685 
686  const base_matrix &di = ci.grad_i3();
687  for (size_type l1 = 0; l1 < N; ++l1)
688  for (size_type l2 = 0; l2 < N; ++l2)
689  for (size_type l3 = 0; l3 < N; ++l3)
690  for (size_type l4 = 0; l4 < N; ++l4)
691  result(l1, l2, l3, l4) += coeff * di(l1, l2) * di(l3, l4);
692 
693 // GMM_ASSERT1(check_symmetry(result) == 7,
694 // "Fourth order tensor not symmetric : " << result);
695  }
696 
697  Neo_Hookean_hyperelastic_law::Neo_Hookean_hyperelastic_law(bool bonet_)
698  : bonet(bonet_)
699  {
700  nb_params_ = 2;
701  }
702 
703 
704 
705  scalar_type generalized_Blatz_Ko_hyperelastic_law::strain_energy
706  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
707  if (det_trans <= scalar_type(0)) return 1e200;
708  scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
709  scalar_type n = params[4];
710  size_type N = gmm::mat_nrows(E);
711  GMM_ASSERT1(N == 3, "Generalized Blatz Ko hyperelastic law only defined "
712  "on dimension 3, sorry");
713  base_matrix C = E;
714  gmm::scale(C, scalar_type(2));
715  gmm::add(gmm::identity_matrix(), C);
716  compute_invariants ci(C);
717 
718  return pow(a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
719  + c*ci.i2() / ci.i3() + d, n);
720  }
721 
722  void generalized_Blatz_Ko_hyperelastic_law::sigma
723  (const base_matrix &E, base_matrix &result,
724  const base_vector &params, scalar_type det_trans) const {
725  scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
726  scalar_type n = params[4];
727  size_type N = gmm::mat_nrows(E);
728  GMM_ASSERT1(N == 3, "Generalized Blatz Ko hyperelastic law only defined "
729  "on dimension 3, sorry");
730  base_matrix C = E;
731  gmm::scale(C, scalar_type(2));
732  gmm::add(gmm::identity_matrix(), C);
733  compute_invariants ci(C);
734 
735  scalar_type z = a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
736  + c*ci.i2() / ci.i3() + d;
737  scalar_type nz = n * pow(z, n-1.);
738  scalar_type di1 = nz * a;
739  scalar_type di2 = nz * c / ci.i3();
740  scalar_type di3 = nz *
741  (b / (2. * sqrt(gmm::abs(ci.i3()))) - c * ci.i2() / gmm::sqr(ci.i3()));
742 
743  gmm::copy(gmm::scaled(ci.grad_i1(), di1 * 2.0), result);
744  gmm::add(gmm::scaled(ci.grad_i2(), di2 * 2.0), result);
745  gmm::add(gmm::scaled(ci.grad_i3(), di3 * 2.0), result);
746  if (det_trans <= scalar_type(0))
747  gmm::add(gmm::scaled(C, 1e200), result);
748 
749  }
750 
751  void generalized_Blatz_Ko_hyperelastic_law::grad_sigma
752  (const base_matrix &E, base_tensor &result,
753  const base_vector &params, scalar_type) const {
754  scalar_type a = params[0], b = params[1], c = params[2], d = params[3];
755  scalar_type n = params[4];
756  size_type N = gmm::mat_nrows(E);
757  GMM_ASSERT1(N == 3, "Generalized Blatz Ko hyperelastic law only defined "
758  "on dimension 3, sorry");
759  base_matrix C = E;
760  gmm::scale(C, scalar_type(2));
761  gmm::add(gmm::identity_matrix(), C);
762  compute_invariants ci(C);
763 
764 
765  scalar_type z = a*ci.i1() + b*sqrt(gmm::abs(ci.i3()))
766  + c*ci.i2() / ci.i3() + d;
767  scalar_type nz = n * pow(z, n-1.);
768  scalar_type di1 = nz * a;
769  scalar_type di2 = nz * c / ci.i3();
770  scalar_type y = (b / (2. * sqrt(gmm::abs(ci.i3()))) - c * ci.i2() / gmm::sqr(ci.i3()));
771  scalar_type di3 = nz * y;
772 
773  gmm::copy(gmm::scaled(ci.sym_grad_grad_i1().as_vector(),
774  scalar_type(4)*di1), result.as_vector());
775  gmm::add(gmm::scaled(ci.sym_grad_grad_i2().as_vector(),
776  scalar_type(4)*di2), result.as_vector());
777  gmm::add(gmm::scaled(ci.sym_grad_grad_i3().as_vector(),
778  scalar_type(4)*di3), result.as_vector());
779 
780  scalar_type nnz = n * (n-1.) * pow(z, n-2.);
781  base_matrix A(3, 3); // second derivatives of W with respect to invariants
782  A(0, 0) = nnz * a * a;
783  A(1, 0) = A(0, 1) = nnz * a * c / ci.i3();
784  A(2, 0) = A(0, 2) = nnz * a * y;
785  A(1, 1) = nnz * c * c / gmm::sqr(ci.i3());
786  A(2, 1) = A(1, 2) = nnz * y * c / ci.i3() - nz * c / gmm::sqr(ci.i3());
787  A(2, 2) = nnz * y * y + nz * (2. * c * ci.i2() / pow(ci.i3(), 3.) - b / (4. * pow(ci.i3(), 1.5)));
788 
789  typedef const base_matrix * pointer_base_matrix__;
790  pointer_base_matrix__ di[3];
791  di[0] = &(ci.grad_i1());
792  di[1] = &(ci.grad_i2());
793  di[2] = &(ci.grad_i3());
794 
795  for (size_type j = 0; j < N; ++j)
796  for (size_type k = 0; k < N; ++k) {
797  for (size_type l1 = 0; l1 < N; ++l1)
798  for (size_type l2 = 0; l2 < N; ++l2)
799  for (size_type l3 = 0; l3 < N; ++l3)
800  for (size_type l4 = 0; l4 < N; ++l4)
801  result(l1, l2, l3, l4)
802  += 4. * A(j, k) * (*di[j])(l1, l2) * (*di[k])(l3, l4);
803  }
804 
805 // GMM_ASSERT1(check_symmetry(result) == 7,
806 // "Fourth order tensor not symmetric : " << result);
807  }
808 
809  generalized_Blatz_Ko_hyperelastic_law::generalized_Blatz_Ko_hyperelastic_law() {
810  nb_params_ = 5;
811  base_vector V(5);
812  V[0] = 1.0; V[1] = 1.0, V[2] = 1.5; V[3] = -0.5; V[4] = 1.5;
813  }
814 
815 
816  scalar_type Ciarlet_Geymonat_hyperelastic_law::strain_energy
817  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
818  if (det_trans <= scalar_type(0)) return 1e200;
819  size_type N = gmm::mat_nrows(E);
820  scalar_type a = params[2];
821  scalar_type b = params[1]/scalar_type(2) - params[2];
822  scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
823  + params[2];
824  scalar_type d = params[0]/scalar_type(2) + params[1];
825  scalar_type e = -(scalar_type(3)*(a+b) + c);
826  base_matrix C(N, N);
827  gmm::copy(gmm::scaled(E, scalar_type(2)), C);
828  gmm::add(gmm::identity_matrix(), C);
829  scalar_type det = bgeot::lu_det(&(*(C.begin())), N);
830  return a * gmm::mat_trace(C)
831  + b * (gmm::sqr(gmm::mat_trace(C)) -
832  gmm::mat_euclidean_norm_sqr(C))/scalar_type(2)
833  + c * det - d * log(det) / scalar_type(2) + e;
834  }
835 
836  void Ciarlet_Geymonat_hyperelastic_law::sigma
837  (const base_matrix &E, base_matrix &result, const base_vector &params, scalar_type det_trans) const {
838  size_type N = gmm::mat_nrows(E);
839  scalar_type a = params[2];
840  scalar_type b = params[1]/scalar_type(2) - params[2];
841  scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
842  + params[2];
843  scalar_type d = params[0]/scalar_type(2) + params[1];
844  base_matrix C(N, N);
845  if (a > params[1]/scalar_type(2)
846  || a < params[1]/scalar_type(2) - params[0]/scalar_type(4) || a < 0)
847  GMM_WARNING1("Inconsistent third parameter for Ciarlet-Geymonat "
848  "hyperelastic law");
849  gmm::copy(gmm::scaled(E, scalar_type(2)), C);
850  gmm::add(gmm::identity_matrix(), C);
851  gmm::copy(gmm::identity_matrix(), result);
852  gmm::scale(result, scalar_type(2) * (a + b * gmm::mat_trace(C)));
853  gmm::add(gmm::scaled(C, -scalar_type(2) * b), result);
854  if (det_trans <= scalar_type(0))
855  gmm::add(gmm::scaled(C, 1e200), result);
856  else {
857  scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
858  gmm::add(gmm::scaled(C, scalar_type(2) * c * det - d), result);
859  }
860  }
861 
862  void Ciarlet_Geymonat_hyperelastic_law::grad_sigma
863  (const base_matrix &E, base_tensor &result,const base_vector &params, scalar_type) const {
864  size_type N = gmm::mat_nrows(E);
865  // scalar_type a = params[2];
866  scalar_type b2 = params[1] - params[2]*scalar_type(2); // b*2
867  scalar_type c = params[0]/scalar_type(4) - params[1]/scalar_type(2)
868  + params[2];
869  scalar_type d = params[0]/scalar_type(2) + params[1];
870  base_matrix C(N, N);
871  gmm::copy(gmm::scaled(E, scalar_type(2)), C);
872  gmm::add(gmm::identity_matrix(), C);
873  scalar_type det = bgeot::lu_inverse(&(*(C.begin())), N);
874  std::fill(result.begin(), result.end(), scalar_type(0));
875  for (size_type i = 0; i < N; ++i)
876  for (size_type j = 0; j < N; ++j) {
877  result(i, i, j, j) += 2*b2;
878  result(i, j, i, j) -= b2;
879  result(i, j, j, i) -= b2;
880  for (size_type k = 0; k < N; ++k)
881  for (size_type l = 0; l < N; ++l)
882  result(i, j, k, l) +=
883  (C(i, k)*C(l, j) + C(i, l)*C(k, j)) * (d-scalar_type(2)*det*c)
884  + (C(i, j) * C(k, l)) * det*c*scalar_type(4);
885  }
886 
887 // GMM_ASSERT1(check_symmetry(result) == 7,
888 // "Fourth order tensor not symmetric : " << result);
889  }
890 
891 
892  int levi_civita(int i, int j, int k) {
893  int ii=i+1;
894  int jj=j+1;
895  int kk=k+1; //i,j,k from 0 to 2 !
896  return static_cast<int>
897  (int(- 1)*(static_cast<int>(pow(double(ii-jj),2.))%3)
898  * (static_cast<int> (pow(double(ii-kk),double(2)))%3 )
899  * (static_cast<int> (pow(double(jj-kk),double(2)))%3)
900  * (pow(double(jj-(ii%3))-double(0.5),double(2))-double(1.25)));
901  }
902 
903 
904 
905  scalar_type plane_strain_hyperelastic_law::strain_energy
906  (const base_matrix &E, const base_vector &params, scalar_type det_trans) const {
907  GMM_ASSERT1(gmm::mat_nrows(E) == 2, "Plane strain law is for 2D only.");
908  base_matrix E3D(3,3);
909  E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
910  return pl->strain_energy(E3D, params, det_trans);
911  }
912 
913  void plane_strain_hyperelastic_law::sigma
914  (const base_matrix &E, base_matrix &result,const base_vector &params, scalar_type det_trans) const {
915  GMM_ASSERT1(gmm::mat_nrows(E) == 2, "Plane strain law is for 2D only.");
916  base_matrix E3D(3,3), result3D(3,3);
917  E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
918  pl->sigma(E3D, result3D, params, det_trans);
919  result(0,0) = result3D(0,0); result(1,0) = result3D(1,0);
920  result(0,1) = result3D(0,1); result(1,1) = result3D(1,1);
921  }
922 
923  void plane_strain_hyperelastic_law::grad_sigma
924  (const base_matrix &E, base_tensor &result,const base_vector &params, scalar_type det_trans) const {
925  GMM_ASSERT1(gmm::mat_nrows(E) == 2, "Plane strain law is for 2D only.");
926  base_matrix E3D(3,3);
927  base_tensor result3D(3,3,3,3);
928  E3D(0,0)=E(0,0); E3D(1,0)=E(1,0); E3D(0,1)=E(0,1); E3D(1,1)=E(1,1);
929  pl->grad_sigma(E3D, result3D, params, det_trans);
930  result(0,0,0,0) = result3D(0,0,0,0); result(1,0,0,0) = result3D(1,0,0,0);
931  result(0,1,0,0) = result3D(0,1,0,0); result(1,1,0,0) = result3D(1,1,0,0);
932  result(0,0,1,0) = result3D(0,0,1,0); result(1,0,1,0) = result3D(1,0,1,0);
933  result(0,1,1,0) = result3D(0,1,1,0); result(1,1,1,0) = result3D(1,1,1,0);
934  result(0,0,0,1) = result3D(0,0,0,1); result(1,0,0,1) = result3D(1,0,0,1);
935  result(0,1,0,1) = result3D(0,1,0,1); result(1,1,0,1) = result3D(1,1,0,1);
936  result(0,0,1,1) = result3D(0,0,1,1); result(1,0,1,1) = result3D(1,0,1,1);
937  result(0,1,1,1) = result3D(0,1,1,1); result(1,1,1,1) = result3D(1,1,1,1);
938  }
939 
940 
941 
942 
943 
944 
945 
946  //=========================================================================
947  //
948  // Nonlinear elasticity (old) Brick
949  //
950  //=========================================================================
951 
952  struct nonlinear_elasticity_brick : public virtual_brick {
953 
954  phyperelastic_law AHL;
955 
956  virtual void asm_real_tangent_terms(const model &md, size_type /* ib */,
957  const model::varnamelist &vl,
958  const model::varnamelist &dl,
959  const model::mimlist &mims,
960  model::real_matlist &matl,
961  model::real_veclist &vecl,
962  model::real_veclist &,
963  size_type region,
964  build_version version) const {
965  GMM_ASSERT1(mims.size() == 1,
966  "Nonlinear elasticity brick need a single mesh_im");
967  GMM_ASSERT1(vl.size() == 1,
968  "Nonlinear elasticity brick need a single variable");
969  GMM_ASSERT1(dl.size() == 1,
970  "Wrong number of data for nonlinear elasticity brick, "
971  << dl.size() << " should be 1 (vector).");
972  GMM_ASSERT1(matl.size() == 1, "Wrong number of terms for nonlinear "
973  "elasticity brick");
974 
975  const model_real_plain_vector &u = md.real_variable(vl[0]);
976  const mesh_fem &mf_u = *(md.pmesh_fem_of_variable(vl[0]));
977 
978  const mesh_fem *mf_params = md.pmesh_fem_of_variable(dl[0]);
979  const model_real_plain_vector &params = md.real_variable(dl[0]);
980  const mesh_im &mim = *mims[0];
981 
982  size_type sl = gmm::vect_size(params);
983  if (mf_params) sl = sl * mf_params->get_qdim() / mf_params->nb_dof();
984  GMM_ASSERT1(sl == AHL->nb_params(), "Wrong number of coefficients for the "
985  "nonlinear constitutive elastic law");
986 
987  mesh_region rg(region);
988  mf_u.linked_mesh().intersect_with_mpi_region(rg);
989 
990  if (version & model::BUILD_MATRIX) {
991  gmm::clear(matl[0]);
992  GMM_TRACE2("Nonlinear elasticity stiffness matrix assembly");
994  (matl[0], mim, mf_u, u, mf_params, params, *AHL, rg);
995  }
996 
997 
998  if (version & model::BUILD_RHS) {
999  asm_nonlinear_elasticity_rhs(vecl[0], mim,
1000  mf_u, u, mf_params, params, *AHL, rg);
1001  gmm::scale(vecl[0], scalar_type(-1));
1002  }
1003 
1004  }
1005 
1006  nonlinear_elasticity_brick(const phyperelastic_law &AHL_)
1007  : AHL(AHL_) {
1008  set_flags("Nonlinear elasticity brick", false /* is linear*/,
1009  true /* is symmetric */, true /* is coercive */,
1010  true /* is real */, false /* is complex */);
1011  }
1012 
1013  };
1014 
1015  //=========================================================================
1016  // Add a nonlinear elasticity brick.
1017  //=========================================================================
1018 
1019  // Deprecated brick
1021  (model &md, const mesh_im &mim, const std::string &varname,
1022  const phyperelastic_law &AHL, const std::string &dataname,
1023  size_type region) {
1024  pbrick pbr = std::make_shared<nonlinear_elasticity_brick>(AHL);
1025 
1026  model::termlist tl;
1027  tl.push_back(model::term_description(varname, varname, true));
1028  model::varnamelist dl(1, dataname);
1029  model::varnamelist vl(1, varname);
1030  return md.add_brick(pbr, vl, dl, tl, model::mimlist(1,&mim), region);
1031  }
1032 
1033  //=========================================================================
1034  // Von Mises or Tresca stress computation.
1035  //=========================================================================
1036 
1037  void compute_Von_Mises_or_Tresca(model &md,
1038  const std::string &varname,
1039  const phyperelastic_law &AHL,
1040  const std::string &dataname,
1041  const mesh_fem &mf_vm,
1042  model_real_plain_vector &VM,
1043  bool tresca) {
1044  GMM_ASSERT1(gmm::vect_size(VM) == mf_vm.nb_dof(),
1045  "The vector has not the good size");
1046  const mesh_fem &mf_u = md.mesh_fem_of_variable(varname);
1047  const model_real_plain_vector &u = md.real_variable(varname);
1048  const mesh_fem *mf_params = md.pmesh_fem_of_variable(dataname);
1049  const model_real_plain_vector &params = md.real_variable(dataname);
1050 
1051  size_type sl = gmm::vect_size(params);
1052  if (mf_params) sl = sl * mf_params->get_qdim() / mf_params->nb_dof();
1053  GMM_ASSERT1(sl == AHL->nb_params(), "Wrong number of coefficients for "
1054  "the nonlinear constitutive elastic law");
1055 
1056  unsigned N = unsigned(mf_u.linked_mesh().dim());
1057  unsigned NP = unsigned(AHL->nb_params()), NFem = mf_u.get_qdim();
1058  model_real_plain_vector GRAD(mf_vm.nb_dof()*NFem*N);
1059  model_real_plain_vector PARAMS(mf_vm.nb_dof()*NP);
1060  if (mf_params) interpolation(*mf_params, mf_vm, params, PARAMS);
1061  compute_gradient(mf_u, mf_vm, u, GRAD);
1062  base_matrix E(N, N), gradphi(NFem,N),gradphit(N,NFem), Id(N, N),
1063  sigmahathat(N,N),aux(NFem,N), sigma(NFem,NFem),
1064  IdNFem(NFem, NFem);
1065  base_vector p(NP);
1066  if (!mf_params) gmm::copy(params, p);
1067  base_vector eig(NFem);
1068  base_vector ez(NFem); // vector normal at deformed surface, (ex X ey)
1069  double normEz(0); //norm of ez
1070  gmm::copy(gmm::identity_matrix(), Id);
1071  gmm::copy(gmm::identity_matrix(), IdNFem);
1072  for (size_type i = 0; i < mf_vm.nb_dof(); ++i) {
1073  gmm::resize(gradphi,NFem,N);
1074  std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
1075  gradphit.begin());
1076  gmm::copy(gmm::transposed(gradphit),gradphi);
1077  for (unsigned int alpha = 0; alpha <N; ++alpha)
1078  gradphi(alpha, alpha)+=1;
1079  gmm::mult(gmm::transposed(gradphi), gradphi, E);
1080  gmm::add(gmm::scaled(Id, -scalar_type(1)), E);
1081  gmm::scale(E, scalar_type(1)/scalar_type(2));
1082  if (mf_params)
1083  gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*NP,NP)), p);
1084  AHL->sigma(E, sigmahathat, p, scalar_type(1));
1085  if (NFem == 3 && N == 2) {
1086  //jyh : compute ez, normal on deformed surface
1087  for (unsigned int l = 0; l <NFem; ++l) {
1088  ez[l]=0;
1089  for (unsigned int m = 0; m <NFem; ++m)
1090  for (unsigned int n = 0; n <NFem; ++n){
1091  ez[l]+=levi_civita(l,m,n)*gradphi(m,0)*gradphi(n,1);
1092  }
1093  normEz= gmm::vect_norm2(ez);
1094  }
1095  //jyh : end compute ez
1096  }
1097  gmm::mult(gradphi, sigmahathat, aux);
1098  gmm::mult(aux, gmm::transposed(gradphi), sigma);
1099 
1100  /* jyh : complete gradphi for virtual 3rd dim (perpendicular to
1101  deformed surface, same thickness) */
1102  if (NFem == 3 && N == 2) {
1103  gmm::resize(gradphi,NFem,NFem);
1104  for (unsigned int ll = 0; ll <NFem; ++ll)
1105  for (unsigned int ii = 0; ii <NFem; ++ii)
1106  for (unsigned int jj = 0; jj <NFem; ++jj)
1107  gradphi(ll,2)+=(levi_civita(ll,ii,jj)*gradphi(ii,0)
1108  *gradphi(jj,1))/normEz;
1109  //jyh : end complete graphi
1110  }
1111 
1112  gmm::scale(sigma, scalar_type(1) / bgeot::lu_det(&(*(gradphi.begin())), NFem));
1113 
1114  if (!tresca) {
1115  /* von mises: norm(deviator(sigma)) */
1116  gmm::add(gmm::scaled(IdNFem, -gmm::mat_trace(sigma) / NFem), sigma);
1117 
1118  //jyh : von mises stress=sqrt(3/2)* norm(sigma) ?
1119  VM[i] = sqrt(3.0/2)*gmm::mat_euclidean_norm(sigma);
1120  } else {
1121  /* else compute the tresca criterion */
1122  //jyh : to be adapted for membrane if necessary
1123  gmm::symmetric_qr_algorithm(sigma, eig);
1124  std::sort(eig.begin(), eig.end());
1125  VM[i] = eig.back() - eig.front();
1126  }
1127  }
1128  }
1129 
1130 
1131  void compute_sigmahathat(model &md,
1132  const std::string &varname,
1133  const phyperelastic_law &AHL,
1134  const std::string &dataname,
1135  const mesh_fem &mf_sigma,
1136  model_real_plain_vector &SIGMA) {
1137  const mesh_fem &mf_u = md.mesh_fem_of_variable(varname);
1138  const model_real_plain_vector &u = md.real_variable(varname);
1139  const mesh_fem *mf_params = md.pmesh_fem_of_variable(dataname);
1140  const model_real_plain_vector &params = md.real_variable(dataname);
1141 
1142  size_type sl = gmm::vect_size(params);
1143  if (mf_params) sl = sl * mf_params->get_qdim() / mf_params->nb_dof();
1144  GMM_ASSERT1(sl == AHL->nb_params(), "Wrong number of coefficients for "
1145  "the nonlinear constitutive elastic law");
1146 
1147  unsigned N = unsigned(mf_u.linked_mesh().dim());
1148  unsigned NP = unsigned(AHL->nb_params()), NFem = mf_u.get_qdim();
1149  GMM_ASSERT1(mf_sigma.nb_dof() > 0, "Bad mf_sigma");
1150  size_type qqdim = mf_sigma.get_qdim();
1151  size_type ratio = N*N / qqdim;
1152 
1153  GMM_ASSERT1(((ratio == 1) || (ratio == N*N)) &&
1154  (gmm::vect_size(SIGMA) == mf_sigma.nb_dof()*ratio),
1155  "The vector has not the good size");
1156 
1157  model_real_plain_vector GRAD(mf_sigma.nb_dof()*ratio*NFem/N);
1158  model_real_plain_vector PARAMS(mf_sigma.nb_dof()*NP);
1159 
1160 
1161  getfem::mesh_trans_inv mti(mf_sigma.linked_mesh());
1162  if (mf_params) {
1163  for (size_type i = 0; i < mf_sigma.nb_dof(); ++i)
1164  mti.add_point(mf_sigma.point_of_basic_dof(i));
1165  interpolation(*mf_params, mti, params, PARAMS);
1166  }
1167 
1168  compute_gradient(mf_u, mf_sigma, u, GRAD);
1169  base_matrix E(N, N), gradphi(NFem,N),gradphit(N,NFem), Id(N, N),
1170  sigmahathat(N,N),aux(NFem,N), sigma(NFem,NFem),
1171  IdNFem(NFem, NFem);
1172 
1173 
1174  base_vector p(NP);
1175  if (!mf_params) gmm::copy(params, p);
1176  base_vector eig(NFem);
1177  base_vector ez(NFem); // vector normal at deformed surface, (ex X ey)
1178  gmm::copy(gmm::identity_matrix(), Id);
1179  gmm::copy(gmm::identity_matrix(), IdNFem);
1180 
1181  // cout << "GRAD = " << GRAD << endl;
1182  // GMM_ASSERT1(false, "oops");
1183 
1184  for (size_type i = 0; i < mf_sigma.nb_dof()/qqdim; ++i) {
1185 // cout << "GRAD size = " << gmm::vect_size(GRAD) << endl;
1186 // cout << "i = " << i << endl;
1187 // cout << "i*NFem*N = " << i*NFem*N << endl;
1188 
1189  gmm::resize(gradphi,NFem,N);
1190  std::copy(GRAD.begin()+i*NFem*N, GRAD.begin()+(i+1)*NFem*N,
1191  gradphit.begin());
1192  // cout << "GRAD = " << gradphit << endl;
1193  gmm::copy(gmm::transposed(gradphit),gradphi);
1194  for (unsigned int alpha = 0; alpha <N; ++alpha)
1195  gradphi(alpha, alpha) += scalar_type(1);
1196  gmm::mult(gmm::transposed(gradphi), gradphi, E);
1197  gmm::add(gmm::scaled(Id, -scalar_type(1)), E);
1198  gmm::scale(E, scalar_type(1)/scalar_type(2));
1199  if (mf_params)
1200  gmm::copy(gmm::sub_vector(PARAMS, gmm::sub_interval(i*ratio*NP,NP)),p);
1201  // cout << "E = " << E << endl;
1202  AHL->sigma(E, sigmahathat, p, scalar_type(1));
1203  // cout << "ok" << endl;
1204  std::copy(sigmahathat.begin(), sigmahathat.end(), SIGMA.begin()+i*N*N);
1205  }
1206  }
1207 
1208 
1209 
1210  // ----------------------------------------------------------------------
1211  //
1212  // Nonlinear incompressibility brick
1213  //
1214  // ----------------------------------------------------------------------
1215 
1216  struct nonlinear_incompressibility_brick : public virtual_brick {
1217 
1218  virtual void asm_real_tangent_terms(const model &md, size_type,
1219  const model::varnamelist &vl,
1220  const model::varnamelist &dl,
1221  const model::mimlist &mims,
1222  model::real_matlist &matl,
1223  model::real_veclist &vecl,
1224  model::real_veclist &veclsym,
1225  size_type region,
1226  build_version version) const {
1227 
1228  GMM_ASSERT1(matl.size() == 2, "Wrong number of terms for nonlinear "
1229  "incompressibility brick");
1230  GMM_ASSERT1(dl.size() == 0, "Nonlinear incompressibility brick need no "
1231  "data");
1232  GMM_ASSERT1(mims.size() == 1, "Nonlinear incompressibility brick need a "
1233  "single mesh_im");
1234  GMM_ASSERT1(vl.size() == 2, "Wrong number of variables for nonlinear "
1235  "incompressibility brick");
1236 
1237  const mesh_fem &mf_u = md.mesh_fem_of_variable(vl[0]);
1238  const mesh_fem &mf_p = md.mesh_fem_of_variable(vl[1]);
1239  const model_real_plain_vector &u = md.real_variable(vl[0]);
1240  const model_real_plain_vector &p = md.real_variable(vl[1]);
1241  const mesh_im &mim = *mims[0];
1242  mesh_region rg(region);
1243  mim.linked_mesh().intersect_with_mpi_region(rg);
1244 
1245  if (version & model::BUILD_MATRIX) {
1246  gmm::clear(matl[0]);
1247  gmm::clear(matl[1]);
1248  asm_nonlinear_incomp_tangent_matrix(matl[0], matl[1],
1249  mim, mf_u, mf_p, u, p, rg);
1250  }
1251 
1252  if (version & model::BUILD_RHS) {
1253  asm_nonlinear_incomp_rhs(vecl[0], veclsym[1], mim, mf_u, mf_p,u,p, rg);
1254  gmm::scale(vecl[0], scalar_type(-1));
1255  gmm::scale(veclsym[1], scalar_type(-1));
1256  }
1257  }
1258 
1259 
1260  nonlinear_incompressibility_brick() {
1261  set_flags("Nonlinear incompressibility brick",
1262  false /* is linear*/,
1263  true /* is symmetric */, false /* is coercive */,
1264  true /* is real */, false /* is complex */);
1265  }
1266 
1267 
1268  };
1269 
1271  (model &md, const mesh_im &mim, const std::string &varname,
1272  const std::string &multname, size_type region) {
1273  pbrick pbr = std::make_shared<nonlinear_incompressibility_brick>();
1274  model::termlist tl;
1275  tl.push_back(model::term_description(varname, varname, true));
1276  tl.push_back(model::term_description(varname, multname, true));
1277  model::varnamelist vl(1, varname);
1278  vl.push_back(multname);
1279  model::varnamelist dl;
1280  return md.add_brick(pbr, vl, dl, tl, model::mimlist(1, &mim), region);
1281  }
1282 
1283 
1284 
1285 
1286 
1287  // ----------------------------------------------------------------------
1288  //
1289  // High-level Generic assembly operators
1290  //
1291  // ----------------------------------------------------------------------
1292 
1293 
1294  static void ga_init_scalar_(bgeot::multi_index &mi) { mi.resize(0); }
1295  static void ga_init_square_matrix_(bgeot::multi_index &mi, size_type N)
1296  { mi.resize(2); mi[0] = mi[1] = N; }
1297 
1298 
1299  // Matrix_i2 Operator (second invariant of square matrix of size >= 2)
1300  // (For 2x2 matrices, it is equivalent to det(M))
1301  struct matrix_i2_operator : public ga_nonlinear_operator {
1302  bool result_size(const arg_list &args, bgeot::multi_index &sizes) const {
1303  if (args.size() != 1 || args[0]->sizes().size() != 2
1304  || args[0]->sizes()[0] != args[0]->sizes()[1]) return false;
1305  ga_init_scalar_(sizes);
1306  return true;
1307  }
1308 
1309  // Value : (Trace(M))^2 - Trace(M^2))/2
1310  void value(const arg_list &args, base_tensor &result) const {
1311  size_type N = args[0]->sizes()[0];
1312  const base_tensor &t = *args[0];
1313  scalar_type tr = scalar_type(0);
1314  for (size_type i = 0; i < N; ++i) tr += t[i*N+i];
1315  scalar_type tr2 = scalar_type(0);
1316  for (size_type i = 0; i < N; ++i)
1317  for (size_type j = 0; j < N; ++j)
1318  tr2 += t[i+ j*N] * t[j + i*N];
1319  result[0] = (tr*tr - tr2)/2;
1320  }
1321 
1322  // Derivative : Trace(M)I - M^T
1323  void derivative(const arg_list &args, size_type,
1324  base_tensor &result) const { // to be verified
1325  size_type N = args[0]->sizes()[0];
1326  const base_tensor &t = *args[0];
1327  scalar_type tr = scalar_type(0);
1328  for (size_type i = 0; i < N; ++i) tr += t[i*N+i];
1329  base_tensor::iterator it = result.begin();
1330  for (size_type j = 0; j < N; ++j)
1331  for (size_type i = 0; i < N; ++i, ++it)
1332  *it = ((i == j) ? tr : scalar_type(0)) - t[i*N+j];
1333  GMM_ASSERT1(it == result.end(), "Internal error");
1334  }
1335 
1336  // Second derivative : I@I - \delta_{il}\delta_{jk}
1337  void second_derivative(const arg_list &args, size_type, size_type,
1338  base_tensor &result) const { // To be verified
1339  size_type N = args[0]->sizes()[0];
1340  gmm::clear(result.as_vector());
1341  for (size_type i = 0; i < N; ++i)
1342  for (size_type j = 0; j < N; ++j) {
1343  result[(N+1)*(i+N*N*j)] += scalar_type(1);
1344  result[(N+1)*N*j + i*(N*N*N + 1)] -= scalar_type(1);
1345  }
1346  }
1347  };
1348 
1349 
1350  // Matrix_j1 Operator
1351  struct matrix_j1_operator : public ga_nonlinear_operator {
1352  bool result_size(const arg_list &args, bgeot::multi_index &sizes) const {
1353  if (args.size() != 1 || args[0]->sizes().size() != 2
1354  || args[0]->sizes()[0] != args[0]->sizes()[1]) return false;
1355  ga_init_scalar_(sizes);
1356  return true;
1357  }
1358 
1359  // Value : Trace(M)/(det(M)^1/3)
1360  void value(const arg_list &args, base_tensor &result) const {
1361  size_type N = args[0]->sizes()[0];
1362  base_matrix M(N, N);
1363  gmm::copy(args[0]->as_vector(), M.as_vector());
1364  scalar_type det = bgeot::lu_det(&(*(M.begin())), N);
1365  scalar_type tr = scalar_type(0);
1366  for (size_type i = 0; i < N; ++i) tr += M(i,i);
1367  if (det > 0)
1368  result[0] = tr / pow(det, scalar_type(1)/scalar_type(3));
1369  else
1370  result[0] = 1.E200;
1371  }
1372 
1373  // Derivative : (I-Trace(M)*M^{-T}/3)/(det(M)^1/3)
1374  void derivative(const arg_list &args, size_type,
1375  base_tensor &result) const { // to be verified
1376  size_type N = args[0]->sizes()[0];
1377  base_matrix M(N, N);
1378  gmm::copy(args[0]->as_vector(), M.as_vector());
1379  scalar_type tr = scalar_type(0);
1380  for (size_type i = 0; i < N; ++i) tr += M(i,i);
1381  scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1382  if (det > 0) {
1383  base_tensor::iterator it = result.begin();
1384  for (size_type j = 0; j < N; ++j)
1385  for (size_type i = 0; i < N; ++i, ++it)
1386  *it = (((i == j) ? scalar_type(1) : scalar_type(0))
1387  - tr*M(j,i)/scalar_type(3))
1388  / pow(det, scalar_type(1)/scalar_type(3));
1389  GMM_ASSERT1(it == result.end(), "Internal error");
1390  } else
1391  std::fill(result.begin(), result.end(), 1.E200);
1392  }
1393 
1394  // Second derivative : (-M^{-T}@I + Trace(M)*M^{-T}_{ik}M^{-T}_{lj}
1395  // -I@M^{-T} + Trace(M)*M^{-T}@M^{-T}/3)/(3det(M)^1/3)
1396  void second_derivative(const arg_list &args, size_type, size_type,
1397  base_tensor &result) const { // To be verified
1398  size_type N = args[0]->sizes()[0];
1399  base_matrix M(N, N);
1400  gmm::copy(args[0]->as_vector(), M.as_vector());
1401  scalar_type tr = scalar_type(0);
1402  for (size_type i = 0; i < N; ++i) tr += M(i,i);
1403  scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1404  if (det > 0) {
1405  base_tensor::iterator it = result.begin();
1406  for (size_type l = 0; l < N; ++l)
1407  for (size_type k = 0; k < N; ++k)
1408  for (size_type j = 0; j < N; ++j)
1409  for (size_type i = 0; i < N; ++i, ++it)
1410  *it = (- ((k == l) ? M(j, i) : scalar_type(0))
1411  + tr*M(i,k)*M(l,j)
1412  - ((i == j) ? M(l, k) : scalar_type(0))
1413  + tr*M(j,i)*M(k,l)/ scalar_type(3))
1414  / (scalar_type(3)*pow(det, scalar_type(1)/scalar_type(3)));
1415  GMM_ASSERT1(it == result.end(), "Internal error");
1416  } else
1417  std::fill(result.begin(), result.end(), 1.E200);
1418  }
1419  };
1420 
1421 
1422  // Matrix_j2 Operator
1423  struct matrix_j2_operator : public ga_nonlinear_operator {
1424  bool result_size(const arg_list &args, bgeot::multi_index &sizes) const {
1425  if (args.size() != 1 || args[0]->sizes().size() != 2
1426  || args[0]->sizes()[0] != args[0]->sizes()[1]) return false;
1427  ga_init_scalar_(sizes);
1428  return true;
1429  }
1430 
1431  // Value : i2(M)/(det(M)^2/3)
1432  void value(const arg_list &args, base_tensor &result) const {
1433  size_type N = args[0]->sizes()[0];
1434  base_matrix M(N, N);
1435  gmm::copy(args[0]->as_vector(), M.as_vector());
1436  scalar_type tr = scalar_type(0);
1437  for (size_type i = 0; i < N; ++i) tr += M(i,i);
1438  scalar_type tr2 = scalar_type(0);
1439  for (size_type i = 0; i < N; ++i)
1440  for (size_type j = 0; j < N; ++j)
1441  tr2 += M(i,j)*M(j,i);
1442  scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
1443  scalar_type det = bgeot::lu_det(&(*(M.begin())), N);
1444 
1445  if (det > 0)
1446  result[0] = i2 / pow(det, scalar_type(2)/scalar_type(3));
1447  else
1448  result[0] = 1.E200;
1449  }
1450 
1451  // Derivative : (Trace(M)*I-M^T-2i2(M)M^{-T}/3)/(det(M)^2/3)
1452  void derivative(const arg_list &args, size_type,
1453  base_tensor &result) const { // to be verified
1454  size_type N = args[0]->sizes()[0];
1455  const base_tensor &t = *args[0];
1456  base_matrix M(N, N);
1457  gmm::copy(t.as_vector(), M.as_vector());
1458  scalar_type tr = scalar_type(0);
1459  for (size_type i = 0; i < N; ++i) tr += M(i,i);
1460  scalar_type tr2 = scalar_type(0);
1461  for (size_type i = 0; i < N; ++i)
1462  for (size_type j = 0; j < N; ++j)
1463  tr2 += M(i,j)*M(j,i);
1464  scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
1465  scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1466  base_tensor::iterator it = result.begin();
1467  for (size_type j = 0; j < N; ++j)
1468  for (size_type i = 0; i < N; ++i, ++it)
1469  *it = (((i == j) ? tr : scalar_type(0)) - t[j+N*i]
1470  - scalar_type(2)*i2*M(j,i)/scalar_type(3))
1471  / pow(det, scalar_type(2)/scalar_type(3));
1472  GMM_ASSERT1(it == result.end(), "Internal error");
1473  }
1474 
1475  // Second derivative
1476  void second_derivative(const arg_list &args, size_type, size_type,
1477  base_tensor &result) const { // To be verified
1478  size_type N = args[0]->sizes()[0];
1479  const base_tensor &t = *args[0];
1480  base_matrix M(N, N);
1481  gmm::copy(t.as_vector(), M.as_vector());
1482  scalar_type tr = scalar_type(0);
1483  for (size_type i = 0; i < N; ++i) tr += M(i,i);
1484  scalar_type tr2 = scalar_type(0);
1485  for (size_type i = 0; i < N; ++i)
1486  for (size_type j = 0; j < N; ++j)
1487  tr2 += M(i,j)*M(j,i);
1488  scalar_type i2 = (tr*tr-tr2)/scalar_type(2);
1489  scalar_type det = bgeot::lu_inverse(&(*(M.begin())), N);
1490  base_tensor::iterator it = result.begin();
1491  for (size_type l = 0; l < N; ++l)
1492  for (size_type k = 0; k < N; ++k)
1493  for (size_type j = 0; j < N; ++j)
1494  for (size_type i = 0; i < N; ++i, ++it)
1495  *it = ( ((i==j) ? 1. : 0.) * ((k==l) ? 1. : 0.)
1496  - ((i==l) ? 1. : 0.) * ((k==j) ? 1. : 0.)
1497  - 2.*tr*M(j,i)*((k==l) ? 1. : 0.)/3.
1498  + 2.*tr*M(j,i)*M(l,k)/3.
1499  - 2.*i2*M(i,k)*M(l,j)/3.
1500  - 2.*((tr*((i==j) ? 1. : 0.))-t[j+N*i]
1501  - 2.*i2*M(j,i)/3)*M(l,k)/3.)
1502  / pow(det, scalar_type(2)/scalar_type(3));
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 { // to be verified
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 { // to be verified
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_potential",
2049  std::make_shared<AHL_wrapper_potential>
2050  (std::make_shared<SaintVenant_Kirchhoff_hyperelastic_law>()));
2051  PREDEF_OPERATORS.add_method
2052  ("Plane_Strain_Saint_Venant_Kirchhoff_sigma",
2053  std::make_shared<Saint_Venant_Kirchhoff_sigma>());
2054  PREDEF_OPERATORS.add_method
2055  ("Plane_Strain_Saint_Venant_Kirchhoff_potential",
2056  std::make_shared<AHL_wrapper_potential>
2057  (std::make_shared<SaintVenant_Kirchhoff_hyperelastic_law>()));
2058 
2059  phyperelastic_law gbklaw
2060  = std::make_shared<generalized_Blatz_Ko_hyperelastic_law>();
2061  PREDEF_OPERATORS.add_method
2062  ("Generalized_Blatz_Ko_sigma",
2063  std::make_shared<AHL_wrapper_sigma>(gbklaw));
2064  PREDEF_OPERATORS.add_method
2065  ("Generalized_Blatz_Ko_potential",
2066  std::make_shared<AHL_wrapper_potential>
2067  (std::make_shared<generalized_Blatz_Ko_hyperelastic_law>()));
2068  PREDEF_OPERATORS.add_method
2069  ("Plane_Strain_Generalized_Blatz_Ko_sigma",
2070  std::make_shared<AHL_wrapper_sigma>
2071  (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
2072  PREDEF_OPERATORS.add_method
2073  ("Plane_Strain_Generalized_Blatz_Ko_potential",
2074  std::make_shared<AHL_wrapper_potential>
2075  (std::make_shared<plane_strain_hyperelastic_law>(gbklaw)));
2076 
2077  phyperelastic_law cigelaw
2078  = std::make_shared<Ciarlet_Geymonat_hyperelastic_law>();
2079  PREDEF_OPERATORS.add_method
2080  ("Ciarlet_Geymonat_sigma", std::make_shared<AHL_wrapper_sigma>(cigelaw));
2081  PREDEF_OPERATORS.add_method
2082  ("Ciarlet_Geymonat_potential",
2083  std::make_shared<AHL_wrapper_potential>
2084  (std::make_shared<Ciarlet_Geymonat_hyperelastic_law>()));
2085  PREDEF_OPERATORS.add_method
2086  ("Plane_Strain_Ciarlet_Geymonat_sigma",
2087  std::make_shared<AHL_wrapper_sigma>
2088  (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
2089  PREDEF_OPERATORS.add_method
2090  ("Plane_Strain_Ciarlet_Geymonat_potential",
2091  std::make_shared<AHL_wrapper_potential>
2092  (std::make_shared<plane_strain_hyperelastic_law>(cigelaw)));
2093 
2094  phyperelastic_law morilaw
2095  = std::make_shared<Mooney_Rivlin_hyperelastic_law>();
2096  PREDEF_OPERATORS.add_method
2097  ("Incompressible_Mooney_Rivlin_sigma",
2098  std::make_shared<AHL_wrapper_sigma>(morilaw));
2099  PREDEF_OPERATORS.add_method
2100  ("Incompressible_Mooney_Rivlin_potential",
2101  std::make_shared<AHL_wrapper_potential>
2102  (std::make_shared<Mooney_Rivlin_hyperelastic_law>()));
2103  PREDEF_OPERATORS.add_method
2104  ("Plane_Strain_Incompressible_Mooney_Rivlin_sigma",
2105  std::make_shared<AHL_wrapper_sigma>
2106  (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
2107  PREDEF_OPERATORS.add_method
2108  ("Plane_Strain_Incompressible_Mooney_Rivlin_potential",
2109  std::make_shared<AHL_wrapper_potential>
2110  (std::make_shared<plane_strain_hyperelastic_law>(morilaw)));
2111 
2112  phyperelastic_law cmorilaw
2113  = std::make_shared<Mooney_Rivlin_hyperelastic_law>(true);
2114  PREDEF_OPERATORS.add_method
2115  ("Compressible_Mooney_Rivlin_sigma",
2116  std::make_shared<AHL_wrapper_sigma>(cmorilaw));
2117  PREDEF_OPERATORS.add_method
2118  ("Compressible_Mooney_Rivlin_potential",
2119  std::make_shared<AHL_wrapper_potential>
2120  (std::make_shared<Mooney_Rivlin_hyperelastic_law>(true)));
2121  PREDEF_OPERATORS.add_method
2122  ("Plane_Strain_Compressible_Mooney_Rivlin_sigma",
2123  std::make_shared<AHL_wrapper_sigma>
2124  (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
2125  PREDEF_OPERATORS.add_method
2126  ("Plane_Strain_Compressible_Mooney_Rivlin_potential",
2127  std::make_shared<AHL_wrapper_potential>
2128  (std::make_shared<plane_strain_hyperelastic_law>(cmorilaw)));
2129 
2130  phyperelastic_law ineolaw
2131  = std::make_shared<Mooney_Rivlin_hyperelastic_law>(false, true);
2132  PREDEF_OPERATORS.add_method
2133  ("Incompressible_Neo_Hookean_sigma",
2134  std::make_shared<AHL_wrapper_sigma>(ineolaw));
2135  PREDEF_OPERATORS.add_method
2136  ("Incompressible_Neo_Hookean_potential",
2137  std::make_shared<AHL_wrapper_potential>
2138  (std::make_shared<Mooney_Rivlin_hyperelastic_law>(false,true)));
2139  PREDEF_OPERATORS.add_method
2140  ("Plane_Strain_Incompressible_Neo_Hookean_sigma",
2141  std::make_shared<AHL_wrapper_sigma>
2142  (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
2143  PREDEF_OPERATORS.add_method
2144  ("Plane_Strain_Incompressible_Neo_Hookean_potential",
2145  std::make_shared<AHL_wrapper_potential>
2146  (std::make_shared<plane_strain_hyperelastic_law>(ineolaw)));
2147 
2148  phyperelastic_law cneolaw
2149  = std::make_shared<Mooney_Rivlin_hyperelastic_law>(true, true);
2150  PREDEF_OPERATORS.add_method
2151  ("Compressible_Neo_Hookean_sigma",
2152  std::make_shared<AHL_wrapper_sigma>(cneolaw));
2153  PREDEF_OPERATORS.add_method
2154  ("Compressible_Neo_Hookean_potential",
2155  std::make_shared<AHL_wrapper_potential>
2156  (std::make_shared<Mooney_Rivlin_hyperelastic_law>(true,true)));
2157  PREDEF_OPERATORS.add_method
2158  ("Plane_Strain_Compressible_Neo_Hookean_sigma",
2159  std::make_shared<AHL_wrapper_sigma>
2160  (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
2161  PREDEF_OPERATORS.add_method
2162  ("Plane_Strain_Compressible_Neo_Hookean_potential",
2163  std::make_shared<AHL_wrapper_potential>
2164  (std::make_shared<plane_strain_hyperelastic_law>(cneolaw)));
2165 
2166  phyperelastic_law cneobolaw
2167  = std::make_shared<Neo_Hookean_hyperelastic_law>(true);
2168  PREDEF_OPERATORS.add_method
2169  ("Compressible_Neo_Hookean_Bonet_sigma",
2170  std::make_shared<AHL_wrapper_sigma>(cneobolaw));
2171  PREDEF_OPERATORS.add_method
2172  ("Compressible_Neo_Hookean_Bonet_potential",
2173  std::make_shared<AHL_wrapper_potential>
2174  (std::make_shared<Neo_Hookean_hyperelastic_law>(true)));
2175  PREDEF_OPERATORS.add_method
2176  ("Plane_Strain_Compressible_Neo_Hookean_Bonet_sigma",
2177  std::make_shared<AHL_wrapper_sigma>
2178  (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
2179  PREDEF_OPERATORS.add_method
2180  ("Plane_Strain_Compressible_Neo_Hookean_Bonet_potential",
2181  std::make_shared<AHL_wrapper_potential>
2182  (std::make_shared<plane_strain_hyperelastic_law>(cneobolaw)));
2183 
2184  phyperelastic_law cneocilaw
2185  = std::make_shared<Neo_Hookean_hyperelastic_law>(false);
2186  PREDEF_OPERATORS.add_method
2187  ("Compressible_Neo_Hookean_Ciarlet_sigma",
2188  std::make_shared<AHL_wrapper_sigma>(cneocilaw));
2189  PREDEF_OPERATORS.add_method
2190  ("Compressible_Neo_Hookean_Ciarlet_potential",
2191  std::make_shared<AHL_wrapper_potential>
2192  (std::make_shared<Neo_Hookean_hyperelastic_law>(false)));
2193  PREDEF_OPERATORS.add_method
2194  ("Plane_Strain_Compressible_Neo_Hookean_Ciarlet_sigma",
2195  std::make_shared<AHL_wrapper_sigma>
2196  (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
2197  PREDEF_OPERATORS.add_method
2198  ("Plane_Strain_Compressible_Neo_Hookean_Ciarlet_potential",
2199  std::make_shared<AHL_wrapper_potential>
2200  (std::make_shared<plane_strain_hyperelastic_law>(cneocilaw)));
2201 
2202  return true;
2203  }
2204 
2205  // declared in getfem_generic_assembly.cc
2206  bool predef_operators_nonlinear_elasticity_initialized
2207  = init_predef_operators();
2208 
2209 
2210  std::string adapt_law_name(const std::string &lawname, size_type N) {
2211  std::string adapted_lawname = lawname;
2212 
2213  for (size_type i = 0; i < lawname.size(); ++i)
2214  if (adapted_lawname[i] == ' ') adapted_lawname[i] = '_';
2215 
2216  if (adapted_lawname.compare("SaintVenant_Kirchhoff") == 0) {
2217  adapted_lawname = "Saint_Venant_Kirchhoff";
2218  } else if (adapted_lawname.compare("Saint_Venant_Kirchhoff") == 0) {
2219 
2220  } else if (adapted_lawname.compare("Generalized_Blatz_Ko") == 0) {
2221  if (N == 2) adapted_lawname = "Plane_Strain_" + adapted_lawname;
2222  } else if (adapted_lawname.compare("Ciarlet_Geymonat") == 0) {
2223  if (N == 2) adapted_lawname = "Plane_Strain_" + adapted_lawname;
2224  } else if (adapted_lawname.compare("Incompressible_Mooney_Rivlin") == 0) {
2225  if (N == 2) adapted_lawname = "Plane_Strain_" + adapted_lawname;
2226  } else if (adapted_lawname.compare("Compressible_Mooney_Rivlin") == 0) {
2227  if (N == 2) adapted_lawname = "Plane_Strain_" + adapted_lawname;
2228  } else if (adapted_lawname.compare("Incompressible_Neo_Hookean") == 0) {
2229  if (N == 2) adapted_lawname = "Plane_Strain_" + adapted_lawname;
2230  } else if (adapted_lawname.compare("Compressible_Neo_Hookean") == 0 ||
2231  adapted_lawname.compare("Compressible_Neo_Hookean_Bonet") == 0 ||
2232  adapted_lawname.compare("Compressible_Neo_Hookean_Ciarlet") == 0 ) {
2233  if (N == 2) adapted_lawname = "Plane_Strain_" + adapted_lawname;
2234  } else
2235  GMM_ASSERT1(false, lawname << " is not a known hyperelastic law");
2236 
2237  return adapted_lawname;
2238  }
2239 
2240 
2242  (model &md, const mesh_im &mim, const std::string &lawname,
2243  const std::string &varname, const std::string &params,
2244  size_type region) {
2245  std::string test_varname = "Test_" + sup_previous_and_dot_to_varname(varname);
2246  size_type N = mim.linked_mesh().dim();
2247  GMM_ASSERT1(N >= 2 && N <= 3,
2248  "Finite strain elasticity brick works only in 2D or 3D");
2249 
2250  const mesh_fem *mf = md.pmesh_fem_of_variable(varname);
2251  GMM_ASSERT1(mf, "Finite strain elasticity brick can only be applied on "
2252  "fem variables");
2253  size_type Q = mf->get_qdim();
2254  GMM_ASSERT1(Q == N, "Finite strain elasticity brick can only be applied "
2255  "on a fem variable having the same dimension than the mesh");
2256 
2257  std::string adapted_lawname = adapt_law_name(lawname, N);
2258 
2259  std::string expr = "((Id(meshdim)+Grad_"+varname+")*(" + adapted_lawname
2260  + "_sigma(Grad_"+varname+","+params+"))):Grad_" + test_varname;
2261 
2262  return add_nonlinear_generic_assembly_brick
2263  (md, mim, expr, region, true, false,
2264  "Finite strain elasticity brick for " + adapted_lawname + " law");
2265  }
2266 
2268  (model &md, const mesh_im &mim, const std::string &varname,
2269  const std::string &multname, size_type region) {
2270  std::string test_varname = "Test_" + sup_previous_and_dot_to_varname(varname);
2271  std::string test_multname = "Test_" + sup_previous_and_dot_to_varname(multname);
2272 
2273  std::string expr
2274  = "(" + test_multname+ ")*(1-Det(Id(meshdim)+Grad_" + varname + "))"
2275  + "-(" + multname + ")*(Det(Id(meshdim)+Grad_" + varname + ")"
2276  + "*((Inv(Id(meshdim)+Grad_" + varname + "))':Grad_"
2277  + test_varname + "))" ;
2278  return add_nonlinear_generic_assembly_brick
2279  (md, mim, expr, region, true, false,
2280  "Finite strain incompressibility brick");
2281  }
2282 
2284  (model &md, const std::string &lawname, const std::string &varname,
2285  const std::string &params, const mesh_fem &mf_vm,
2286  model_real_plain_vector &VM, const mesh_region &rg) {
2287  size_type N = mf_vm.linked_mesh().dim();
2288 
2289  std::string adapted_lawname = adapt_law_name(lawname, N);
2290 
2291  std::string expr = "sqrt(3/2)*Norm(Deviator(Cauchy_stress_from_PK2("
2292  + adapted_lawname + "_sigma(Grad_" + varname + "," + params + "),Grad_"
2293  + varname + ")))";
2294  ga_interpolation_Lagrange_fem(md, expr, mf_vm, VM, rg);
2295  }
2296 
2297 
2298 
2299 } /* end of namespace getfem. */
2300 
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
structure used to hold a set of convexes and/or convex faces.
size_type add_nonlinear_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a nonlinear incompressibility term (for large strain elasticity) to the model with respect to the...
Non-linear elasticty and incompressibility bricks.
virtual void cauchy_updated_lagrangian(const base_matrix &F, const base_matrix &E, base_matrix &cauchy_stress, const base_vector &params, scalar_type det_trans) const
True Cauchy stress (for Updated Lagrangian formulation)
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
void compute_gradient(const mesh_fem &mf, const mesh_fem &mf_target, const VECT1 &UU, VECT2 &VV)
Compute the gradient of a field on a getfem::mesh_fem.
Describe an integration method linked to a mesh.
void fill_random(L &l, double cfill)
*/
Definition: gmm_blas.h:154
static T & instance()
Instance from the current thread.
size_type add_finite_strain_elasticity_brick(model &md, const mesh_im &mim, const std::string &lawname, const std::string &varname, const std::string &params, size_type region=size_type(-1))
Add a finite strain elasticity brick to the model with respect to the variable varname (the displacem...
const mesh_fem & mesh_fem_of_variable(const std::string &name) const
Gives the access to the mesh_fem of a variable if any.
void compute_finite_strain_elasticity_Von_Mises(model &md, const std::string &lawname, const std::string &varname, const std::string &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...
``Model&#39;&#39; variables store the variables, the data and the description of a model. ...
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
Model representation in Getfem.
size_type add_brick(pbrick pbr, const varnamelist &varnames, const varnamelist &datanames, const termlist &terms, const mimlist &mims, size_type region)
Add a brick to the model.
virtual dim_type get_qdim() const
Return the Q dimension.
virtual base_node point_of_basic_dof(size_type cv, size_type i) const
Return the geometrical location of a degree of freedom.
A langage for generic assembly of pde boundary value problems.
GEneric Tool for Finite Element Methods.
Describe a finite element method linked to a mesh.
The virtual brick has to be derived to describe real model bricks.
void asm_nonlinear_elasticity_tangent_matrix(const MAT &K_, const mesh_im &mim, const getfem::mesh_fem &mf, const VECT1 &U, const getfem::mesh_fem *mf_data, const VECT2 &PARAMS, const abstract_hyperelastic_law &AHL, const mesh_region &rg=mesh_region::all_convexes())
Tangent matrix for the non-linear elasticity.
linalg_traits< M >::value_type mat_trace(const M &m)
Trace of a matrix.
Definition: gmm_blas.h:528
size_type add_nonlinear_elasticity_brick(model &md, const mesh_im &mim, const std::string &varname, const phyperelastic_law &AHL, const std::string &dataname, size_type region=size_type(-1))
Add a nonlinear (large strain) elasticity term to the model with respect to the variable varname (dep...
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
const mesh & linked_mesh() const
Give a reference to the linked mesh of type mesh.
const mesh_fem * pmesh_fem_of_variable(const std::string &name) const
Gives a pointer to the mesh_fem of a variable if any.
const model_real_plain_vector & real_variable(const std::string &name, size_type niter) const
Gives the access to the vector value of a variable.
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector &params, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian
void resize(M &v, size_type m, size_type n)
*/
Definition: gmm_blas.h:231
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
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...