38 #ifndef GETFEM_NONLINEAR_ELASTICITY_H__ 39 #define GETFEM_NONLINEAR_ELASTICITY_H__ 51 int check_symmetry(
const base_tensor &t);
53 class abstract_hyperelastic_law;
54 typedef std::shared_ptr<const abstract_hyperelastic_law> phyperelastic_law;
64 void reset_unvalid_flag()
const { uvflag = 0; }
65 void inc_unvalid_flag()
const { uvflag++; }
66 int get_unvalid_flag()
const {
return uvflag; }
67 std::string adapted_tangent_term_assembly_fem_data;
68 std::string adapted_tangent_term_assembly_cte_data;
71 virtual scalar_type strain_energy(
const base_matrix &E,
72 const base_vector ¶ms,
73 scalar_type det_trans)
const = 0;
77 base_matrix &cauchy_stress,
78 const base_vector ¶ms,
79 scalar_type det_trans)
const;
80 virtual void sigma(
const base_matrix &E, base_matrix &result,
81 const base_vector ¶ms,
82 scalar_type det_trans)
const = 0;
84 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
85 const base_vector ¶ms, scalar_type det_trans)
const = 0;
90 const base_vector ¶ms,
91 scalar_type det_trans,
92 base_tensor &grad_sigma_ul)
const;
94 size_type nb_params()
const {
return nb_params_; }
97 static void random_E(base_matrix &E);
98 void test_derivatives(
size_type N, scalar_type h,
99 const base_vector& param)
const;
110 virtual scalar_type strain_energy(
const base_matrix &E,
111 const base_vector ¶ms,
112 scalar_type det_trans)
const;
114 virtual void sigma(
const base_matrix &E, base_matrix &result,
115 const base_vector ¶ms, scalar_type det_trans)
const;
116 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
117 const base_vector ¶ms, scalar_type det_trans)
const;
119 const base_matrix& E,
120 const base_vector ¶ms,
121 scalar_type det_trans,
122 base_tensor &grad_sigma_ul)
const;
133 virtual scalar_type strain_energy(
const base_matrix &E,
134 const base_vector ¶ms,
135 scalar_type det_trans)
const;
136 virtual void sigma(
const base_matrix &E, base_matrix &result,
137 const base_vector ¶ms, scalar_type det_trans)
const;
138 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
139 const base_vector ¶ms, scalar_type det_trans)
const;
154 const bool compressible, neohookean;
155 virtual scalar_type strain_energy(
const base_matrix &E,
156 const base_vector ¶ms, scalar_type det_trans)
const;
157 virtual void sigma(
const base_matrix &E, base_matrix &result,
158 const base_vector ¶ms, scalar_type det_trans)
const;
159 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
160 const base_vector ¶ms, scalar_type det_trans)
const;
162 bool neohookean_=
false);
172 virtual scalar_type strain_energy(
const base_matrix &E,
173 const base_vector ¶ms, scalar_type det_trans)
const;
174 virtual void sigma(
const base_matrix &E, base_matrix &result,
175 const base_vector ¶ms, scalar_type det_trans)
const;
176 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
177 const base_vector ¶ms, scalar_type det_trans)
const;
186 virtual scalar_type strain_energy(
const base_matrix &E,
187 const base_vector ¶ms, scalar_type det_trans)
const;
188 virtual void sigma(
const base_matrix &E, base_matrix &result,
189 const base_vector ¶ms, scalar_type det_trans)
const;
190 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
191 const base_vector ¶ms, scalar_type det_trans)
const;
201 virtual scalar_type strain_energy(
const base_matrix &E,
202 const base_vector ¶ms, scalar_type det_trans)
const;
203 virtual void sigma(
const base_matrix &E, base_matrix &result,
204 const base_vector ¶ms, scalar_type det_trans)
const;
205 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
206 const base_vector ¶ms, scalar_type det_trans)
const;
213 virtual scalar_type strain_energy(
const base_matrix &E,
214 const base_vector ¶ms, scalar_type det_trans)
const;
215 virtual void sigma(
const base_matrix &E, base_matrix &result,
216 const base_vector ¶ms, scalar_type det_trans)
const;
217 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
218 const base_vector ¶ms, scalar_type det_trans)
const;
220 { pl = pl_; nb_params_ = pl->nb_params(); }
226 template<
typename VECT1,
typename VECT2>
class elasticity_nonlinear_term
229 std::vector<scalar_type> U;
235 base_vector params, coeff;
236 base_matrix E, Sigma, gradU;
238 bgeot::multi_index sizes_;
241 elasticity_nonlinear_term(
const mesh_fem &mf_,
const VECT1 &U_,
242 const mesh_fem *mf_data_,
const VECT2 &PARAMS_,
245 : mf(mf_), U(mf_.
nb_basic_dof()), mf_data(mf_data_), PARAMS(PARAMS_),
247 params(AHL_.nb_params()), E(N, N), Sigma(N, N), gradU(NFem, N),
248 tt(N, N, N, N), sizes_(NFem, N, NFem, N),
252 case 1 : sizes_.resize(2);
break;
253 case 2 : sizes_.resize(1); sizes_[0] = 1;
break;
254 case 3 : sizes_.resize(2);
break;
257 mf.extend_vector(U_, U);
258 if (gmm::vect_size(PARAMS) == AHL_.nb_params())
259 gmm::copy(PARAMS, params);
262 const bgeot::multi_index &sizes(
size_type)
const {
return sizes_; }
265 bgeot::base_tensor &t) {
268 ctx.
pf()->interpolation_grad(ctx, coeff, gradU, mf.
get_qdim());
270 for (
unsigned int alpha = 0; alpha < N; ++alpha)
271 gradU(alpha, alpha) += scalar_type(1);
280 gmm::mult(gmm::transposed(gradU), gradU, E);
281 for (
unsigned int alpha = 0; alpha < N; ++alpha)
282 E(alpha, alpha) -= scalar_type(1);
283 gmm::scale(E, scalar_type(0.5));
285 scalar_type det_trans = gmm::lu_det(gradU);
288 t[0] = AHL.strain_energy(E, params, det_trans);
292 AHL.sigma(E, Sigma, params, det_trans);
295 AHL.grad_sigma(E, tt, params, det_trans);
300 scalar_type aux = (k == n) ? Sigma(m,l) : 0.0;
303 aux += gradU(n ,j) * gradU(k, i) * tt(j, m, i, l);
308 if (det_trans < scalar_type(0)) AHL.inc_unvalid_flag();
313 aux += gradU(i, k) * Sigma(k, j);
328 ctx.
pf()->interpolation(ctx, coeff, params, dim_type(nb));
339 template<
typename MAT,
typename VECT1,
typename VECT2>
345 MAT &K =
const_cast<MAT &
>(K_);
347 "wrong qdim for the mesh_fem");
349 elasticity_nonlinear_term<VECT1, VECT2>
350 nterm(mf, U, mf_data, PARAMS, AHL, 0);
351 elasticity_nonlinear_term<VECT1, VECT2>
352 nterm2(mf, U, mf_data, PARAMS, AHL, 3);
356 if (AHL.adapted_tangent_term_assembly_fem_data.size() > 0)
357 assem.set(AHL.adapted_tangent_term_assembly_cte_data);
359 assem.set(
"M(#1,#1)+=sym(comp(NonLin$1(#1,#2)(i,j,k,l).vGrad(#1)(:,i,j).vGrad(#1)(:,k,l)))");
361 if (AHL.adapted_tangent_term_assembly_cte_data.size() > 0)
362 assem.set(AHL.adapted_tangent_term_assembly_cte_data);
364 assem.set(
"M(#1,#1)+=sym(comp(NonLin$1(#1)(i,j,k,l).vGrad(#1)(:,i,j).vGrad(#1)(:,k,l)))");
367 if (mf_data) assem.
push_mf(*mf_data);
379 template<
typename VECT1,
typename VECT2,
typename VECT3>
380 void asm_nonlinear_elasticity_rhs
385 VECT1 &R =
const_cast<VECT1 &
>(R_);
387 "wrong qdim for the mesh_fem");
389 elasticity_nonlinear_term<VECT2, VECT3>
390 nterm(mf, U, mf_data, PARAMS, AHL, 1);
394 assem.set(
"t=comp(NonLin(#1,#2).vGrad(#1)); V(#1) += t(i,j,:,i,j)");
396 assem.set(
"t=comp(NonLin(#1).vGrad(#1)); V(#1) += t(i,j,:,i,j)");
400 if (mf_data) assem.
push_mf(*mf_data);
407 int levi_civita(
int i,
int j,
int k);
412 template<
typename VECT2,
typename VECT3>
413 scalar_type asm_elastic_strain_energy
420 "wrong qdim for the mesh_fem");
422 elasticity_nonlinear_term<VECT2, VECT3>
423 nterm(mf, U, mf_data, PARAMS, AHL, 2);
424 std::vector<scalar_type> V(1);
428 assem.set(
"V() += comp(NonLin(#1,#2))(i)");
430 assem.set(
"V() += comp(NonLin(#1))(i)");
434 if (mf_data) assem.
push_mf(*mf_data);
454 template<
typename VECT1>
class incomp_nonlinear_term
458 std::vector<scalar_type> U;
462 bgeot::multi_index sizes_;
466 incomp_nonlinear_term(
const mesh_fem &mf_,
const VECT1 &U_,
470 gradPhi(N, N), sizes_(N, N),
472 if (version == 1) { sizes_.resize(1); sizes_[0] = 1; }
473 mf.extend_vector(U_, U);
476 const bgeot::multi_index &sizes(
size_type)
const {
return sizes_; }
479 bgeot::base_tensor &t) {
482 ctx.
pf()->interpolation_grad(ctx, coeff, gradPhi, mf.
get_qdim());
483 gmm::add(gmm::identity_matrix(), gradPhi);
484 scalar_type det = gmm::lu_inverse(gradPhi);
487 if (version == 2) det = sqrt(gmm::abs(det));
490 t(i,j) = - det * gradPhi(j,i);
493 else t[0] = scalar_type(1) - det;
499 template<
typename MAT1,
typename MAT2,
typename VECT1,
typename VECT2>
500 void asm_nonlinear_incomp_tangent_matrix(
const MAT1 &K_,
const MAT2 &B_,
504 const VECT1 &U,
const VECT2 &P,
506 MAT1 &K =
const_cast<MAT1 &
>(K_);
507 MAT2 &B =
const_cast<MAT2 &
>(B_);
509 "wrong qdim for the mesh_fem");
511 incomp_nonlinear_term<VECT1> ntermk(mf_u, U, 0);
512 incomp_nonlinear_term<VECT1> ntermb(mf_u, U, 2);
515 "t=comp(NonLin$1(#1).vGrad(#1).Base(#2));" 516 "M$2(#1,#2)+= t(i,j,:,i,j,:);" 525 "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));" 526 "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));" 544 template<
typename VECT1,
typename VECT2,
typename VECT3>
545 void asm_nonlinear_incomp_rhs
546 (
const VECT1 &R_U_,
const VECT1 &R_P_,
const mesh_im &mim,
548 const VECT2 &U,
const VECT3 &P,
550 VECT1 &R_U =
const_cast<VECT1 &
>(R_U_);
551 VECT1 &R_P =
const_cast<VECT1 &
>(R_P_);
553 "wrong qdim for the mesh_fem");
555 incomp_nonlinear_term<VECT2> nterm_tg(mf_u, U, 0);
556 incomp_nonlinear_term<VECT2> nterm(mf_u, U, 1);
560 "t=comp(NonLin$1(#1).vGrad(#1).Base(#2));" 561 "V$1(#1) += t(i,j,:,i,j,k).P(k);" 562 "w=comp(NonLin$2(#1).Base(#2)); V$2(#2) += w(1,:)");
595 (
model &md,
const mesh_im &mim,
const std::string &varname,
596 const phyperelastic_law &AHL,
const std::string &dataname,
601 void compute_Von_Mises_or_Tresca
602 (
model &md,
const std::string &varname,
const phyperelastic_law &AHL,
603 const std::string &dataname,
const mesh_fem &mf_vm,
604 model_real_plain_vector &VM,
bool tresca);
607 void compute_sigmahathat(
model &md,
608 const std::string &varname,
609 const phyperelastic_law &AHL,
610 const std::string &dataname,
612 model_real_plain_vector &SIGMA);
619 template <
class VECTVM>
void compute_Von_Mises_or_Tresca
620 (
model &md,
const std::string &varname,
const phyperelastic_law &AHL,
621 const std::string &dataname,
const mesh_fem &mf_vm,
622 VECTVM &VM,
bool tresca) {
623 model_real_plain_vector VMM(mf_vm.
nb_dof());
624 compute_Von_Mises_or_Tresca
625 (md, varname, AHL, dataname, mf_vm, VMM, tresca);
634 (
model &md,
const mesh_im &mim,
const std::string &varname,
648 (
model &md,
const mesh_im &mim,
const std::string &lawname,
649 const std::string &varname,
const std::string ¶ms,
659 (
model &md,
const mesh_im &mim,
const std::string &varname,
668 (
model &md,
const std::string &lawname,
const std::string &varname,
669 const std::string ¶ms,
const mesh_fem &mf_vm,
670 model_real_plain_vector &VM,
673 IS_DEPRECATED
inline void finite_strain_elasticity_Von_Mises
674 (
model &md,
const std::string &varname,
const std::string &lawname,
675 const std::string ¶ms,
const mesh_fem &mf_vm,
676 model_real_plain_vector &VM,
virtual size_type nb_basic_dof_of_element(size_type cv) const
Return the number of degrees of freedom attached to a given convex.
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.
Mooney-Rivlin hyperelastic law.
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.
void push_data(const VEC &d)
Add a new data (dense array)
static mesh_region all_convexes()
provide a default value for the mesh_region parameters of assembly procedures etc.
size_type convex_num() const
get the current convex number
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 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.
virtual void cauchy_updated_lagrangian(const base_matrix &F, const base_matrix &E, base_matrix &cauchy_stress, const base_vector ¶ms, scalar_type det_trans) const
True Cauchy stress (for Updated Lagrangian formulation)
Describe an integration method linked to a mesh.
Base class for material law.
void assembly(const mesh_region ®ion=mesh_region::all_convexes())
do the assembly on the specified region (boundary or set of convexes)
size_type add_finite_strain_elasticity_brick(model &md, const mesh_im &mim, const std::string &lawname, const std::string &varname, const std::string ¶ms, size_type region=size_type(-1))
Add a finite strain elasticity brick to the model with respect to the variable varname (the displacem...
void push_mf(const mesh_fem &mf_)
Add a new mesh_fem.
void compute_finite_strain_elasticity_Von_Mises(model &md, const std::string &lawname, const std::string &varname, const std::string ¶ms, const mesh_fem &mf_vm, model_real_plain_vector &VM, const mesh_region &rg=mesh_region::all_convexes())
Interpolate the Von-Mises stress of a field varname with respect to the nonlinear elasticity constitu...
Generic assembly of vectors, matrices.
``Model'' variables store the variables, the data and the description of a model. ...
size_t size_type
used as the common size type in the library
virtual size_type nb_basic_dof() const
Return the total number of basic degrees of freedom (before the optional reduction).
Model representation in Getfem.
Ciarlet-Geymonat hyperelastic law.
Input/output on sparse matrices.
Interpolation of fields from a mesh_fem onto another.
virtual dim_type get_qdim() const
Return the Q dimension.
structure passed as the argument of fem interpolation functions.
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.
Neo-Hookean hyperelastic law variants.
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 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 pfem pf() const
get the current FEM descriptor
Generic assembly implementation.
Linear law for a membrane (plane stress), orthotropic material caracterized by Ex, Ey, vYX and G, with optional orthotropic prestresses.
void push_mat(const MAT &m)
Add a new output matrix (fake const version..)
abstract class for integration of non-linear terms into the mat_elem computations the nonlinear term ...
virtual ind_dof_ct ind_basic_dof_of_element(size_type cv) const
Give an array of the dof numbers a of convex.
Blatz_Ko hyperelastic law.
Compute the gradient of a field on a getfem::mesh_fem.
void push_vec(VEC &v)
Add a new output vector.
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector ¶ms, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian
Plane strain hyperelastic law (takes another law as a parameter)
Saint-Venant Kirchhoff hyperelastic law.
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...