40 #ifndef _NEOHOOKEANMATERIAL_H_ 41 #define _NEOHOOKEANMATERIAL_H_ 44 #include <lifev/structure/solver/isotropic/StructuralIsotropicConstitutiveLaw.hpp> 49 template <
typename MeshType>
50 class NeoHookeanMaterialNonLinear :
51 public StructuralIsotropicConstitutiveLaw<MeshType>
57 typedef StructuralIsotropicConstitutiveLaw<MeshType> super;
61 typedef typename super::vector_Type vector_Type;
62 typedef typename super::matrix_Type matrix_Type;
64 typedef typename super::matrixPtr_Type matrixPtr_Type;
65 typedef typename super::vectorPtr_Type vectorPtr_Type;
66 typedef typename super::dataPtr_Type dataPtr_Type;
67 typedef typename super::displayerPtr_Type displayerPtr_Type;
69 typedef typename super::mapMarkerVolumesPtr_Type mapMarkerVolumesPtr_Type;
70 typedef typename super::mapMarkerVolumes_Type mapMarkerVolumes_Type;
71 typedef typename mapMarkerVolumes_Type::const_iterator mapIterator_Type;
73 typedef typename super::vectorVolumes_Type vectorVolumes_Type;
74 typedef std::shared_ptr<vectorVolumes_Type> vectorVolumesPtr_Type;
76 typedef std::vector<UInt> vectorIndexes_Type;
77 typedef std::shared_ptr<vectorIndexes_Type> vectorIndexesPtr_Type;
79 typedef typename super::mapMarkerIndexesPtr_Type mapMarkerIndexesPtr_Type;
80 typedef typename super::mapMarkerIndexes_Type mapMarkerIndexes_Type;
81 typedef typename mapMarkerIndexes_Type::const_iterator mapIteratorIndex_Type;
83 typedef typename super::FESpacePtr_Type FESpacePtr_Type;
84 typedef typename super::ETFESpacePtr_Type ETFESpacePtr_Type;
87 typedef typename super::vectorsParameters_Type vectorsParameters_Type;
88 typedef typename super::vectorsParametersPtr_Type vectorsParametersPtr_Type;
90 typedef MatrixSmall<3, 3> matrixSmall_Type;
93 typedef typename super::tensorF_Type tensorF_Type;
94 typedef typename super::determinantF_Type determinantF_Type;
95 typedef typename super::tensorC_Type tensorC_Type;
96 typedef typename super::minusT_Type minusT_Type;
97 typedef typename super::traceTensor_Type traceTensor_Type;
105 NeoHookeanMaterialNonLinear();
107 virtual ~NeoHookeanMaterialNonLinear();
120 void setup (
const FESpacePtr_Type& dFESpace,
121 const ETFESpacePtr_Type& dETFESpace,
122 const std::shared_ptr<
const MapEpetra>& monolithicMap,
123 const UInt offset,
const dataPtr_Type& dataMaterial);
130 void computeLinearStiff ( dataPtr_Type& ,
const mapMarkerVolumesPtr_Type ,
131 const mapMarkerIndexesPtr_Type );
141 void updateJacobianMatrix (
const vector_Type& disp,
142 const dataPtr_Type& dataMaterial,
143 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
144 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
145 const displayerPtr_Type& displayer);
156 void updateNonLinearJacobianTerms ( matrixPtr_Type& jacobian,
157 const vector_Type& disp,
158 const dataPtr_Type& dataMaterial,
159 const mapMarkerVolumesPtr_Type ,
160 const mapMarkerIndexesPtr_Type ,
161 const displayerPtr_Type& displayer);
172 void computeStiffness (
const vector_Type& disp,
Real factor,
const dataPtr_Type& dataMaterial,
173 const mapMarkerVolumesPtr_Type ,
174 const mapMarkerIndexesPtr_Type ,
175 const displayerPtr_Type& displayer );
201 void computeKinematicsVariables (
const VectorElemental& dk_loc ) {}
204 void showMe ( std::string
const& fileNameVectStiff,
205 std::string
const& fileNameJacobain);
216 void computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
217 const Epetra_SerialDenseMatrix& tensorF,
218 const Epetra_SerialDenseMatrix& cofactorF,
219 const std::vector<Real>& invariants,
229 void computeCauchyStressTensor (
const vectorPtr_Type disp,
231 vectorPtr_Type sigma_1,
232 vectorPtr_Type sigma_2,
233 vectorPtr_Type sigma_3);
241 matrixPtr_Type
const stiffMatrix()
const 243 return super::M_jacobian;
247 vectorPtr_Type
const stiffVector()
const 252 void apply (
const vector_Type& sol, vector_Type& res,
253 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
254 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
255 const displayerPtr_Type displayer);
268 void setupVectorsParameters (
void );
271 vectorPtr_Type M_stiff;
274 matrixSmall_Type M_identity;
278 template <
typename MeshType>
279 NeoHookeanMaterialNonLinear<MeshType>::NeoHookeanMaterialNonLinear() :
290 template <
typename MeshType>
291 NeoHookeanMaterialNonLinear<MeshType>::~NeoHookeanMaterialNonLinear()
298 template <
typename MeshType>
300 NeoHookeanMaterialNonLinear<MeshType>::setup (
const FESpacePtr_Type& dFESpace,
301 const ETFESpacePtr_Type& dETFESpace,
302 const std::shared_ptr<
const MapEpetra>& monolithicMap,
304 const dataPtr_Type& dataMaterial)
306 this->M_dataMaterial = dataMaterial;
307 this->M_dispFESpace = dFESpace;
308 this->M_dispETFESpace = dETFESpace;
309 this->M_localMap = monolithicMap;
310 this->M_offset = offset;
311 M_stiff.reset (
new vector_Type (*
this->M_localMap) );
313 M_identity (0, 0) = 1.0;
314 M_identity (0, 1) = 0.0;
315 M_identity (0, 2) = 0.0;
316 M_identity (1, 0) = 0.0;
317 M_identity (1, 1) = 1.0;
318 M_identity (1, 2) = 0.0;
319 M_identity (2, 0) = 0.0;
320 M_identity (2, 1) = 0.0;
321 M_identity (2, 2) = 1.0;
326 this->M_vectorsParameters.reset (
new vectorsParameters_Type ( 2 ) );
328 this->setupVectorsParameters();
331 template <
typename MeshType>
333 NeoHookeanMaterialNonLinear<MeshType>::setupVectorsParameters (
void )
342 UInt nbElements =
this->M_dispFESpace->mesh()->numVolumes();
346 (* (
this->M_vectorsParameters) ) [0].resize ( nbElements );
349 (* (
this->M_vectorsParameters) ) [1].resize ( nbElements );
351 for (
UInt i (0); i < nbElements; i++ )
354 UInt markerID =
this->M_dispFESpace->mesh()->element ( i ).markerID();
356 Real mu =
this->M_dataMaterial->mu ( markerID );
357 Real bulk =
this->M_dataMaterial->bulk ( markerID );
359 ( (* (
this->M_vectorsParameters) ) [0]) [ i ] = mu;
360 ( (* (
this->M_vectorsParameters) ) [1]) [ i ] = bulk;
365 template <
typename MeshType>
366 void NeoHookeanMaterialNonLinear<MeshType>::computeLinearStiff (dataPtr_Type& ,
367 const mapMarkerVolumesPtr_Type ,
368 const mapMarkerIndexesPtr_Type )
374 template <
typename MeshType>
375 void NeoHookeanMaterialNonLinear<MeshType>::updateJacobianMatrix (
const vector_Type& disp,
376 const dataPtr_Type& dataMaterial,
377 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
378 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
379 const displayerPtr_Type& displayer )
381 this->M_jacobian.reset (
new matrix_Type (*
this->M_localMap) );
383 displayer->leaderPrint (
" \n*********************************\n ");
384 updateNonLinearJacobianTerms (
this->M_jacobian, disp,
this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
385 displayer->leaderPrint (
" \n*********************************\n ");
389 template <
typename MeshType>
390 void NeoHookeanMaterialNonLinear<MeshType>::updateNonLinearJacobianTerms ( matrixPtr_Type& jacobian,
391 const vector_Type& disp,
392 const dataPtr_Type& dataMaterial,
393 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
394 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
395 const displayerPtr_Type& displayer )
400 displayer->leaderPrint (
" Non-Linear S- updating non linear terms in the Jacobian Matrix (Neo-Hookean)");
427 tensorF_Type F = ExpressionDefinitions::deformationGradient(
this->M_dispETFESpace, disp,
this->M_offset,
this->M_identity );
442 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
443 this->M_dispFESpace->qr(),
444 this->M_dispETFESpace,
445 this->M_dispETFESpace,
449 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
450 this->M_dispFESpace->qr(),
451 this->M_dispETFESpace,
452 this->M_dispETFESpace,
453 value ( - 1.0 / 2.0 ) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * ( pow (J, 2.0) - J + log (J) ) * dot ( F_T * transpose (grad (phi_j) ) * F_T, grad (phi_i) )
459 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
460 this->M_dispFESpace->qr(),
461 this->M_dispETFESpace,
462 this->M_dispETFESpace,
463 value (-2.0 / 3.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) * pow (J, - (2.0 / 3.0) ) * dot ( F_T ,
grad (phi_j) ) * dot ( F ,
grad (phi_i) )
468 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
469 this->M_dispFESpace->qr(),
470 this->M_dispETFESpace,
471 this->M_dispETFESpace,
472 value (2.0 / 9.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) *
477 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
478 this->M_dispFESpace->qr(),
479 this->M_dispETFESpace,
480 this->M_dispETFESpace,
481 parameter ( (* (
this->M_vectorsParameters) ) [0] ) * pow (J, - (2.0 / 3.0) ) * dot ( grad (phi_j), grad (phi_i) )
485 integrate ( elements (
this->M_dispETFESpace->mesh() ),
486 this->M_dispFESpace->qr(),
487 this->M_dispETFESpace,
488 this->M_dispETFESpace,
489 value (-2.0 / 3.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) *
496 integrate ( elements (
this->M_dispETFESpace->mesh() ),
497 this->M_dispFESpace->qr(),
498 this->M_dispETFESpace,
499 this->M_dispETFESpace,
500 value (1.0 / 3.0) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) *
501 pow( J, -2.0/3.0) * I_C * dot ( ( F_T * transpose (grad (phi_j) ) * F_T ), grad (phi_i) )
506 jacobian->globalAssemble();
510 template <
typename MeshType>
511 void NeoHookeanMaterialNonLinear<MeshType>::apply (
const vector_Type& sol, vector_Type& res,
512 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
513 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
514 const displayerPtr_Type displayer)
516 computeStiffness (sol, 0.,
this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
521 template <
typename MeshType>
522 void NeoHookeanMaterialNonLinear<MeshType>::computeStiffness (
const vector_Type& disp,
524 const dataPtr_Type& dataMaterial,
525 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
526 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
527 const displayerPtr_Type& displayer )
531 this->M_stiff.reset (
new vector_Type (*
this->M_localMap) );
533 displayer->leaderPrint (
" \n******************************************************************\n ");
534 displayer->leaderPrint (
" Non-Linear S- Computing the Neo-Hookean nonlinear stiffness vector" );
535 displayer->leaderPrint (
" \n******************************************************************\n ");
537 M_stiff.reset (
new vector_Type (*
this->M_localMap) );
563 tensorF_Type F = ExpressionDefinitions::deformationGradient(
this->M_dispETFESpace, disp,
this->M_offset,
this->M_identity );
579 integrate ( elements (
this->M_dispETFESpace->mesh() ),
580 this->M_dispFESpace->qr(),
581 this->M_dispETFESpace,
582 value (1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * ( pow ( J , 2.0) - J + log (J) ) * dot ( F_T,
grad (phi_i) )
586 integrate ( elements (
this->M_dispETFESpace->mesh() ),
587 this->M_dispFESpace->qr(),
588 this->M_dispETFESpace,
589 parameter ( (* (
this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) * (dot ( F -
value (1.0 / 3.0
) * I_C * F_T,
grad (phi_i) ) )
593 this->M_stiff->globalAssemble();
596 template <
typename MeshType>
604 template <
typename MeshType>
605 void NeoHookeanMaterialNonLinear<MeshType>::computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
606 const Epetra_SerialDenseMatrix& tensorF,
607 const Epetra_SerialDenseMatrix& cofactorF,
608 const std::vector<Real>& invariants,
613 Real mu =
this->M_dataMaterial->mu (marker);
614 Real bulk =
this->M_dataMaterial->bulk (marker);
617 Epetra_SerialDenseMatrix firstTerm (tensorF);
618 Epetra_SerialDenseMatrix copyCofactorF (cofactorF);
620 scale = -1 * (1.0 / 3.0) * invariants[0];
621 copyCofactorF.Scale ( scale );
622 firstTerm += copyCofactorF;
625 coef = mu * std::pow (invariants[3], -2.0 / 3.0);
626 firstTerm.Scale ( coef );
629 Epetra_SerialDenseMatrix secondTerm (cofactorF);
631 sCoef = invariants[3] * (bulk / 2.0) * (invariants[3] - 1 + (1 / invariants[3]) * std::log (invariants[3]) );
632 secondTerm.Scale ( sCoef );
634 firstPiola += firstTerm;
635 firstPiola += secondTerm;
639 template <
typename MeshType>
640 void NeoHookeanMaterialNonLinear<MeshType>::computeCauchyStressTensor (
const vectorPtr_Type disp,
642 vectorPtr_Type sigma_1,
643 vectorPtr_Type sigma_2,
644 vectorPtr_Type sigma_3)
652 tensorF_Type F = ExpressionDefinitions::deformationGradient(
this->M_dispETFESpace, *disp,
this->M_offset,
this->M_identity );
666 evaluateNode( elements (
this->M_dispETFESpace->mesh() ),
668 this->M_dispETFESpace,
669 meas_K * dot ( vectorFromMatrix( ( 1 / J )*
670 (
value (1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * ( pow ( J , 2.0) - J + log (J) ) * F_T +
671 parameter ( (* (
this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) * ( F -
value (1.0 / 3.0
) * I_C * F_T )
673 transpose(F), 0 ),
phi_i)
675 sigma_1->globalAssemble();
678 evaluateNode( elements (
this->M_dispETFESpace->mesh() ),
680 this->M_dispETFESpace,
681 meas_K * dot ( vectorFromMatrix( ( 1 / J )*
682 (
value (1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * ( pow ( J , 2.0) - J + log (J) ) * F_T +
683 parameter ( (* (
this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) * ( F -
value (1.0 / 3.0
) * I_C * F_T ) ) *
684 transpose( F ) , 1 ),
phi_i)
686 sigma_2->globalAssemble();
688 evaluateNode( elements (
this->M_dispETFESpace->mesh() ),
690 this->M_dispETFESpace,
691 meas_K * dot ( vectorFromMatrix( ( 1 / J )*
692 (
value (1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * ( pow ( J , 2.0) - J + log (J) ) * F_T +
693 parameter ( (* (
this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) * ( F -
value (1.0 / 3.0
) * I_C * F_T )
695 transpose( F ) , 2 ),
phi_i)
697 sigma_3->globalAssemble();
702 template <
typename MeshType>
703 inline StructuralIsotropicConstitutiveLaw<MeshType>* createNeoHookeanMaterialNonLinear()
705 return new NeoHookeanMaterialNonLinear<MeshType >();
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
rightCauchyGreenTensor_Type tensorC(const ExpressionTranspose< deformationGradient_Type > tF, const deformationGradient_Type F)
determinantTensorF_Type determinantF(const deformationGradient_Type F)
static const LifeV::UInt elm_nodes_num[]
const ExpressionPhiJ phi_j
Simple function to be used in the construction of an expression.
minusTransposedTensor_Type minusT(const deformationGradient_Type F)
traceTensor_Type traceTensor(const rightCauchyGreenTensor_Type C)
double Real
Generic real data.
ExpressionScalar value(const Real &myValue)
Simple function to be used in the construction of an expression.
const ExpressionMeas meas_K
Instance to be used in the expressions.
QuadratureRule - The basis class for storing and accessing quadrature rules.
ExpressionDphiJ grad(const ExpressionPhiJ &)
Simple function to be used in the construction of an expression.
DataElasticStructure - Data container for solid problems with elastic structure.
const ExpressionPhiI phi_i
Simple function to be used in the construction of an expression.
uint32_type UInt
generic unsigned integer (used mainly for addressing)
ExpressionDphiI grad(const ExpressionPhiI &)
Simple function to be used in the construction of an expression.