37 #ifndef GETFEM_NONLINEAR_ELASTICITY_H__
38 #define GETFEM_NONLINEAR_ELASTICITY_H__
50 int check_symmetry(
const base_tensor &t);
52 class abstract_hyperelastic_law;
53 typedef std::shared_ptr<const abstract_hyperelastic_law> phyperelastic_law;
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;
67 std::string adapted_tangent_term_assembly_cte_data;
70 virtual scalar_type strain_energy(
const base_matrix &E,
71 const base_vector ¶ms,
72 scalar_type det_trans)
const = 0;
76 base_matrix &cauchy_stress,
77 const base_vector ¶ms,
78 scalar_type det_trans)
const;
79 virtual void sigma(
const base_matrix &E, base_matrix &result,
80 const base_vector ¶ms,
81 scalar_type det_trans)
const = 0;
83 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
84 const base_vector ¶ms, scalar_type det_trans)
const = 0;
89 const base_vector ¶ms,
90 scalar_type det_trans,
91 base_tensor &grad_sigma_ul)
const;
93 size_type nb_params()
const {
return nb_params_; }
96 static void random_E(base_matrix &E);
97 void test_derivatives(
size_type N, scalar_type h,
98 const base_vector& param)
const;
109 virtual scalar_type strain_energy(
const base_matrix &E,
110 const base_vector ¶ms,
111 scalar_type det_trans)
const;
113 virtual void sigma(
const base_matrix &E, base_matrix &result,
114 const base_vector ¶ms, scalar_type det_trans)
const;
115 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
116 const base_vector ¶ms, scalar_type det_trans)
const;
118 const base_matrix& E,
119 const base_vector ¶ms,
120 scalar_type det_trans,
121 base_tensor &grad_sigma_ul)
const;
132 virtual scalar_type strain_energy(
const base_matrix &E,
133 const base_vector ¶ms,
134 scalar_type det_trans)
const;
135 virtual void sigma(
const base_matrix &E, base_matrix &result,
136 const base_vector ¶ms, scalar_type det_trans)
const;
137 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
138 const base_vector ¶ms, scalar_type det_trans)
const;
153 const bool compressible, neohookean;
154 virtual scalar_type strain_energy(
const base_matrix &E,
155 const base_vector ¶ms, scalar_type det_trans)
const;
156 virtual void sigma(
const base_matrix &E, base_matrix &result,
157 const base_vector ¶ms, scalar_type det_trans)
const;
158 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
159 const base_vector ¶ms, scalar_type det_trans)
const;
161 bool neohookean_=
false);
171 virtual scalar_type strain_energy(
const base_matrix &E,
172 const base_vector ¶ms, scalar_type det_trans)
const;
173 virtual void sigma(
const base_matrix &E, base_matrix &result,
174 const base_vector ¶ms, scalar_type det_trans)
const;
175 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
176 const base_vector ¶ms, scalar_type det_trans)
const;
185 virtual scalar_type strain_energy(
const base_matrix &E,
186 const base_vector ¶ms, scalar_type det_trans)
const;
187 virtual void sigma(
const base_matrix &E, base_matrix &result,
188 const base_vector ¶ms, scalar_type det_trans)
const;
189 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
190 const base_vector ¶ms, scalar_type det_trans)
const;
200 virtual scalar_type strain_energy(
const base_matrix &E,
201 const base_vector ¶ms, scalar_type det_trans)
const;
202 virtual void sigma(
const base_matrix &E, base_matrix &result,
203 const base_vector ¶ms, scalar_type det_trans)
const;
204 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
205 const base_vector ¶ms, scalar_type det_trans)
const;
212 virtual scalar_type strain_energy(
const base_matrix &E,
213 const base_vector ¶ms, scalar_type det_trans)
const;
214 virtual void sigma(
const base_matrix &E, base_matrix &result,
215 const base_vector ¶ms, scalar_type det_trans)
const;
216 virtual void grad_sigma(
const base_matrix &E, base_tensor &result,
217 const base_vector ¶ms, scalar_type det_trans)
const;
219 { pl = pl_; nb_params_ = pl->nb_params(); }
225 template<
typename VECT1,
typename VECT2>
class elasticity_nonlinear_term
228 std::vector<scalar_type> U;
234 base_vector params, coeff;
235 base_matrix E, Sigma, gradU;
237 bgeot::multi_index sizes_;
240 elasticity_nonlinear_term(
const mesh_fem &mf_,
const VECT1 &U_,
241 const mesh_fem *mf_data_,
const VECT2 &PARAMS_,
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),
251 case 1 : sizes_.resize(2);
break;
252 case 2 : sizes_.resize(1); sizes_[0] = 1;
break;
253 case 3 : sizes_.resize(2);
break;
256 mf.extend_vector(U_, U);
257 if (gmm::vect_size(PARAMS) == AHL_.nb_params())
258 gmm::copy(PARAMS, params);
261 const bgeot::multi_index &sizes(
size_type)
const {
return sizes_; }
264 bgeot::base_tensor &t) {
267 ctx.
pf()->interpolation_grad(ctx, coeff, gradU, mf.
get_qdim());
269 for (
unsigned int alpha = 0;
alpha < N; ++
alpha)
270 gradU(alpha, alpha) += scalar_type(1);
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));
287 t[0] = AHL.strain_energy(E, params, det_trans);
291 AHL.sigma(E, Sigma, params, det_trans);
294 AHL.grad_sigma(E, tt, params, det_trans);
299 scalar_type aux = (k == n) ? Sigma(m,l) : 0.0;
302 aux += gradU(n ,j) * gradU(k, i) * tt(j, m, i, l);
307 if (det_trans < scalar_type(0)) AHL.inc_unvalid_flag();
312 aux += gradU(i, k) * Sigma(k, j);
318 virtual void prepare(fem_interpolation_context& ctx,
size_type ) {
327 ctx.pf()->interpolation(ctx, coeff, params, dim_type(nb));
338 template<
typename MAT,
typename VECT1,
typename VECT2>
344 MAT &K =
const_cast<MAT &
>(K_);
346 "wrong qdim for the mesh_fem");
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);
355 if (AHL.adapted_tangent_term_assembly_fem_data.size() > 0)
356 assem.set(AHL.adapted_tangent_term_assembly_cte_data);
358 assem.set(
"M(#1,#1)+=sym(comp(NonLin$1(#1,#2)(i,j,k,l).vGrad(#1)(:,i,j).vGrad(#1)(:,k,l)))");
360 if (AHL.adapted_tangent_term_assembly_cte_data.size() > 0)
361 assem.set(AHL.adapted_tangent_term_assembly_cte_data);
363 assem.set(
"M(#1,#1)+=sym(comp(NonLin$1(#1)(i,j,k,l).vGrad(#1)(:,i,j).vGrad(#1)(:,k,l)))");
366 if (mf_data) assem.
push_mf(*mf_data);
378 template<
typename VECT1,
typename VECT2,
typename VECT3>
379 void asm_nonlinear_elasticity_rhs
382 const abstract_hyperelastic_law &AHL,
384 VECT1 &R =
const_cast<VECT1 &
>(R_);
386 "wrong qdim for the mesh_fem");
388 elasticity_nonlinear_term<VECT2, VECT3>
389 nterm(mf, U, mf_data, PARAMS, AHL, 1);
393 assem.set(
"t=comp(NonLin(#1,#2).vGrad(#1)); V(#1) += t(i,j,:,i,j)");
395 assem.set(
"t=comp(NonLin(#1).vGrad(#1)); V(#1) += t(i,j,:,i,j)");
399 if (mf_data) assem.
push_mf(*mf_data);
406 int levi_civita(
int i,
int j,
int k);
411 template<
typename VECT2,
typename VECT3>
412 scalar_type asm_elastic_strain_energy
415 const abstract_hyperelastic_law &AHL,
419 "wrong qdim for the mesh_fem");
421 elasticity_nonlinear_term<VECT2, VECT3>
422 nterm(mf, U, mf_data, PARAMS, AHL, 2);
423 std::vector<scalar_type> V(1);
427 assem.set(
"V() += comp(NonLin(#1,#2))(i)");
429 assem.set(
"V() += comp(NonLin(#1))(i)");
433 if (mf_data) assem.
push_mf(*mf_data);
453 template<
typename VECT1>
class incomp_nonlinear_term
457 std::vector<scalar_type> U;
461 bgeot::multi_index sizes_;
465 incomp_nonlinear_term(
const mesh_fem &mf_,
const VECT1 &U_,
467 : mf(mf_), U(mf_.nb_basic_dof()),
469 gradPhi(N, N), sizes_(N, N),
471 if (version == 1) { sizes_.resize(1); sizes_[0] = 1; }
472 mf.extend_vector(U_, U);
475 const bgeot::multi_index &sizes(
size_type)
const {
return sizes_; }
478 bgeot::base_tensor &t) {
481 ctx.
pf()->interpolation_grad(ctx, coeff, gradPhi, mf.get_qdim());
482 gmm::add(gmm::identity_matrix(), gradPhi);
486 if (version == 2) det = sqrt(gmm::abs(det));
489 t(i,j) = - det * gradPhi(j,i);
492 else t[0] = scalar_type(1) - det;
498 template<
typename MAT1,
typename MAT2,
typename VECT1,
typename VECT2>
499 void asm_nonlinear_incomp_tangent_matrix(
const MAT1 &K_,
const MAT2 &B_,
501 const mesh_fem &mf_u,
502 const mesh_fem &mf_p,
503 const VECT1 &U,
const VECT2 &P,
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");
510 incomp_nonlinear_term<VECT1> ntermk(mf_u, U, 0);
511 incomp_nonlinear_term<VECT1> ntermb(mf_u, U, 2);
514 "t=comp(NonLin$1(#1).vGrad(#1).Base(#2));"
515 "M$2(#1,#2)+= t(i,j,:,i,j,:);"
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));"
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,
547 const VECT2 &U,
const VECT3 &P,
549 VECT1 &R_U =
const_cast<VECT1 &
>(R_U_);
550 VECT1 &R_P =
const_cast<VECT1 &
>(R_P_);
552 "wrong qdim for the mesh_fem");
554 incomp_nonlinear_term<VECT2> nterm_tg(mf_u, U, 0);
555 incomp_nonlinear_term<VECT2> nterm(mf_u, U, 1);
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,:)");
594 (model &md,
const mesh_im &mim,
const std::string &varname,
595 const phyperelastic_law &AHL,
const std::string &dataname,
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);
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);
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);
633 (model &md,
const mesh_im &mim,
const std::string &varname,
647 (model &md,
const mesh_im &mim,
const std::string &lawname,
648 const std::string &varname,
const std::string ¶ms,
658 (model &md,
const mesh_im &mim,
const std::string &varname,
667 (model &md,
const std::string &lawname,
const std::string &varname,
668 const std::string ¶ms,
const mesh_fem &mf_vm,
669 model_real_plain_vector &VM,
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 ¶ms,
const mesh_fem &mf_vm,
675 model_real_plain_vector &VM,
Base class for material law.
virtual void grad_sigma_updated_lagrangian(const base_matrix &F, const base_matrix &E, const base_vector ¶ms, scalar_type det_trans, base_tensor &grad_sigma_ul) const
cauchy-truesdel tangent moduli, used in updated lagrangian
virtual void cauchy_updated_lagrangian(const base_matrix &F, const base_matrix &E, base_matrix &cauchy_stress, const base_vector ¶ms, scalar_type det_trans) const
True Cauchy stress (for Updated Lagrangian formulation)
structure passed as the argument of fem interpolation functions.
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 ®ion=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)
*/
void add(const L1 &l1, L2 &l2)
*/
linalg_traits< DenseMatrixLU >::value_type lu_det(const DenseMatrixLU &LU, const Pvector &pvector)
Compute the matrix determinant (via a LU factorization)
void lu_inverse(const DenseMatrixLU &LU, const Pvector &pvector, const DenseMatrix &AInv_)
Given an LU factored matrix, build the inverse of the matrix.
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
const pfem pf() const
get the current FEM descriptor
size_t size_type
used as the common size type in the library
size_type alpha(short_type n, short_type d)
Return the value of which is the number of monomials of a polynomial of variables and degree .
GEneric Tool for Finite Element Methods.
size_type add_finite_strain_elasticity_brick(model &md, const mesh_im &mim, const std::string &lawname, const std::string &varname, const std::string ¶ms, size_type region=size_type(-1))
Add a finite strain elasticity brick to the model with respect to the variable varname (the displacem...
size_type add_finite_strain_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a finite strain incompressibility term (for large strain elasticity) to the model with respect to...
size_type add_nonlinear_elasticity_brick(model &md, const mesh_im &mim, const std::string &varname, const phyperelastic_law &AHL, const std::string &dataname, size_type region=size_type(-1))
Add a nonlinear (large strain) elasticity term to the model with respect to the variable varname (dep...
size_type add_nonlinear_incompressibility_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region=size_type(-1))
Add a nonlinear incompressibility term (for large strain elasticity) to the model with respect to the...
void compute_finite_strain_elasticity_Von_Mises(model &md, const std::string &lawname, const std::string &varname, const std::string ¶ms, const mesh_fem &mf_vm, model_real_plain_vector &VM, const mesh_region &rg=mesh_region::all_convexes())
Interpolate the Von-Mises stress of a field varname with respect to the nonlinear elasticity constitu...
void 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.
Ciarlet-Geymonat hyperelastic law.
Mooney-Rivlin hyperelastic law.
Neo-Hookean hyperelastic law variants.
Saint-Venant Kirchhoff hyperelastic law.
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
Blatz_Ko hyperelastic law.
Linear law for a membrane (plane stress), orthotropic material caracterized by Ex,...
Plane strain hyperelastic law (takes another law as a parameter)