38 #ifndef _NONLINVENANTKIRCHHOFFMATERIAL_H_ 39 #define _NONLINVENANTKIRCHHOFFMATERIAL_H_ 41 #include <lifev/structure/solver/isotropic/StructuralIsotropicConstitutiveLaw.hpp> 45 template <
typename MeshType>
46 class VenantKirchhoffMaterialNonLinear :
47 public StructuralIsotropicConstitutiveLaw<MeshType>
53 typedef StructuralIsotropicConstitutiveLaw<MeshType> super;
55 typedef typename super::data_Type data_Type;
57 typedef typename super::vector_Type vector_Type;
58 typedef typename super::matrix_Type matrix_Type;
60 typedef typename super::matrixPtr_Type matrixPtr_Type;
61 typedef typename super::vectorPtr_Type vectorPtr_Type;
62 typedef typename super::dataPtr_Type dataPtr_Type;
63 typedef typename super::displayerPtr_Type displayerPtr_Type;
65 typedef typename super::mapMarkerVolumesPtr_Type mapMarkerVolumesPtr_Type;
66 typedef typename super::mapMarkerVolumes_Type mapMarkerVolumes_Type;
67 typedef typename mapMarkerVolumes_Type::const_iterator mapIterator_Type;
69 typedef typename super::vectorVolumes_Type vectorVolumes_Type;
70 typedef std::shared_ptr<vectorVolumes_Type> vectorVolumesPtr_Type;
72 typedef std::vector<UInt> vectorIndexes_Type;
73 typedef std::shared_ptr<vectorIndexes_Type> vectorIndexesPtr_Type;
74 typedef std::map< UInt, vectorIndexes_Type> mapMarkerIndexes_Type;
75 typedef std::shared_ptr<mapMarkerIndexes_Type> mapMarkerIndexesPtr_Type;
76 typedef typename mapMarkerIndexes_Type::const_iterator mapIteratorIndex_Type;
78 typedef typename super::FESpacePtr_Type FESpacePtr_Type;
79 typedef typename super::ETFESpacePtr_Type ETFESpacePtr_Type;
82 typedef typename super::vectorsParameters_Type vectorsParameters_Type;
83 typedef typename super::vectorsParametersPtr_Type vectorsParametersPtr_Type;
85 typedef MatrixSmall<3, 3> matrixSmall_Type;
88 typedef typename super::tensorF_Type tensorF_Type;
89 typedef typename super::determinantF_Type determinantF_Type;
90 typedef typename super::tensorC_Type tensorC_Type;
91 typedef typename super::minusT_Type minusT_Type;
92 typedef typename super::traceTensor_Type traceTensor_Type;
100 VenantKirchhoffMaterialNonLinear();
102 virtual ~VenantKirchhoffMaterialNonLinear();
117 void setup (
const FESpacePtr_Type& dFESpace,
118 const ETFESpacePtr_Type& dETFESpace,
119 const std::shared_ptr<
const MapEpetra>& monolithicMap,
120 const UInt offset,
const dataPtr_Type& dataMaterial);
126 void computeLinearStiff ( dataPtr_Type& ,
const mapMarkerVolumesPtr_Type ,
const mapMarkerIndexesPtr_Type );
136 void updateJacobianMatrix (
const vector_Type& disp,
137 const dataPtr_Type& dataMaterial,
138 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
139 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
140 const displayerPtr_Type& displayer );
151 void updateNonLinearJacobianTerms ( matrixPtr_Type& jacobian,
152 const vector_Type& disp,
153 const dataPtr_Type& dataMaterial,
154 const mapMarkerVolumesPtr_Type ,
155 const mapMarkerIndexesPtr_Type ,
156 const displayerPtr_Type& displayer );
168 void computeStiffness (
const vector_Type& disp,
Real factor,
const dataPtr_Type& dataMaterial,
const mapMarkerVolumesPtr_Type ,
169 const mapMarkerIndexesPtr_Type ,
const displayerPtr_Type& displayer );
195 void computeKinematicsVariables (
const VectorElemental& )
200 void showMe ( std::string
const& fileNameVectStiff,
201 std::string
const& fileNameJacobain );
209 matrixPtr_Type
const stiffMatrix()
const 211 return super::M_jacobian;
216 vectorPtr_Type
const stiffVector()
const 221 void apply (
const vector_Type& sol, vector_Type& res,
222 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
223 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
224 const displayerPtr_Type displayer) ;
234 void computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
235 const Epetra_SerialDenseMatrix& tensorF,
236 const Epetra_SerialDenseMatrix& cofactorF,
237 const std::vector<Real>& invariants,
246 void computeCauchyStressTensor (
const vectorPtr_Type disp,
248 vectorPtr_Type sigma_1,
249 vectorPtr_Type sigma_2,
250 vectorPtr_Type sigma_3);
264 void setupVectorsParameters (
void );
268 vectorPtr_Type M_stiff;
271 matrixSmall_Type M_identity;
279 template <
typename MeshType>
280 VenantKirchhoffMaterialNonLinear<MeshType>::VenantKirchhoffMaterialNonLinear() :
291 template <
typename MeshType>
292 VenantKirchhoffMaterialNonLinear<MeshType>::~VenantKirchhoffMaterialNonLinear()
299 template <
typename MeshType>
301 VenantKirchhoffMaterialNonLinear<MeshType>::setup (
const FESpacePtr_Type& dFESpace,
302 const ETFESpacePtr_Type& dETFESpace,
303 const std::shared_ptr<
const MapEpetra>& monolithicMap,
305 const dataPtr_Type& dataMaterial)
307 this->M_dataMaterial = dataMaterial;
308 this->M_dispFESpace = dFESpace;
309 this->M_dispETFESpace = dETFESpace;
310 this->M_localMap = monolithicMap;
311 this->M_offset = offset;
313 M_stiff.reset (
new vector_Type (*
this->M_localMap) );
315 M_identity (0, 0) = 1.0;
316 M_identity (0, 1) = 0.0;
317 M_identity (0, 2) = 0.0;
318 M_identity (1, 0) = 0.0;
319 M_identity (1, 1) = 1.0;
320 M_identity (1, 2) = 0.0;
321 M_identity (2, 0) = 0.0;
322 M_identity (2, 1) = 0.0;
323 M_identity (2, 2) = 1.0;
328 this->M_vectorsParameters.reset (
new vectorsParameters_Type ( 2 ) );
330 this->setupVectorsParameters( );
333 template <
typename MeshType>
335 VenantKirchhoffMaterialNonLinear<MeshType>::setupVectorsParameters (
void )
344 UInt nbElements =
this->M_dispFESpace->mesh()->numVolumes();
348 (* (
this->M_vectorsParameters) ) [0].resize ( nbElements );
351 (* (
this->M_vectorsParameters) ) [1].resize ( nbElements );
353 for (
UInt i (0); i < nbElements; i++ )
356 UInt markerID =
this->M_dispFESpace->mesh()->element ( i ).markerID();
358 Real lambda =
this->M_dataMaterial->lambda ( markerID );
359 Real mu =
this->M_dataMaterial->mu ( markerID );
362 ( (* (
this->M_vectorsParameters) ) [0]) [ i ] = lambda;
363 ( (* (
this->M_vectorsParameters) ) [1]) [ i ] = mu;
368 template <
typename MeshType>
369 void VenantKirchhoffMaterialNonLinear<MeshType>::computeLinearStiff (dataPtr_Type& ,
370 const mapMarkerVolumesPtr_Type ,
371 const mapMarkerIndexesPtr_Type )
380 template <
typename MeshType>
381 void VenantKirchhoffMaterialNonLinear<MeshType>::updateJacobianMatrix (
const vector_Type& disp,
382 const dataPtr_Type& dataMaterial,
383 const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
384 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
385 const displayerPtr_Type& displayer )
388 this->M_jacobian.reset (
new matrix_Type (*
this->M_localMap) );
390 displayer->leaderPrint (
" \n*********************************\n ");
391 updateNonLinearJacobianTerms (
this->M_jacobian, disp,
this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
392 displayer->leaderPrint (
" \n*********************************\n ");
400 template <
typename MeshType>
401 void VenantKirchhoffMaterialNonLinear<MeshType>::updateNonLinearJacobianTerms ( matrixPtr_Type& jacobian,
402 const vector_Type& disp,
403 const dataPtr_Type& ,
404 const mapMarkerVolumesPtr_Type ,
405 const mapMarkerIndexesPtr_Type ,
406 const displayerPtr_Type& displayer )
411 displayer->leaderPrint (
" Non-Linear S- updating non linear terms in the Jacobian Matrix (St. Venant-Kirchhoff)");
440 tensorF_Type F = ExpressionDefinitions::deformationGradient(
this->M_dispETFESpace, disp,
this->M_offset,
this->M_identity );
450 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
451 this->M_dispFESpace->qr(),
452 this->M_dispETFESpace,
453 this->M_dispETFESpace,
458 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
459 this->M_dispFESpace->qr(),
460 this->M_dispETFESpace,
461 this->M_dispETFESpace,
462 value ( 1.0 / 2.0 ) *parameter ( (* (
this->M_vectorsParameters) ) [0] ) * ( I_C - value (3.0) ) * dot ( grad (phi_j), grad (phi_i) )
466 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
467 this->M_dispFESpace->qr(),
468 this->M_dispETFESpace,
469 this->M_dispETFESpace,
470 value (-1.0) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * dot ( grad (phi_j), grad (phi_i) )
474 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
475 this->M_dispFESpace->qr(),
476 this->M_dispETFESpace,
477 this->M_dispETFESpace,
482 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
483 this->M_dispFESpace->qr(),
484 this->M_dispETFESpace,
485 this->M_dispETFESpace,
486 parameter ( (* (
this->M_vectorsParameters) ) [1]) * dot ( F * transpose (grad (phi_j) ) * F , grad (phi_i) )
490 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
491 this->M_dispFESpace->qr(),
492 this->M_dispETFESpace,
493 this->M_dispETFESpace,
494 parameter ( (* (
this->M_vectorsParameters) ) [1] ) * dot ( F * transpose (F) *
grad (phi_j) ,
grad (phi_i) )
497 jacobian->globalAssemble();
504 template <
typename MeshType>
505 void VenantKirchhoffMaterialNonLinear<MeshType>::computeStiffness (
const vector_Type& disp,
507 const dataPtr_Type& ,
508 const mapMarkerVolumesPtr_Type ,
509 const mapMarkerIndexesPtr_Type ,
510 const displayerPtr_Type& displayer )
515 M_stiff.reset (
new vector_Type (*
this->M_localMap) );
518 displayer->leaderPrint (
" \n*********************************\n ");
519 displayer->leaderPrint (
" Non-Linear S- Computing the Exponential nonlinear stiffness vector ");
520 displayer->leaderPrint (
" \n*********************************\n ");
545 tensorF_Type F = ExpressionDefinitions::deformationGradient(
this->M_dispETFESpace, disp,
this->M_offset,
this->M_identity );
555 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
556 this->M_dispFESpace->qr(),
557 this->M_dispETFESpace,
558 value ( 1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) * ( I_C - 3.0 ) * dot ( F,
grad (phi_i) )
563 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
564 this->M_dispFESpace->qr(),
565 this->M_dispETFESpace,
566 value (-1.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * dot ( F ,
grad (phi_i) )
571 integrate ( elements (
this->M_dispETFESpace->mesh() ) ,
572 this->M_dispFESpace->qr(),
573 this->M_dispETFESpace,
574 parameter ( (* (
this->M_vectorsParameters) ) [1] ) * dot ( F * C ,
grad (phi_i) )
578 this->M_stiff->globalAssemble();
582 template <
typename MeshType>
591 template <
typename MeshType>
592 void VenantKirchhoffMaterialNonLinear<MeshType>::apply (
const vector_Type& sol,
593 vector_Type& res,
const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
594 const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
const displayerPtr_Type displayer )
596 computeStiffness (sol, 0,
this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
600 template <
typename MeshType>
601 void VenantKirchhoffMaterialNonLinear<MeshType>::computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
602 const Epetra_SerialDenseMatrix& tensorF,
603 const Epetra_SerialDenseMatrix& ,
604 const std::vector<Real>& invariants,
609 Real lambda =
this->M_dataMaterial->lambda (marker);
610 Real mu =
this->M_dataMaterial->mu (marker);
611 Real coef = ( lambda / 2.0 ) * ( invariants[0] - 3.0 );
613 Epetra_SerialDenseMatrix firstTerm (tensorF);
614 firstTerm.Scale (coef);
616 Epetra_SerialDenseMatrix secondTerm (tensorF);
617 Real coeff = -1.0 * mu;
618 secondTerm.Scale ( coeff );
620 Epetra_SerialDenseMatrix thirdTerm (
this->M_dispFESpace->fieldDim(),
this->M_dispFESpace->fieldDim() );
621 Epetra_SerialDenseMatrix rightCauchyC (
this->M_dispFESpace->fieldDim(),
this->M_dispFESpace->fieldDim() );
622 rightCauchyC.Scale (0.0);
626 rightCauchyC.Multiply (
'T',
'N', 1.0, tensorF, tensorF, 0.0);
628 thirdTerm.Multiply (
'N',
'N', 1.0, tensorF, rightCauchyC, 0.0);
629 thirdTerm.Scale (mu);
631 firstPiola.Scale (0.0);
633 firstPiola += firstTerm;
634 firstPiola += secondTerm;
635 firstPiola += thirdTerm;
638 template <
typename MeshType>
639 void VenantKirchhoffMaterialNonLinear<MeshType>::computeCauchyStressTensor (
const vectorPtr_Type disp,
641 vectorPtr_Type sigma_1,
642 vectorPtr_Type sigma_2,
643 vectorPtr_Type sigma_3)
650 tensorF_Type F = ExpressionDefinitions::deformationGradient(
this->M_dispETFESpace, *disp,
this->M_offset,
this->M_identity );
665 evaluateNode( elements (
this->M_dispETFESpace->mesh() ),
667 this->M_dispETFESpace,
668 meas_K * dot ( vectorFromMatrix( ( 1 / J )*
669 (
value ( 1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) * ( I_C - 3.0 ) * F +
670 value (-1.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * F +
671 parameter ( (* (
this->M_vectorsParameters) ) [1] ) * F * C
673 transpose( F ), 0 ),
phi_i)
675 sigma_1->globalAssemble();
677 evaluateNode( elements (
this->M_dispETFESpace->mesh() ),
679 this->M_dispETFESpace,
680 meas_K * dot ( vectorFromMatrix( ( 1 / J )*
681 (
value ( 1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) * ( I_C - 3.0 ) * F +
682 value (-1.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * F +
683 parameter ( (* (
this->M_vectorsParameters) ) [1] ) * F * C
685 transpose( F ) , 1 ),
phi_i)
687 sigma_2->globalAssemble();
689 evaluateNode( elements (
this->M_dispETFESpace->mesh() ),
691 this->M_dispETFESpace,
692 meas_K * dot ( vectorFromMatrix( ( 1 / J )*
693 (
value ( 1.0 / 2.0
) * parameter ( (* (
this->M_vectorsParameters) ) [0] ) * ( I_C - 3.0 ) * F +
694 value (-1.0
) * parameter ( (* (
this->M_vectorsParameters) ) [1] ) * F +
695 parameter ( (* (
this->M_vectorsParameters) ) [1] ) * F * C
697 transpose( F ) , 2 ),
phi_i)
699 sigma_3->globalAssemble();
704 template <
typename MeshType>
705 inline StructuralIsotropicConstitutiveLaw<MeshType>* createVenantKirchhoffNonLinear()
707 return new VenantKirchhoffMaterialNonLinear<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.
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.