GetFEM  5.5
getfem_nonlinear_elasticity.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2000-2026 Yves Renard, Julien Pommier
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 /**@file getfem_nonlinear_elasticity.h
32  @author Yves Renard <Yves.Renard@insa-lyon.fr>,
33  @author Julien Pommier <Julien.Pommier@insa-toulouse.fr>
34  @date July 6, 2004.
35  @brief Non-linear elasticty and incompressibility bricks.
36 */
37 #ifndef GETFEM_NONLINEAR_ELASTICITY_H__
38 #define GETFEM_NONLINEAR_ELASTICITY_H__
39 
40 #include "getfem_models.h"
42 #include "getfem_derivatives.h"
43 #include "getfem_interpolation.h"
45 #include "gmm/gmm_inoutput.h"
46 
47 namespace getfem {
48 
49 
50  int check_symmetry(const base_tensor &t);
51 
52  class abstract_hyperelastic_law;
53  typedef std::shared_ptr<const abstract_hyperelastic_law> phyperelastic_law;
54 
55  /** Base class for material law.
56  Inherit from this class to define a new law.
57  */
59  public:
60  mutable int uvflag;
61  size_type nb_params_;
62  phyperelastic_law pl; /* optional reference */
63  void reset_unvalid_flag() const { uvflag = 0; }
64  void inc_unvalid_flag() const { uvflag++; }
65  int get_unvalid_flag() const { return uvflag; }
66  std::string adapted_tangent_term_assembly_fem_data; // should be filled
67  std::string adapted_tangent_term_assembly_cte_data; // to replace the
68  // default assembly
69 
70  virtual scalar_type strain_energy(const base_matrix &E,
71  const base_vector &params,
72  scalar_type det_trans) const = 0;
73  /**True Cauchy stress (for Updated Lagrangian formulation)*/
74  virtual void cauchy_updated_lagrangian(const base_matrix& F,
75  const base_matrix &E,
76  base_matrix &cauchy_stress,
77  const base_vector &params,
78  scalar_type det_trans) const;
79  virtual void sigma(const base_matrix &E, base_matrix &result,
80  const base_vector &params,
81  scalar_type det_trans) const = 0;
82  // the result of grad_sigma has to be completely symmetric.
83  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
84  const base_vector &params, scalar_type det_trans) const = 0;
85 
86  /**cauchy-truesdel tangent moduli, used in updated lagrangian*/
87  virtual void grad_sigma_updated_lagrangian(const base_matrix& F,
88  const base_matrix& E,
89  const base_vector &params,
90  scalar_type det_trans,
91  base_tensor &grad_sigma_ul)const;
92 
93  size_type nb_params() const { return nb_params_; }
94  abstract_hyperelastic_law() { nb_params_ = 0; }
95  virtual ~abstract_hyperelastic_law() {}
96  static void random_E(base_matrix &E);
97  void test_derivatives(size_type N, scalar_type h,
98  const base_vector& param) const;
99  };
100 
101  /** Saint-Venant Kirchhoff hyperelastic law.
102 
103  This is the linear law used in linear elasticity, it is not well
104  suited to large strain. (the convexes may become flat)
105  */
108  /* W = lambda*0.5*trace(E)^2 + mu*tr(E^2) */
109  virtual scalar_type strain_energy(const base_matrix &E,
110  const base_vector &params,
111  scalar_type det_trans) const;
112  /* sigma = lambda*trace(E) + 2 mu * E */
113  virtual void sigma(const base_matrix &E, base_matrix &result,
114  const base_vector &params, scalar_type det_trans) const;
115  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
116  const base_vector &params, scalar_type det_trans) const;
117  virtual void grad_sigma_updated_lagrangian(const base_matrix& F,
118  const base_matrix& E,
119  const base_vector &params,
120  scalar_type det_trans,
121  base_tensor &grad_sigma_ul)const;
123  };
124 
125 
126  /** Linear law for a membrane (plane stress), orthotropic material
127  caracterized by Ex, Ey, vYX and G, with optional orthotropic prestresses.
128  due to Jean-Yves Heddebaut <jyhed@alpino.be>
129  */
132  virtual scalar_type strain_energy(const base_matrix &E,
133  const base_vector &params,
134  scalar_type det_trans) const;
135  virtual void sigma(const base_matrix &E, base_matrix &result,
136  const base_vector &params, scalar_type det_trans) const;
137  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
138  const base_vector &params, scalar_type det_trans) const;
139  membrane_elastic_law() { nb_params_ = 6; }
140  };
141 
142 
143  /** Mooney-Rivlin hyperelastic law
144 
145  To be used for compressible and incompressible problems.
146  Following combinations are possible:
147  not compressible, not neohookean (default): 2 parameters (C1,C2)
148  not compressible, neohookean: 1 parameter (C1)
149  compressible, not neohookean: 3 parameters (C1,C2,D1)
150  compressible, neohookean: 2 parameters (C1,D1)
151  */
153  const bool compressible, neohookean;
154  virtual scalar_type strain_energy(const base_matrix &E,
155  const base_vector &params, scalar_type det_trans) const;
156  virtual void sigma(const base_matrix &E, base_matrix &result,
157  const base_vector &params, scalar_type det_trans) const;
158  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
159  const base_vector &params, scalar_type det_trans) const;
160  explicit Mooney_Rivlin_hyperelastic_law(bool compressible_=false,
161  bool neohookean_=false);
162  };
163 
164  /** Neo-Hookean hyperelastic law variants
165 
166  To be used for compressible problems.
167  For incompressible ones Mooney-Rivlin hyperelastic law can be used.
168  */
170  const bool bonet;
171  virtual scalar_type strain_energy(const base_matrix &E,
172  const base_vector &params, scalar_type det_trans) const;
173  virtual void sigma(const base_matrix &E, base_matrix &result,
174  const base_vector &params, scalar_type det_trans) const;
175  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
176  const base_vector &params, scalar_type det_trans) const;
177  explicit Neo_Hookean_hyperelastic_law(bool bonet_=true);
178  };
179 
180  /** Blatz_Ko hyperelastic law
181 
182  To be used for compressible problems.
183  */
185  virtual scalar_type strain_energy(const base_matrix &E,
186  const base_vector &params, scalar_type det_trans) const;
187  virtual void sigma(const base_matrix &E, base_matrix &result,
188  const base_vector &params, scalar_type det_trans) const;
189  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
190  const base_vector &params, scalar_type det_trans) const;
192  };
193 
194 
195  /** Ciarlet-Geymonat hyperelastic law
196  */
198  // parameters are lambda=params[0], mu=params[1], a=params[2]
199  // The parameter a has to verify a in ]0, mu/2[
200  virtual scalar_type strain_energy(const base_matrix &E,
201  const base_vector &params, scalar_type det_trans) const;
202  virtual void sigma(const base_matrix &E, base_matrix &result,
203  const base_vector &params, scalar_type det_trans) const;
204  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
205  const base_vector &params, scalar_type det_trans) const;
206  Ciarlet_Geymonat_hyperelastic_law() { nb_params_ = 3; }
207  };
208 
209  /** Plane strain hyperelastic law (takes another law as a parameter)
210  */
212  virtual scalar_type strain_energy(const base_matrix &E,
213  const base_vector &params, scalar_type det_trans) const;
214  virtual void sigma(const base_matrix &E, base_matrix &result,
215  const base_vector &params, scalar_type det_trans) const;
216  virtual void grad_sigma(const base_matrix &E, base_tensor &result,
217  const base_vector &params, scalar_type det_trans) const;
218  plane_strain_hyperelastic_law(const phyperelastic_law &pl_)
219  { pl = pl_; nb_params_ = pl->nb_params(); }
220  };
221 
222 
223 
224 
225  template<typename VECT1, typename VECT2> class elasticity_nonlinear_term
226  : public getfem::nonlinear_elem_term {
227  const mesh_fem &mf;
228  std::vector<scalar_type> U;
229  const mesh_fem *mf_data;
230  const VECT2 &PARAMS;
231  size_type N;
232  size_type NFem;
233  const abstract_hyperelastic_law &AHL;
234  base_vector params, coeff;
235  base_matrix E, Sigma, gradU;
236  base_tensor tt;
237  bgeot::multi_index sizes_;
238  int version;
239  public:
240  elasticity_nonlinear_term(const mesh_fem &mf_, const VECT1 &U_,
241  const mesh_fem *mf_data_, const VECT2 &PARAMS_,
242  const abstract_hyperelastic_law &AHL_,
243  int version_)
244  : mf(mf_), U(mf_.nb_basic_dof()), mf_data(mf_data_), PARAMS(PARAMS_),
245  N(mf_.linked_mesh().dim()), NFem(mf_.get_qdim()), AHL(AHL_),
246  params(AHL_.nb_params()), E(N, N), Sigma(N, N), gradU(NFem, N),
247  tt(N, N, N, N), sizes_(NFem, N, NFem, N),
248  version(version_) {
249  switch (version) {
250  case 0 : break; // tangent term
251  case 1 : sizes_.resize(2); break; // rhs
252  case 2 : sizes_.resize(1); sizes_[0] = 1; break; // strain energy
253  case 3 : sizes_.resize(2); break; // Id + grad(u)
254  }
255 
256  mf.extend_vector(U_, U);
257  if (gmm::vect_size(PARAMS) == AHL_.nb_params())
258  gmm::copy(PARAMS, params);
259  }
260 
261  const bgeot::multi_index &sizes(size_type) const { return sizes_; }
262 
263  virtual void compute(getfem::fem_interpolation_context& ctx,
264  bgeot::base_tensor &t) {
265  size_type cv = ctx.convex_num();
266  slice_vector_on_basic_dof_of_element(mf, U, cv, coeff);
267  ctx.pf()->interpolation_grad(ctx, coeff, gradU, mf.get_qdim());
268 
269  for (unsigned int alpha = 0; alpha < N; ++alpha)
270  gradU(alpha, alpha) += scalar_type(1);
271 
272  if (version == 3) {
273  for (size_type n = 0; n < NFem; ++n)
274  for (size_type m = 0; m < N; ++m)
275  t(n,m) = gradU(n,m);
276  return;
277  }
278 
279  gmm::mult(gmm::transposed(gradU), gradU, E);
280  for (unsigned int alpha = 0; alpha < N; ++alpha)
281  E(alpha, alpha) -= scalar_type(1);
282  gmm::scale(E, scalar_type(0.5));
283 
284  scalar_type det_trans = gmm::lu_det(gradU);
285 
286  if (version == 2) {
287  t[0] = AHL.strain_energy(E, params, det_trans);
288  return;
289  }
290 
291  AHL.sigma(E, Sigma, params, det_trans);
292 
293  if (version == 0) {
294  AHL.grad_sigma(E, tt, params, det_trans);
295  for (size_type n = 0; n < NFem; ++n)
296  for (size_type m = 0; m < N; ++m)
297  for (size_type l = 0; l < N; ++l)
298  for (size_type k = 0; k < NFem; ++k) {
299  scalar_type aux = (k == n) ? Sigma(m,l) : 0.0;
300  for (size_type j = 0; j < N; ++j)
301  for (size_type i = 0; i < N; ++i) {
302  aux += gradU(n ,j) * gradU(k, i) * tt(j, m, i, l);
303  }
304  t(n, m, k, l) = aux;
305  }
306  } else {
307  if (det_trans < scalar_type(0)) AHL.inc_unvalid_flag();
308  for (size_type i = 0; i < NFem; ++i)
309  for (size_type j = 0; j < N; ++j) {
310  scalar_type aux(0);
311  for (size_type k = 0; k < N; ++k)
312  aux += gradU(i, k) * Sigma(k, j);
313  t(i,j) = aux;
314  }
315  }
316  }
317 
318  virtual void prepare(fem_interpolation_context& ctx, size_type ) {
319  if (mf_data) {
320  size_type cv = ctx.convex_num();
321  size_type nb = AHL.nb_params();
322  coeff.resize(mf_data->nb_basic_dof_of_element(cv)*nb);
323  for (size_type i = 0; i < mf_data->nb_basic_dof_of_element(cv); ++i)
324  for (size_type k = 0; k < nb; ++k)
325  coeff[i * nb + k]
326  = PARAMS[mf_data->ind_basic_dof_of_element(cv)[i]*nb+k];
327  ctx.pf()->interpolation(ctx, coeff, params, dim_type(nb));
328  }
329  }
330 
331  };
332 
333 
334  /**
335  Tangent matrix for the non-linear elasticity
336  @ingroup asm
337  */
338  template<typename MAT, typename VECT1, typename VECT2>
340  (const MAT &K_, const mesh_im &mim, const getfem::mesh_fem &mf,
341  const VECT1 &U, const getfem::mesh_fem *mf_data, const VECT2 &PARAMS,
342  const abstract_hyperelastic_law &AHL,
343  const mesh_region &rg = mesh_region::all_convexes()) {
344  MAT &K = const_cast<MAT &>(K_);
345  GMM_ASSERT1(mf.get_qdim() >= mf.linked_mesh().dim(),
346  "wrong qdim for the mesh_fem");
347 
348  elasticity_nonlinear_term<VECT1, VECT2>
349  nterm(mf, U, mf_data, PARAMS, AHL, 0);
350  elasticity_nonlinear_term<VECT1, VECT2>
351  nterm2(mf, U, mf_data, PARAMS, AHL, 3);
352 
354  if (mf_data)
355  if (AHL.adapted_tangent_term_assembly_fem_data.size() > 0)
356  assem.set(AHL.adapted_tangent_term_assembly_cte_data);
357  else
358  assem.set("M(#1,#1)+=sym(comp(NonLin$1(#1,#2)(i,j,k,l).vGrad(#1)(:,i,j).vGrad(#1)(:,k,l)))");
359  else
360  if (AHL.adapted_tangent_term_assembly_cte_data.size() > 0)
361  assem.set(AHL.adapted_tangent_term_assembly_cte_data);
362  else
363  assem.set("M(#1,#1)+=sym(comp(NonLin$1(#1)(i,j,k,l).vGrad(#1)(:,i,j).vGrad(#1)(:,k,l)))");
364  assem.push_mi(mim);
365  assem.push_mf(mf);
366  if (mf_data) assem.push_mf(*mf_data);
367  assem.push_data(PARAMS);
368  assem.push_nonlinear_term(&nterm);
369  assem.push_nonlinear_term(&nterm2);
370  assem.push_mat(K);
371  assem.assembly(rg);
372  }
373 
374 
375  /**
376  @ingroup asm
377  */
378  template<typename VECT1, typename VECT2, typename VECT3>
379  void asm_nonlinear_elasticity_rhs
380  (const VECT1 &R_, const mesh_im &mim, const getfem::mesh_fem &mf,
381  const VECT2 &U, const getfem::mesh_fem *mf_data, const VECT3 &PARAMS,
382  const abstract_hyperelastic_law &AHL,
383  const mesh_region &rg = mesh_region::all_convexes()) {
384  VECT1 &R = const_cast<VECT1 &>(R_);
385  GMM_ASSERT1(mf.get_qdim() >= mf.linked_mesh().dim(),
386  "wrong qdim for the mesh_fem");
387 
388  elasticity_nonlinear_term<VECT2, VECT3>
389  nterm(mf, U, mf_data, PARAMS, AHL, 1);
390 
392  if (mf_data)
393  assem.set("t=comp(NonLin(#1,#2).vGrad(#1)); V(#1) += t(i,j,:,i,j)");
394  else
395  assem.set("t=comp(NonLin(#1).vGrad(#1)); V(#1) += t(i,j,:,i,j)");
396  // comp() to be optimized ?
397  assem.push_mi(mim);
398  assem.push_mf(mf);
399  if (mf_data) assem.push_mf(*mf_data);
400  assem.push_nonlinear_term(&nterm);
401  assem.push_vec(R);
402  assem.assembly(rg);
403  }
404 
405  // added by Jean-Yves Heddebaut <jyhed@alpino.be>
406  int levi_civita(int i,int j,int k);
407 
408 
409  /**@ingroup asm
410  */
411  template<typename VECT2, typename VECT3>
412  scalar_type asm_elastic_strain_energy
413  (const mesh_im &mim, const getfem::mesh_fem &mf,
414  const VECT2 &U, const getfem::mesh_fem *mf_data, const VECT3 &PARAMS,
415  const abstract_hyperelastic_law &AHL,
416  const mesh_region &rg = mesh_region::all_convexes()) {
417 
418  GMM_ASSERT1(mf.get_qdim() >= mf.linked_mesh().dim(),
419  "wrong qdim for the mesh_fem");
420 
421  elasticity_nonlinear_term<VECT2, VECT3>
422  nterm(mf, U, mf_data, PARAMS, AHL, 2);
423  std::vector<scalar_type> V(1);
424 
426  if (mf_data)
427  assem.set("V() += comp(NonLin(#1,#2))(i)");
428  else
429  assem.set("V() += comp(NonLin(#1))(i)");
430 
431  assem.push_mi(mim);
432  assem.push_mf(mf);
433  if (mf_data) assem.push_mf(*mf_data);
434  assem.push_nonlinear_term(&nterm);
435  assem.push_vec(V);
436  assem.assembly(rg);
437 
438  return V[0];
439  }
440 
441 
442 
443 
444 
445 
446 
447 
448 
449  /* ******************************************************************** */
450  /* Mixed nonlinear incompressibility assembly procedure */
451  /* ******************************************************************** */
452 
453  template<typename VECT1> class incomp_nonlinear_term
454  : public getfem::nonlinear_elem_term {
455 
456  const mesh_fem &mf;
457  std::vector<scalar_type> U;
458  size_type N;
459  base_vector coeff;
460  base_matrix gradPhi;
461  bgeot::multi_index sizes_;
462  int version;
463 
464  public:
465  incomp_nonlinear_term(const mesh_fem &mf_, const VECT1 &U_,
466  int version_)
467  : mf(mf_), U(mf_.nb_basic_dof()),
468  N(mf_.get_qdim()),
469  gradPhi(N, N), sizes_(N, N),
470  version(version_) {
471  if (version == 1) { sizes_.resize(1); sizes_[0] = 1; }
472  mf.extend_vector(U_, U);
473  }
474 
475  const bgeot::multi_index &sizes(size_type) const { return sizes_; }
476 
477  virtual void compute(getfem::fem_interpolation_context& ctx,
478  bgeot::base_tensor &t) {
479  size_type cv = ctx.convex_num();
480  slice_vector_on_basic_dof_of_element(mf, U, cv, coeff);
481  ctx.pf()->interpolation_grad(ctx, coeff, gradPhi, mf.get_qdim());
482  gmm::add(gmm::identity_matrix(), gradPhi);
483  scalar_type det = gmm::lu_inverse(gradPhi);
484 
485  if (version != 1) {
486  if (version == 2) det = sqrt(gmm::abs(det));
487  for (size_type i = 0; i < N; ++i)
488  for (size_type j = 0; j < N; ++j) {
489  t(i,j) = - det * gradPhi(j,i);
490  }
491  }
492  else t[0] = scalar_type(1) - det;
493 
494  }
495  };
496 
497  /**@ingroup asm*/
498  template<typename MAT1, typename MAT2, typename VECT1, typename VECT2>
499  void asm_nonlinear_incomp_tangent_matrix(const MAT1 &K_, const MAT2 &B_,
500  const mesh_im &mim,
501  const mesh_fem &mf_u,
502  const mesh_fem &mf_p,
503  const VECT1 &U, const VECT2 &P,
504  const mesh_region &rg = mesh_region::all_convexes()) {
505  MAT1 &K = const_cast<MAT1 &>(K_);
506  MAT2 &B = const_cast<MAT2 &>(B_);
507  GMM_ASSERT1(mf_u.get_qdim() == mf_u.linked_mesh().dim(),
508  "wrong qdim for the mesh_fem");
509 
510  incomp_nonlinear_term<VECT1> ntermk(mf_u, U, 0);
511  incomp_nonlinear_term<VECT1> ntermb(mf_u, U, 2);
513  assem("P=data(#2);"
514  "t=comp(NonLin$1(#1).vGrad(#1).Base(#2));"
515  "M$2(#1,#2)+= t(i,j,:,i,j,:);"
516  /*"w=comp(NonLin$2(#1).vGrad(#1).NonLin$2(#1).vGrad(#1).Base(#2));"
517  "M$1(#1,#1)+= w(j,i,:,j,k, m,k,:,m,i,p).P(p)"
518  "-w(i,j,:,i,j, k,l,:,k,l,p).P(p)"*/
519  /*
520  "w=comp(vGrad(#1).NonLin$2(#1).vGrad(#1).NonLin$2(#1).Base(#2));"
521  "M$1(#1,#1)+= w(:,j,k, j,i, :,m,i, m,k, p).P(p)"
522  "-w(:,j,i, j,i, :,m,l, m,l, p).P(p)"
523  */
524  "w1=comp(vGrad(#1)(:,j,k).NonLin$2(#1)(j,i).vGrad(#1)(:,m,i).NonLin$2(#1)(m,k).Base(#2)(p).P(p));"
525  "w2=comp(vGrad(#1)(:,j,i).NonLin$2(#1)(j,i).vGrad(#1)(:,m,l).NonLin$2(#1)(m,l).Base(#2)(p).P(p));"
526  "M$1(#1,#1)+= w1-w2"
527  );
528 
529  assem.push_mi(mim);
530  assem.push_mf(mf_u);
531  assem.push_mf(mf_p);
532  assem.push_nonlinear_term(&ntermk);
533  assem.push_nonlinear_term(&ntermb);
534  assem.push_mat(K);
535  assem.push_mat(B);
536  assem.push_data(P);
537  assem.assembly(rg);
538  }
539 
540 
541  /**@ingroup asm
542  */
543  template<typename VECT1, typename VECT2, typename VECT3>
544  void asm_nonlinear_incomp_rhs
545  (const VECT1 &R_U_, const VECT1 &R_P_, const mesh_im &mim,
546  const getfem::mesh_fem &mf_u, const getfem::mesh_fem &mf_p,
547  const VECT2 &U, const VECT3 &P,
548  const mesh_region &rg = mesh_region::all_convexes()) {
549  VECT1 &R_U = const_cast<VECT1 &>(R_U_);
550  VECT1 &R_P = const_cast<VECT1 &>(R_P_);
551  GMM_ASSERT1(mf_u.get_qdim() == mf_u.linked_mesh().dim(),
552  "wrong qdim for the mesh_fem");
553 
554  incomp_nonlinear_term<VECT2> nterm_tg(mf_u, U, 0);
555  incomp_nonlinear_term<VECT2> nterm(mf_u, U, 1);
556 
558  assem("P=data(#2); "
559  "t=comp(NonLin$1(#1).vGrad(#1).Base(#2));"
560  "V$1(#1) += t(i,j,:,i,j,k).P(k);"
561  "w=comp(NonLin$2(#1).Base(#2)); V$2(#2) += w(1,:)");
562  // assem() to be optimized ?
563 
564  assem.push_mi(mim);
565  assem.push_mf(mf_u);
566  assem.push_mf(mf_p);
567  assem.push_nonlinear_term(&nterm_tg);
568  assem.push_nonlinear_term(&nterm);
569  assem.push_vec(R_U);
570  assem.push_vec(R_P);
571  assem.push_data(P);
572  assem.assembly(rg);
573  }
574 
575 
576 
577  //===========================================================================
578  //
579  // Bricks
580  //
581  //===========================================================================
582 
583 
584  /** Add a nonlinear (large strain) elasticity term to the model with
585  respect to the variable
586  `varname` (deprecated brick, use add_finite_strain_elaticity instead).
587  Note that the constitutive law is described by `AHL` which
588  should not be freed while the model is used. `dataname` described the
589  parameters of the constitutive laws. It could be a vector of value
590  of length the number of parameter of the constitutive law or a vector
591  field described on a finite element method.
592  */
594  (model &md, const mesh_im &mim, const std::string &varname,
595  const phyperelastic_law &AHL, const std::string &dataname,
596  size_type region = size_type(-1));
597 
598 
599 
600  void compute_Von_Mises_or_Tresca
601  (model &md, const std::string &varname, const phyperelastic_law &AHL,
602  const std::string &dataname, const mesh_fem &mf_vm,
603  model_real_plain_vector &VM, bool tresca);
604 
605 
606  void compute_sigmahathat(model &md,
607  const std::string &varname,
608  const phyperelastic_law &AHL,
609  const std::string &dataname,
610  const mesh_fem &mf_sigma,
611  model_real_plain_vector &SIGMA);
612 
613 
614  /**
615  Compute the Von-Mises stress or the Tresca stress of a field
616  with respect to the constitutive elasticity law AHL (only valid in 3D).
617  */
618  template <class VECTVM> void compute_Von_Mises_or_Tresca
619  (model &md, const std::string &varname, const phyperelastic_law &AHL,
620  const std::string &dataname, const mesh_fem &mf_vm,
621  VECTVM &VM, bool tresca) {
622  model_real_plain_vector VMM(mf_vm.nb_dof());
623  compute_Von_Mises_or_Tresca
624  (md, varname, AHL, dataname, mf_vm, VMM, tresca);
625  gmm::copy(VMM, VM);
626  }
627 
628  /** Add a nonlinear incompressibility term (for large strain elasticity)
629  to the model with respect to the variable
630  `varname` (the displacement) and `multname` (the pressure).
631  */
633  (model &md, const mesh_im &mim, const std::string &varname,
634  const std::string &multname, size_type region = size_type(-1));
635 
636  //===========================================================================
637  // High-level generic assembly bricks
638  //===========================================================================
639 
640  /** Add a finite strain elasticity brick
641  to the model with respect to the variable
642  `varname` (the displacement).
643  For 2D meshes, switch automatically to plane strain elasticity.
644  High-level generic assembly version.
645  */
647  (model &md, const mesh_im &mim, const std::string &lawname,
648  const std::string &varname, const std::string &params,
649  size_type region = size_type(-1));
650 
651 
652  /** Add a finite strain incompressibility term (for large strain elasticity)
653  to the model with respect to the variable
654  `varname` (the displacement) and `multname` (the pressure).
655  High-level generic assembly version.
656  */
658  (model &md, const mesh_im &mim, const std::string &varname,
659  const std::string &multname, size_type region = size_type(-1));
660 
661  /**
662  Interpolate the Von-Mises stress of a field `varname`
663  with respect to the nonlinear elasticity constitutive law `lawname`
664  with parameters `params` (only valid in 3D).
665  */
667  (model &md, const std::string &lawname, const std::string &varname,
668  const std::string &params, const mesh_fem &mf_vm,
669  model_real_plain_vector &VM,
670  const mesh_region &rg=mesh_region::all_convexes());
671 
672  IS_DEPRECATED inline void finite_strain_elasticity_Von_Mises
673  (model &md, const std::string &varname, const std::string &lawname,
674  const std::string &params, const mesh_fem &mf_vm,
675  model_real_plain_vector &VM,
676  const mesh_region &rg=mesh_region::all_convexes()) {
677  compute_finite_strain_elasticity_Von_Mises(md, varname, lawname, params,
678  mf_vm, VM, rg);
679  }
680 
681 } /* end of namespace getfem. */
682 
683 
684 #endif /* GETFEM_NONLINEAR_ELASTICITY_H__ */
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector &params, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian
virtual void cauchy_updated_lagrangian(const base_matrix &F, const base_matrix &E, base_matrix &cauchy_stress, const base_vector &params, scalar_type det_trans) const
True Cauchy stress (for Updated Lagrangian formulation)
structure passed as the argument of fem interpolation functions.
Definition: getfem_fem.h:749
Generic assembly of vectors, matrices.
void push_vec(VEC &v)
Add a new output vector.
void push_data(const VEC &d)
Add a new data (dense array)
void push_mat(const MAT &m)
Add a new output matrix (fake const version..)
void assembly(const mesh_region &region=mesh_region::all_convexes())
do the assembly on the specified region (boundary or set of convexes)
void push_nonlinear_term(pnonlinear_elem_term net)
Add a new non-linear term.
void push_mi(const mesh_im &im_)
Add a new mesh_im.
void push_mf(const mesh_fem &mf_)
Add a new mesh_fem.
Describe a finite element method linked to a mesh.
virtual ind_dof_ct ind_basic_dof_of_element(size_type cv) const
Give an array of the dof numbers a of convex.
virtual dim_type get_qdim() const
Return the Q dimension.
virtual size_type nb_dof() const
Return the total number of degrees of freedom.
const mesh & linked_mesh() const
Return a reference to the underlying mesh.
virtual size_type nb_basic_dof_of_element(size_type cv) const
Return the number of degrees of freedom attached to a given convex.
Describe an integration method linked to a mesh.
structure used to hold a set of convexes and/or convex faces.
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
`‘Model’' variables store the variables, the data and the description of a model.
abstract class for integration of non-linear terms into the mat_elem computations the nonlinear term ...
Generic assembly implementation.
Compute the gradient of a field on a getfem::mesh_fem.
A language for generic assembly of pde boundary value problems.
Interpolation of fields from a mesh_fem onto another.
Model representation in Getfem.
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*‍/
Definition: gmm_blas.h:1663
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
linalg_traits< DenseMatrixLU >::value_type lu_det(const DenseMatrixLU &LU, const Pvector &pvector)
Compute the matrix determinant (via a LU factorization)
Definition: gmm_dense_lu.h:240
void lu_inverse(const DenseMatrixLU &LU, const Pvector &pvector, const DenseMatrix &AInv_)
Given an LU factored matrix, build the inverse of the matrix.
Definition: gmm_dense_lu.h:211
Input/output on sparse matrices.
void asm_nonlinear_elasticity_tangent_matrix(const MAT &K_, const mesh_im &mim, const getfem::mesh_fem &mf, const VECT1 &U, const getfem::mesh_fem *mf_data, const VECT2 &PARAMS, const abstract_hyperelastic_law &AHL, const mesh_region &rg=mesh_region::all_convexes())
Tangent matrix for the non-linear elasticity.
size_type convex_num() const
get the current convex number
Definition: getfem_fem.cc:52
const pfem pf() const
get the current FEM descriptor
Definition: getfem_fem.h:781
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
size_type alpha(short_type n, short_type d)
Return the value of which is the number of monomials of a polynomial of variables and degree .
Definition: bgeot_poly.cc:46
GEneric Tool for Finite Element Methods.
size_type add_finite_strain_elasticity_brick(model &md, const mesh_im &mim, const std::string &lawname, const std::string &varname, const std::string &params, size_type region=size_type(-1))
Add a finite strain elasticity brick to the model with respect to the variable varname (the displacem...
size_type add_finite_strain_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a finite strain incompressibility term (for large strain elasticity) to the model with respect to...
size_type add_nonlinear_elasticity_brick(model &md, const mesh_im &mim, const std::string &varname, const phyperelastic_law &AHL, const std::string &dataname, size_type region=size_type(-1))
Add a nonlinear (large strain) elasticity term to the model with respect to the variable varname (dep...
size_type add_nonlinear_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a nonlinear incompressibility term (for large strain elasticity) to the model with respect to the...
void compute_finite_strain_elasticity_Von_Mises(model &md, const std::string &lawname, const std::string &varname, const std::string &params, const mesh_fem &mf_vm, model_real_plain_vector &VM, const mesh_region &rg=mesh_region::all_convexes())
Interpolate the Von-Mises stress of a field varname with respect to the nonlinear elasticity constitu...
void slice_vector_on_basic_dof_of_element(const mesh_fem &mf, const VEC1 &vec, size_type cv, VEC2 &coeff, size_type qmult1=size_type(-1), size_type qmult2=size_type(-1))
Given a mesh_fem.
Neo-Hookean hyperelastic law variants.
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
Linear law for a membrane (plane stress), orthotropic material caracterized by Ex,...
Plane strain hyperelastic law (takes another law as a parameter)