LifeV
DistributedHolzapfelMaterialNonLinear.hpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published
12  by the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  * @file
29  * @brief This file contains the definition for the St. Venant Kirchhoff linear material
30  *
31  * @version 1.0
32  * @date 29-07-2010
33  * @author Paolo Tricerri
34  *
35  * @maintainer Paolo Tricerri <paolo.tricerri@epfl.ch>
36  */
37 
38 #ifndef _DISTRIBUTEDHOLZAPFELMATERIAL_H_
39 #define _DISTRIBUTEDHOLZAPFELMATERIAL_H_
40 
41 #pragma GCC diagnostic ignored "-Wunused-variable"
42 #pragma GCC diagnostic ignored "-Wunused-parameter"
43 
44 
45 #include <lifev/structure/solver/anisotropic/StructuralAnisotropicConstitutiveLaw.hpp>
46 
47 #define PI 3.14159265359
48 
49 namespace LifeV
50 {
51 template <typename MeshType>
52 class DistributedHolzapfelMaterialNonLinear : public StructuralAnisotropicConstitutiveLaw<MeshType>
53 {
54  //!@name Type definitions
55  //@{
56 
57 public:
58  typedef StructuralAnisotropicConstitutiveLaw<MeshType> super;
59 
60  typedef typename super::data_Type data_Type;
61 
62  typedef typename super::vector_Type vector_Type;
63  typedef typename super::matrix_Type matrix_Type;
64 
65  typedef typename super::matrixPtr_Type matrixPtr_Type;
66  typedef typename super::vectorPtr_Type vectorPtr_Type;
67  typedef typename super::dataPtr_Type dataPtr_Type;
69 
73 
77 
80 
83 
84  typedef typename super::FESpace_Type FESpace_Type;
87 
88  //Vector for vector parameters
91 
92  typedef MatrixSmall<3, 3> matrixSmall_Type;
93 
96 
99 
100 
101  // Typedefs for expression definitions
102  typedef typename super::tensorF_Type tensorF_Type;
104  typedef typename super::tensorC_Type tensorC_Type;
105  typedef typename super::minusT_Type minusT_Type;
109 
110  // Anisotropic typedefs
113  typedef typename super::stretch_Type stretch_Type;
115 
116  // typedefs specific for this model that are defined in ExpressionDefinitions.hpp
117  typedef ExpressionDistributedModel::distributedStretch_Type distributedStretch_Type;
118  typedef ExpressionDistributedModel::tensorialPart_distrType tensorialPart_distrType;
120  //@}
121 
122 
123 
124  //! @name Constructor & Destructor
125  //@{
126 
128 
130 
131  //@}
132 
133 
134 
135  //!@name Methods
136  //@{
137 
138  //! Setup the created object of the class StructuralIsotropicConstitutiveLaw
139  /*!
140  \param dFespace: the FiniteElement Space
141  \param monolithicMap: the MapEpetra
142  \param offset: the offset parameter used assembling the matrices
143  */
144  void setup ( const FESpacePtr_Type& dFESpace,
145  const ETFESpacePtr_Type& dETFESpace,
146  const std::shared_ptr<const MapEpetra>& monolithicMap,
147  const UInt offset,const dataPtr_Type& dataMaterial);
148 
149 
150  //! Compute the Stiffness matrix in StructuralSolver::buildSystem()
151  /*!
152  \param dataMaterial the class with Material properties data
153  */
154  void computeLinearStiff ( dataPtr_Type& /*dataMaterial*/,
155  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
156  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/);
157 
158 
159  //! Updates the Jacobian matrix in StructualSolver::updateJacobian
160  /*!
161  \param disp: solution at the k-th iteration of NonLinearRichardson Method
162  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get
163  the material coefficients (e.g. Young modulus, Poisson ratio..)
164  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
165  */
166  void updateJacobianMatrix ( const vector_Type& disp,
167  const dataPtr_Type& dataMaterial,
168  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
169  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
170  const displayerPtr_Type& displayer );
171 
172 
173  //! Updates the nonlinear terms in the Jacobian matrix in StructualSolver::updateJacobian
174  /*!
175  \param stiff: stiffness matrix provided from outside
176  \param disp: solution at the k-th iteration of NonLinearRichardson Method
177  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get
178  the material coefficients (e.g. Young modulus, Poisson ratio..)
179  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
180  */
182  const vector_Type& disp,
183  const dataPtr_Type& dataMaterial,
184  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
185  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/,
186  const displayerPtr_Type& displayer );
187 
188 
189  //! Interface method to compute the new Stiffness matrix in StructuralSolver::evalResidual and in
190  //! StructuralSolver::updateSystem since the matrix is the expression of the matrix is the same.
191  /*!
192  \param sol: the solution vector
193  \param factor: scaling factor used in FSI
194  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get the
195  material coefficients (e.g. Young modulus, Poisson ratio..)
196  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
197  */
198  void computeStiffness ( const vector_Type& disp,
199  const UInt iter,
200  Real factor,
201  const dataPtr_Type& dataMaterial,
202  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
203  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/,
204  const displayerPtr_Type& displayer );
205 
207  const dataPtr_Type& dataMaterial,
208  const displayerPtr_Type& displayer ){};
209 
210  //! Computes the new Stiffness vector for Neo-Hookean and Holzapfel materials in StructuralSolver
211  //! given a certain displacement field.
212  //! This function is used both in StructuralSolver::evalResidual and in StructuralSolver::updateSystem
213  //! since the matrix is the expression of the matrix is the same.
214  /*!
215  \param sol: the solution vector
216  \param factor: scaling factor used in FSI
217  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get
218  the material coefficients (e.g. Young modulus, Poisson ratio..)
219  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
220  */
221 
222  //! Computes the deformation gradient F, the cofactor matrix Cof(F),
223  //! the determinant of F (J = det(F)), the trace of right Cauchy-Green tensor tr(C)
224  //! This function is used in StructuralIsotropicConstitutiveLaw::computeStiffness
225  /*!
226  \param dk_loc: the elemental displacement
227  */
228  void computeKinematicsVariables ( const VectorElemental& dk_loc )
229  {
230  }
231 
232  //! Computes the deformation Gradient F, the cofactor of F Cof(F), the determinant of F J = det(F), the trace of C Tr(C).
233  /*!
234  \param dk_loc: local displacement vector
235  */
236  //void computeStress( const vector_Type& sol );
237 
238  //! ShowMe method of the class (saved on a file the stiffness vector and the jacobian)
239  void showMe ( std::string const& fileNameVectStiff,
240  std::string const& fileNameJacobain );
241 
242  //! Compute the First Piola Kirchhoff Tensor
243  /*!
244  \param firstPiola Epetra_SerialDenseMatrix that has to be filled
245  \param tensorF Epetra_SerialDenseMatrix the deformation gradient
246  \param cofactorF Epetra_SerialDenseMatrix cofactor of F
247  \param invariants std::vector with the invariants of C and the detF
248  \param material UInt number to get the material parameteres form the VenantElasticData class
249  */
250  void computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
251  const Epetra_SerialDenseMatrix& tensorF,
252  const Epetra_SerialDenseMatrix& cofactorF,
253  const std::vector<Real>& invariants,
254  const UInt marker);
255 
256  //! Compute the First Piola Kirchhoff Tensor
257  /*!
258  \param disp the displacement field from which we compute the fisrt piola-Kirchhoff tensor
259  \param sigma_1 the first column of the Cauchy stress tensor
260  \param sigma_2 the second column of the Cauchy stress tensor
261  \param sigma_3 the third column of the Cauchy stress tensor
262  */
263  void computeCauchyStressTensor ( const vectorPtr_Type disp,
264  const QuadratureRule& evalQuad,
265  vectorPtr_Type sigma_1,
266  vectorPtr_Type sigma_2,
267  vectorPtr_Type sigma_3);
268 
269 
270  void setupFiberDirections( vectorFiberFunctionPtr_Type vectorOfFibers );
271 
272 
273  void evaluateActivationFibers( const vector_Type& displacement,
274  vector_Type& fourthInvariant){}
275  //@}
276 
277  //! @name Get Methods
278  //@{
279 
280  //! Get the Stiffness matrix
282  {
283  return super::M_jacobian;
284  }
285 
286 
287  //! Get the stiffness vector
289  {
290  return M_stiff;
291  }
292 
293  //! Specific for multimechanism
294  vectorPtr_Type const selectionCriterion( const UInt /*i*/ ) const
295  {
296  vectorPtr_Type empty( new vector_Type( this->M_dispFESpace->map() ) );
297  return empty;
298  }
299 
300 
301  vectorPtr_Type const activationDisplacement( const UInt /*i*/ ) const
302  {
303  vectorPtr_Type empty( new vector_Type( this->M_dispFESpace->map() ) );
304  return empty;
305  }
306 
307  vectorPtr_Type const activatedUnitFiber( const UInt /*i*/ ) const
308  {
309  vectorPtr_Type empty( new vector_Type( this->M_dispFESpace->map() ) );
310  return empty;
311  }
312 
313  vectorPtr_Type const activatedDeterminant( const UInt /*i*/ ) const
314  {
315  vectorPtr_Type empty( new vector_Type( this->M_dispFESpace->map() ) );
316  return empty;
317  }
318 
319 
320 
321  void apply ( const vector_Type& sol, vector_Type& res,
322  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
323  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
324  const displayerPtr_Type displayer) ;
325 
326  //@}
327 
328 
329 
330 protected:
331 
332  //! construct the vectors for the parameters
333  /*!
334  \param VOID
335  \return VOID
336  */
337  // void setupVectorsParameters ( void );
338 
339 
340  //! Vector: stiffness non-linear
342 
343  //Create the indentity for F
345 
346 };
347 
348 
349 
350 
351 
352 template <typename MeshType>
354  super ( ),
355  M_stiff ( ),
356  M_identity ( )
357 {
358 }
359 
360 
361 
362 
363 
364 template <typename MeshType>
366 {}
367 
368 
369 
370 
371 
372 template <typename MeshType>
373 void
375  const ETFESpacePtr_Type& dETFESpace,
376  const std::shared_ptr<const MapEpetra>& monolithicMap,
377  const UInt offset,
378  const dataPtr_Type& dataMaterial)
379 {
380 
381  this->M_dataMaterial = dataMaterial;
382  this->M_dispFESpace = dFESpace;
383  this->M_dispETFESpace = dETFESpace;
384  this->M_localMap = monolithicMap;
385  this->M_offset = offset;
386 
387  M_stiff.reset (new vector_Type (*this->M_localMap) );
388 
389  M_identity (0, 0) = 1.0;
390  M_identity (0, 1) = 0.0;
391  M_identity (0, 2) = 0.0;
392  M_identity (1, 0) = 0.0;
393  M_identity (1, 1) = 1.0;
394  M_identity (1, 2) = 0.0;
395  M_identity (2, 0) = 0.0;
396  M_identity (2, 1) = 0.0;
397  M_identity (2, 2) = 1.0;
398 
399  this->M_epsilon = this->M_dataMaterial->smoothness();
400  //this->setupVectorsParameters( );
401 }
402 
403 
404 // template <typename MeshType>
405 // void
406 // DistributedHolzapfelMaterialNonLinear<MeshType>::setupVectorsParameters ( void )
407 // {
408 // // Paolo Tricerri: February, 20th
409 // // In each class, the name of the parameters has to inserted in the law
410 // // TODO: move the saving of the material parameters to more abstract objects
411 // // such that in the class of the material we do not need to call explicitly
412 // // the name of the parameter.
413 
414 // // Number of volume on the local part of the mesh
415 // UInt nbElements = this->M_dispFESpace->mesh()->numVolumes();
416 
417 // // Parameter alpha
418 // // 1. resize the vector in the first element of the vector.
419 // (* (this->M_vectorsParameters) ) [0].resize ( nbElements );
420 
421 // // Parameter gamma
422 // (* (this->M_vectorsParameters) ) [1].resize ( nbElements );
423 
424 // // Parameter bulk
425 // (* (this->M_vectorsParameters) ) [2].resize ( nbElements );
426 
427 // for (UInt i (0); i < nbElements; i++ )
428 // {
429 // // Extracting the marker
430 // UInt markerID = this->M_dispFESpace->mesh()->element ( i ).markerID();
431 
432 // Real alpha = this->M_dataMaterial->alpha ( markerID );
433 // Real gamma = this->M_dataMaterial->gamma ( markerID );
434 // Real bulk = this->M_dataMaterial->bulk ( markerID );
435 
436 // ( (* (this->M_vectorsParameters) ) [0]) [ i ] = alpha;
437 // ( (* (this->M_vectorsParameters) ) [1]) [ i ] = gamma;
438 // ( (* (this->M_vectorsParameters) ) [2]) [ i ] = bulk;
439 // }
440 // }
441 
442 
443 template <typename MeshType>
444 void
446 {
447  // Allocating the right space for the vector of fiber function
448  this->M_vectorOfFibers.reset( new vectorFiberFunction_Type( ) );
449 
450  // In this method, the vector of fiber functions has to be properly set in the main
451  // of the test. The functions are then projected on a FESpace for the integration
452 
453  // Number of volume on the local part of the mesh
454  UInt nbFamilies = (*vectorOfFibers).size();
455 
456  ASSERT( nbFamilies == this->M_dataMaterial->numberFibersFamilies(),
457  " The number of families set in the test is different from the one in the data" );
458 
459  this->M_vectorOfFibers->resize( nbFamilies );
460 
461  for( UInt k(0); k < nbFamilies; k++ )
462  {
463  ( *(this->M_vectorOfFibers) )[ k ] = ( *vectorOfFibers )[ k ];
464  }
465 
466  // Setting the vectors that will be used
467  this->M_vectorInterpolated.resize( nbFamilies );
468 
469  for( UInt k(0); k < nbFamilies; k++ )
470  {
471  this->M_vectorInterpolated[ k ].reset( new vector_Type(*this->M_localMap) );
472  this->M_dispFESpace->interpolate ( *( ( *(this->M_vectorOfFibers) )[ k ] ) ,
473  * ( ( this->M_vectorInterpolated )[ k ] ),
474  0.0 );
475  }
476 
477 }
478 
479 
480 template <typename MeshType>
482  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
483  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/)
484 {
485  //! Empty method for exponential material
486 }
487 
488 
489 
490 
491 
492 template <typename MeshType>
494  const dataPtr_Type& dataMaterial,
495  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
496  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
497  const displayerPtr_Type& displayer )
498 {
499 
500  this->M_jacobian.reset (new matrix_Type (*this->M_localMap) );
501 
502  displayer->leaderPrint (" \n*********************************\n ");
503  updateNonLinearJacobianTerms (this->M_jacobian, disp, this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
504  displayer->leaderPrint (" \n*********************************\n ");
505 
506 }
507 
508 template <typename MeshType>
510  const vector_Type& disp,
511  const dataPtr_Type& dataMaterial,
512  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
513  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
514  const displayerPtr_Type& displayer )
515 {
516  using namespace ExpressionAssembly;
517 
518  // In the following definitions, the critical template argument is MapEpetra
519  // When other maps will be available in LifeV, the Holzapfel class and its mother
520  // should have as template the MeshType and the MapType.
521 
522  // Definition of F
523  tensorF_Type F = ExpressionDefinitions::deformationGradient( this->M_dispETFESpace, disp, this->M_offset, this->M_identity );
524 
525  // Definition of J
526  determinantF_Type J = ExpressionDefinitions::determinantF( F );
527 
528  // Definition of J^-(2/3) = det( C ) using the isochoric/volumetric decomposition
529  powerExpression_Type Jel = ExpressionDefinitions::powerExpression( J , (-2.0/3.0) );
530 
531  // Definition of F^-T
532  minusT_Type F_T = ExpressionDefinitions::minusT( F );
533 
534  //Definition of C
535  tensorC_Type C = ExpressionDefinitions::tensorC( transpose(F), F );
536 
537  // Definition of tr( C )
538  traceTensor_Type I_C = ExpressionDefinitions::traceTensor( C );
539 
540  // Definition of J^( -2.0/3.0 ) * tr( C )
541  isochoricTrace_Type trCbar = ExpressionDefinitions::isochoricTrace( Jel, I_C );
542 
543  // Update the heaviside function for the stretch of the fibers
544  * (jacobian) *= 0.0;
545 
546  displayer->leaderPrint (" Non-Linear S - updating non linear terms in the Jacobian Matrix (Distributed Holzapfel): \n");
547 
548  for( UInt i(0); i < this->M_vectorInterpolated.size(); i++ )
549  {
550 
551  displayer->leaderPrint (" ", i + 1,"-th fiber family \n" );
552 
553  // Material parameters
554  Real kappaI = this->M_dataMaterial->ithDistributionFibers( i );
555  Real alphaI = this->M_dataMaterial->ithStiffnessFibers( i );
556  Real gammaI = this->M_dataMaterial->ithNonlinearityFibers( i );
557 
558  // Defining the expression for the i-th fiber
559  // Definitions of the quantities which depend on the fiber directions e.g. I_4^i
560  interpolatedValue_Type fiberIth = ExpressionDefinitions::interpolateFiber( this->M_dispETFESpace, *(this->M_vectorInterpolated[ i ] ) );
561 
562  // Definition of the tensor M = ithFiber \otimes ithFiber
563  // At the moment, it's automatic that the method constructs the expression M = ithFiber \otimes ithFiber
564  // For a more general case, the file ExpressionDefinitions.hpp should be changed
565  outerProduct_Type Mith = ExpressionDefinitions::fiberTensor( fiberIth );
566 
567  // Definition of the fourth invariant : I_4^i = C:Mith
568  stretch_Type IVith = ExpressionDefinitions::fiberStretch( C, Mith );
569 
570  // Definition of the fouth isochoric invariant : J^(-2.0/3.0) * I_4^i
571  isochoricStretch_Type IVithBar = ExpressionDefinitions::isochoricFourthInvariant( Jel, IVith );
572 
573  // Definition of the quantity to measure stretch in the distributed case
574  // \kappaI * \bar{I_C} + ( 1 - 3\kappaI ) * \bar{I_4^i} - 1.0
575  distributedStretch_Type distrStretch = ExpressionDistributedModel::distributedStretch( trCbar, IVithBar, kappaI );
576 
577  // Definition of the matricial ( i.e. tensorial part ) of the Piola - Kirchhoff )
578  tensorialPart_distrType tensorialPart = ExpressionDistributedModel::tensorialPartPiola( kappaI, I_C, IVith, F, F_T, Mith );
579 
580  // Definition of the derivative w.r.t to F of the distributed stretch
581  derivativeDistrStretch_Type derivativeStretch = ExpressionDistributedModel::derivativeDistributedStretch( kappaI, trCbar, IVithBar,
582  Jel, F, F_T, Mith);
583 
584  // first term:
585  integrate( elements ( this->M_dispETFESpace->mesh() ),
586  this->M_dispFESpace->qr(),
587  this->M_dispETFESpace,
588  this->M_dispETFESpace,
589  value( 2.0 * alphaI ) * exp( value( gammaI ) * ( distrStretch ) * ( distrStretch ) ) *
590  Jel * distrStretch * derAtan( distrStretch , this->M_epsilon, ( 1.0 / PI ) ) *
591  derivativeStretch * dot( tensorialPart , grad(phi_i) )
592  ) >> jacobian;
593 
594  // second term:
595  integrate( elements ( this->M_dispETFESpace->mesh() ),
596  this->M_dispFESpace->qr(),
597  this->M_dispETFESpace,
598  this->M_dispETFESpace,
599  atan( distrStretch , this->M_epsilon, ( 1 / PI ), ( 1.0/2.0 ) ) *
600  value( 2.0 * alphaI ) * exp( value( gammaI ) * ( distrStretch ) * ( distrStretch ) ) *
601  Jel * derivativeStretch * dot( tensorialPart , grad(phi_i) )
602  ) >> jacobian;
603 
604  // third term:
605  integrate( elements ( this->M_dispETFESpace->mesh() ),
606  this->M_dispFESpace->qr(),
607  this->M_dispETFESpace,
608  this->M_dispETFESpace,
609  atan( distrStretch , this->M_epsilon, ( 1 / PI ), ( 1.0/2.0 ) ) *
610  value( 2.0 * alphaI ) * Jel * distrStretch *
611  exp( value( gammaI ) * ( distrStretch ) * ( distrStretch ) ) *
612  value( 2.0 * gammaI ) * distrStretch * derivativeStretch * dot( tensorialPart , grad(phi_i) )
613  ) >> jacobian;
614 
615  // fourth term:
616  integrate( elements ( this->M_dispETFESpace->mesh() ),
617  this->M_dispFESpace->qr(),
618  this->M_dispETFESpace,
619  this->M_dispETFESpace,
620  atan( distrStretch , this->M_epsilon, ( 1 / PI ), ( 1.0/2.0 ) ) *
621  value( 2.0 * alphaI ) * distrStretch *
622  exp( value( gammaI ) * ( distrStretch ) * ( distrStretch ) ) *
623  value( -2.0/3.0 ) * Jel * dot( F_T, grad(phi_j) ) * dot( tensorialPart , grad(phi_i) )
624  ) >> jacobian;
625 
626  // fifth term:
627  integrate( elements ( this->M_dispETFESpace->mesh() ),
628  this->M_dispFESpace->qr(),
629  this->M_dispETFESpace,
630  this->M_dispETFESpace,
631  atan( distrStretch , this->M_epsilon, ( 1 / PI ), ( 1.0/2.0 ) ) *
632  value( 2.0 * alphaI ) * distrStretch *
633  exp( value( gammaI ) * ( distrStretch ) * ( distrStretch ) ) * Jel *
634  dot( value( kappaI ) * grad(phi_j) +
635  value( - ( 2.0 * kappaI )/ 3.0 ) * dot( F, grad(phi_j) ) * F_T +
636  value( kappaI / 3.0 ) * I_C * F_T * transpose( grad(phi_j) ) * F_T +
637  value( 1.0 - 3.0 * kappaI ) * grad(phi_j) * Mith +
638  value( ( 1.0 - 3.0*kappaI ) / 3.0 ) * IVith * F_T * transpose( grad(phi_j) ) * F_T +
639  value(-( 1.0 - 3.0*kappaI ) / 3.0 ) * ( dot( transpose(grad(phi_j)) * F , Mith ) + dot( transpose(F)*grad(phi_j), Mith ) ) * F_T , grad(phi_i) )
640  ) >> jacobian;
641 
642 
643 
644  } // closing loop on fibers
645 
646  jacobian->globalAssemble();
647 }
648 
649 
650 
651 
652 
653 template <typename MeshType>
655  const UInt /*iter*/,
656  Real /*factor*/,
657  const dataPtr_Type& dataMaterial,
658  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
659  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/,
660  const displayerPtr_Type& displayer )
661 {
662 
663  using namespace ExpressionAssembly;
664 
665  M_stiff.reset (new vector_Type (*this->M_localMap) );
666  * (M_stiff) *= 0.0;
667 
668  displayer->leaderPrint (" \n*********************************\n ");
669  displayer->leaderPrint (" Non-Linear S- Computing the Distributed Holzapfel nonlinear stiffness vector (Distributed Holzapfel) ");
670  displayer->leaderPrint (" \n*********************************\n ");
671 
672  // For anisotropic part of the Piola-Kirchhoff is assemble summing up the parts of the
673  // Piola-Kirchhoff using the fiber index
674 
675  // Here the fiber vector at the quadrature node is deduced using the method
676  // Interpolate value of the expression template.
677 
678  // Defining quantities which depend of the displacement field and not on the fiber direction
679 
680  // Definition of F
681  tensorF_Type F = ExpressionDefinitions::deformationGradient( this->M_dispETFESpace, disp, this->M_offset, this->M_identity );
682 
683  // Definition of J
684  determinantF_Type J = ExpressionDefinitions::determinantF( F );
685 
686  //Definition of C
687  tensorC_Type C = ExpressionDefinitions::tensorC( transpose(F), F );
688 
689  // Definition I_C
690  traceTensor_Type I_C = ExpressionDefinitions::traceTensor( C );
691 
692  // Definition of J^-(2/3) = det( C ) using the isochoric/volumetric decomposition
693  powerExpression_Type Jel = ExpressionDefinitions::powerExpression( J , (-2.0/3.0) );
694 
695  // Definition of the isochoric trace of C
696  isochoricTrace_Type trCbar = ExpressionDefinitions::isochoricTrace( Jel, I_C );
697 
698  // Definition of F^-T
699  minusT_Type F_T = ExpressionDefinitions::minusT( F );
700 
701  for( UInt i(0); i < this->M_vectorInterpolated.size() ; i++ )
702  {
703  // As in other classes, the specialization of the MapType = MapEpetra makes this expression
704  // not always usable. When other maps will be available in LifeV, the class should be re-templated.
705 
706  // Defining the expression for the i-th fiber
707  // Definitions of the quantities which depend on the fiber directions e.g. I_4^i
708  interpolatedValue_Type fiberIth = ExpressionDefinitions::interpolateFiber( this->M_dispETFESpace, *(this->M_vectorInterpolated[ i ] ) );
709 
710  // Definition of the tensor M = ithFiber \otimes ithFiber
711  // At the moment, it's automatic that the method constructs the expression M = ithFiber \otimes ithFiber
712  // For a more general case, the file ExpressionDefinitions.hpp should be changed
713  outerProduct_Type Mith = ExpressionDefinitions::fiberTensor( fiberIth );
714 
715  // Definition of the fourth invariant : I_4^i = C:Mith
716  stretch_Type IVith = ExpressionDefinitions::fiberStretch( C, Mith );
717 
718  // Definition of the fouth isochoric invariant : J^(-2.0/3.0) * I_4^i
719  isochoricStretch_Type IVithBar = ExpressionDefinitions::isochoricFourthInvariant( Jel, IVith );
720 
721  // Material parameters
722  Real kappaI = this->M_dataMaterial->ithDistributionFibers( i );
723  Real alphaI = this->M_dataMaterial->ithStiffnessFibers( i );
724  Real gammaI = this->M_dataMaterial->ithNonlinearityFibers( i );
725 
726  // Definition of the activation for distributed model
727  distributedStretch_Type distrStretch = ExpressionDistributedModel::distributedStretch( trCbar, IVithBar, kappaI );
728 
729  // For this model, the activation in front of the stress tensor is of the form
730  // atan( k * \bar{I_C} + ( 1 - 3*k ) * \bar(I_4) - 1 )
731 
732  // First term:
733  // 2 alpha_i J^(-2.0/3.0) ( k \bar{I_1} + ( 1 - 3*k) \bar{I_4} - 1 ) ( k )
734  // exp( gamma_i * ( k \bar{I_1} + ( 1 - 3*k) \bar{I_4} - 1)^2 ) * ( F - 1.3 * I_1 * F^-T) : \grad phi_i
735  // where alpha_i and gamma_i are the fiber parameters and M is the 2nd order tensor of type f_i \otimes \ f_i
736 
737  // Second term:
738  // 2 alpha_i J^(-2.0/3.0) ( k \bar{I_1} + ( 1 - 3*k) \bar{I_4} - 1 ) ( 1 - 3k )
739  // exp( gamma_i * ( k \bar{I_1} + ( 1 - 3*k) \bar{I_4} - 1)^2 ) * ( FM - 1.3 * I_4 * F^-T) : \grad phi_i
740  // where alpha_i and gamma_i are the fiber parameters and M is the 2nd order tensor of type f_i \otimes \ f_i
741 
742  // Definition of the matricial ( i.e. tensorial part ) of the Piola - Kirchhoff )
743  tensorialPart_distrType tensorialPart = ExpressionDistributedModel::tensorialPartPiola( kappaI, I_C, IVith, F, F_T, Mith );
744 
745  integrate ( elements ( this->M_dispETFESpace->mesh() ),
746  this->M_dispFESpace->qr(),
747  this->M_dispETFESpace,
748  atan( distrStretch , this->M_epsilon, ( 1 / PI ), ( 1.0/2.0 ) ) *
749  value( 2.0 * alphaI ) * Jel * ( distrStretch ) *
750  exp( value( gammaI ) * ( distrStretch ) * ( distrStretch ) ) *
751  dot( tensorialPart, grad( phi_i ) )
752  ) >> this->M_stiff;
753 
754  }
755 
756  this->M_stiff->globalAssemble();
757 }
758 
759 
760 template <typename MeshType>
762  std::string const& fileNameJacobian )
763 {
764  this->M_stiff->spy (fileNameStiff);
766 }
767 
768 
769 template <typename MeshType>
771  vector_Type& res, const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
772  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
773  const displayerPtr_Type displayer)
774 {
775  computeStiffness (sol, 0, 0, this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
776  res += *M_stiff;
777 }
778 
779 
780 template <typename MeshType>
781 void DistributedHolzapfelMaterialNonLinear<MeshType>::computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
782  const Epetra_SerialDenseMatrix& tensorF,
783  const Epetra_SerialDenseMatrix& cofactorF,
784  const std::vector<Real>& invariants,
785  const UInt marker)
786 {
787  // It can be done using the evaluateNodal framework that has been shown to work
788  // for the isotropic laws!
789 }
790 
791 template <typename MeshType>
793  const QuadratureRule& evalQuad,
794  vectorPtr_Type sigma_1,
795  vectorPtr_Type sigma_2,
796  vectorPtr_Type sigma_3)
797 {
798 
799  ASSERT( 2 < 0 , "For the distributed anisotropic law the computation of the Cauchy stress has to be defined." );
800 
801 }
802 
803 
804 template <typename MeshType>
805 inline StructuralAnisotropicConstitutiveLaw<MeshType>* createDistributedHolzapfelMaterialNonLinear()
806 {
807  return new DistributedHolzapfelMaterialNonLinear<MeshType >();
808 }
809 
810 namespace
811 {
813 }
814 
815 } //Namespace LifeV
816 
817 #endif /* __EXPONENTIALMATERIAL_H */
StructuralAnisotropicConstitutiveLaw< MeshType > * createDistributedHolzapfelMaterialNonLinear()
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
ExpressionAddition< ExpressionAddition< scaledTensorF_Type, scaledTraceTimesMinusTF_Type >, ExpressionAddition< scaledFtimesM_Type, scaledFourthInvariantTimesMinusTF_Type > > tensorialPart_distrType
void updateJacobianMatrix(const vector_Type &disp, const dataPtr_Type &dataMaterial, const mapMarkerVolumesPtr_Type mapsMarkerVolumes, const mapMarkerIndexesPtr_Type mapsMarkerIndexes, const displayerPtr_Type &displayer)
Updates the Jacobian matrix in StructualSolver::updateJacobian.
void computeKinematicsVariables(const VectorElemental &dk_loc)
Computes the new Stiffness vector for Neo-Hookean and Holzapfel materials in StructuralSolver given a...
#define PI
vectorPtr_Type const activationDisplacement(const UInt) const
Get the activation displacement.
void updateNonLinearJacobianTerms(matrixPtr_Type &jacobian, const vector_Type &disp, const dataPtr_Type &dataMaterial, const mapMarkerVolumesPtr_Type, const mapMarkerIndexesPtr_Type, const displayerPtr_Type &displayer)
Updates the nonlinear terms in the Jacobian matrix in StructualSolver::updateJacobian.
matrixPtr_Type const stiffMatrix() const
Get the Stiffness matrix.
vectorPtr_Type const stiffVector() const
Get the stiffness vector.
void showMe(std::string const &fileNameVectStiff, std::string const &fileNameJacobain)
Computes the deformation Gradient F, the cofactor of F Cof(F), the determinant of F J = det(F)...
void evaluateActivationFibers(const vector_Type &displacement, vector_Type &fourthInvariant)
ExpressionDistributedModel::linearizationDistributedStretch_Type derivativeDistrStretch_Type
void computeLinearStiff(dataPtr_Type &, const mapMarkerVolumesPtr_Type, const mapMarkerIndexesPtr_Type)
Compute the Stiffness matrix in StructuralSolver::buildSystem()
void setupFiberDirections(vectorFiberFunctionPtr_Type vectorOfFibers)
vectorPtr_Type const activatedDeterminant(const UInt) const
Get the activation displacement.
#define ASSERT(X, A)
Definition: LifeAssert.hpp:90
StructuralAnisotropicConstitutiveLaw< MeshType > super
void computeLocalFirstPiolaKirchhoffTensor(Epetra_SerialDenseMatrix &firstPiola, const Epetra_SerialDenseMatrix &tensorF, const Epetra_SerialDenseMatrix &cofactorF, const std::vector< Real > &invariants, const UInt marker)
Compute the First Piola Kirchhoff Tensor.
std::vector< typename MeshType::element_Type * > vectorVolumes_Type
void apply(const vector_Type &sol, vector_Type &res, const mapMarkerVolumesPtr_Type mapsMarkerVolumes, const mapMarkerIndexesPtr_Type mapsMarkerIndexes, const displayerPtr_Type displayer)
ExpressionAddition< ExpressionProduct< ExpressionScalar, linearizationFisochoricTrace_Type >, ExpressionProduct< ExpressionScalar, linearizationFisochoricFourthInvariant_Type > > linearizationDistributedStretch_Type
void computeCauchyStressTensor(const vectorPtr_Type disp, const QuadratureRule &evalQuad, vectorPtr_Type sigma_1, vectorPtr_Type sigma_2, vectorPtr_Type sigma_3)
Compute the First Piola Kirchhoff Tensor.
vectorPtr_Type const selectionCriterion(const UInt) const
Specific for multimechanism.
double Real
Generic real data.
Definition: LifeV.hpp:175
void computeReferenceConfigurations(const vector_Type &disp, const dataPtr_Type &dataMaterial, const displayerPtr_Type &displayer)
vectorPtr_Type const activatedUnitFiber(const UInt) const
Get the activation displacement.
void setup(const FESpacePtr_Type &dFESpace, const ETFESpacePtr_Type &dETFESpace, const std::shared_ptr< const MapEpetra > &monolithicMap, const UInt offset, const dataPtr_Type &dataMaterial)
Setup the created object of the class StructuralIsotropicConstitutiveLaw.
QuadratureRule - The basis class for storing and accessing quadrature rules.
vectorPtr_Type M_stiff
construct the vectors for the parameters
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
ExpressionAddition< distributedInvariants_Type, ExpressionScalar > distributedStretch_Type
void computeStiffness(const vector_Type &disp, const UInt iter, Real factor, const dataPtr_Type &dataMaterial, const mapMarkerVolumesPtr_Type, const mapMarkerIndexesPtr_Type, const displayerPtr_Type &displayer)
Interface method to compute the new Stiffness matrix in StructuralSolver::evalResidual and in Structu...