LifeV
ExponentialMaterialNonLinear.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  * @author Gianmarco Mengaldo
35  *
36  * @maintainer Paolo Tricerri <paolo.tricerri@epfl.ch>
37  * @contributor Gianmarco Mengaldo <gianmarco.mengaldo@gmail.com>
38  */
39 
40 #ifndef _EXPONENTIALMATERIAL_H_
41 #define _EXPONENTIALMATERIAL_H_
42 
43 
44 
45 #include <lifev/structure/solver/isotropic/StructuralIsotropicConstitutiveLaw.hpp>
46 
47 namespace LifeV
48 {
49 template <typename MeshType>
50 class ExponentialMaterialNonLinear : public StructuralIsotropicConstitutiveLaw<MeshType>
51 {
52  //!@name Type definitions
53  //@{
54 
55 public:
56  typedef StructuralIsotropicConstitutiveLaw<MeshType> super;
57 
58  typedef typename super::data_Type data_Type;
59 
60  typedef typename super::vector_Type vector_Type;
61  typedef typename super::matrix_Type matrix_Type;
62 
63  typedef typename super::matrixPtr_Type matrixPtr_Type;
64  typedef typename super::vectorPtr_Type vectorPtr_Type;
65  typedef typename super::dataPtr_Type dataPtr_Type;
66  typedef typename super::displayerPtr_Type displayerPtr_Type;
67 
68  typedef typename super::mapMarkerVolumesPtr_Type mapMarkerVolumesPtr_Type;
69  typedef typename super::mapMarkerVolumes_Type mapMarkerVolumes_Type;
70  typedef typename mapMarkerVolumes_Type::const_iterator mapIterator_Type;
71 
72  typedef typename super::mapMarkerIndexesPtr_Type mapMarkerIndexesPtr_Type;
73  typedef typename super::mapMarkerIndexes_Type mapMarkerIndexes_Type;
74  typedef typename mapMarkerIndexes_Type::const_iterator mapIteratorIndex_Type;
75 
76  typedef std::vector<typename MeshType::element_Type*> vectorVolumes_Type;
77  typedef std::shared_ptr<vectorVolumes_Type> vectorVolumesPtr_Type;
78 
79  typedef std::vector<UInt> vectorIndexes_Type;
80  typedef std::shared_ptr<vectorIndexes_Type> vectorIndexesPtr_Type;
81 
82  typedef typename super::FESpacePtr_Type FESpacePtr_Type;
83  typedef typename super::ETFESpacePtr_Type ETFESpacePtr_Type;
84 
85  //Vector for vector parameters
86  typedef typename super::vectorsParameters_Type vectorsParameters_Type;
87  typedef typename super::vectorsParametersPtr_Type vectorsParametersPtr_Type;
88 
89  typedef MatrixSmall<3, 3> matrixSmall_Type;
90 
91  // Typedefs for expression definitions
92  typedef typename super::tensorF_Type tensorF_Type;
93  typedef typename super::determinantF_Type determinantF_Type;
94  typedef typename super::tensorC_Type tensorC_Type;
95  typedef typename super::minusT_Type minusT_Type;
96  typedef typename super::traceTensor_Type traceTensor_Type;
97  //@}
98 
99 
100 
101  //! @name Constructor & Destructor
102  //@{
103 
104  ExponentialMaterialNonLinear();
105 
106  virtual ~ExponentialMaterialNonLinear();
107 
108  //@}
109 
110 
111 
112  //!@name Methods
113  //@{
114 
115  //! Setup the created object of the class StructuralIsotropicConstitutiveLaw
116  /*!
117  \param dFespace: the FiniteElement Space
118  \param monolithicMap: the MapEpetra
119  \param offset: the offset parameter used assembling the matrices
120  */
121  void setup ( const FESpacePtr_Type& dFESpace,
122  const ETFESpacePtr_Type& dETFESpace,
123  const std::shared_ptr<const MapEpetra>& monolithicMap,
124  const UInt offset,const dataPtr_Type& dataMaterial);
125 
126  //! Compute the Stiffness matrix in StructuralSolver::buildSystem()
127  /*!
128  \param dataMaterial the class with Material properties data
129  */
130  void computeLinearStiff ( dataPtr_Type& /*dataMaterial*/,
131  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
132  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/);
133 
134 
135  //! Updates the Jacobian matrix in StructualSolver::updateJacobian
136  /*!
137  \param disp: solution at the k-th iteration of NonLinearRichardson Method
138  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get
139  the material coefficients (e.g. Young modulus, Poisson ratio..)
140  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
141  */
142  void updateJacobianMatrix ( const vector_Type& disp,
143  const dataPtr_Type& dataMaterial,
144  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
145  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
146  const displayerPtr_Type& displayer );
147 
148 
149  //! Updates the nonlinear terms in the Jacobian matrix in StructualSolver::updateJacobian
150  /*!
151  \param stiff: stiffness matrix provided from outside
152  \param disp: solution at the k-th iteration of NonLinearRichardson Method
153  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get
154  the material coefficients (e.g. Young modulus, Poisson ratio..)
155  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
156  */
157  void updateNonLinearJacobianTerms ( matrixPtr_Type& jacobian,
158  const vector_Type& disp,
159  const dataPtr_Type& dataMaterial,
160  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
161  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/,
162  const displayerPtr_Type& displayer );
163 
164 
165  //! Interface method to compute the new Stiffness matrix in StructuralSolver::evalResidual and in
166  //! StructuralSolver::updateSystem since the matrix is the expression of the matrix is the same.
167  /*!
168  \param sol: the solution vector
169  \param factor: scaling factor used in FSI
170  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get the
171  material coefficients (e.g. Young modulus, Poisson ratio..)
172  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
173  */
174  void computeStiffness ( const vector_Type& disp,
175  Real factor,
176  const dataPtr_Type& dataMaterial,
177  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
178  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/,
179  const displayerPtr_Type& displayer );
180 
181 
182  //! Computes the new Stiffness vector for Neo-Hookean and Exponential materials in StructuralSolver
183  //! given a certain displacement field.
184  //! This function is used both in StructuralSolver::evalResidual and in StructuralSolver::updateSystem
185  //! since the matrix is the expression of the matrix is the same.
186  /*!
187  \param sol: the solution vector
188  \param factor: scaling factor used in FSI
189  \param dataMaterial: a pointer to the dataType member in StructuralSolver class to get
190  the material coefficients (e.g. Young modulus, Poisson ratio..)
191  \param displayer: a pointer to the Dysplaier member in the StructuralSolver class
192  */
193  // void computeVector( const vector_Type& sol,
194  // Real factor,
195  // const dataPtr_Type& dataMaterial,
196  // const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
197  // const displayerPtr_Type& displayer );
198 
199  //! Computes the deformation gradient F, the cofactor matrix Cof(F),
200  //! the determinant of F (J = det(F)), the trace of right Cauchy-Green tensor tr(C)
201  //! This function is used in StructuralIsotropicConstitutiveLaw::computeStiffness
202  /*!
203  \param dk_loc: the elemental displacement
204  */
205  void computeKinematicsVariables ( const VectorElemental& dk_loc )
206  {
207  }
208 
209  //! Computes the deformation Gradient F, the cofactor of F Cof(F), the determinant of F J = det(F), the trace of C Tr(C).
210  /*!
211  \param dk_loc: local displacement vector
212  */
213  //void computeStress( const vector_Type& sol );
214 
215  //! ShowMe method of the class (saved on a file the stiffness vector and the jacobian)
216  void showMe ( std::string const& fileNameVectStiff,
217  std::string const& fileNameJacobain );
218 
219  //! Compute the First Piola Kirchhoff Tensor
220  /*!
221  \param firstPiola Epetra_SerialDenseMatrix that has to be filled
222  \param tensorF Epetra_SerialDenseMatrix the deformation gradient
223  \param cofactorF Epetra_SerialDenseMatrix cofactor of F
224  \param invariants std::vector with the invariants of C and the detF
225  \param material UInt number to get the material parameteres form the VenantElasticData class
226  */
227  void computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
228  const Epetra_SerialDenseMatrix& tensorF,
229  const Epetra_SerialDenseMatrix& cofactorF,
230  const std::vector<Real>& invariants,
231  const UInt marker);
232 
233  //! Compute the First Piola Kirchhoff Tensor
234  /*!
235  \param disp the displacement field from which we compute the fisrt piola-Kirchhoff tensor
236  \param sigma_1 the first column of the Cauchy stress tensor
237  \param sigma_2 the second column of the Cauchy stress tensor
238  \param sigma_3 the third column of the Cauchy stress tensor
239  */
240  void computeCauchyStressTensor ( const vectorPtr_Type disp,
241  const QuadratureRule& evalQuad,
242  vectorPtr_Type sigma_1,
243  vectorPtr_Type sigma_2,
244  vectorPtr_Type sigma_3);
245 
246  //@}
247 
248  //! @name Get Methods
249  //@{
250 
251  //! Get the Stiffness matrix
252  matrixPtr_Type const stiffMatrix() const
253  {
254  return super::M_jacobian;
255  }
256 
257 
258  //! Get the stiffness vector
259  vectorPtr_Type const stiffVector() const
260  {
261  return M_stiff;
262  }
263 
264  void apply ( const vector_Type& sol, vector_Type& res,
265  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
266  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
267  const displayerPtr_Type displayer) ;
268 
269  //@}
270 
271 
272 
273 protected:
274 
275  //! construct the vectors for the parameters
276  /*!
277  \param VOID
278  \return VOID
279  */
280  void setupVectorsParameters ( void );
281 
282 
283  //! Vector: stiffness non-linear
284  vectorPtr_Type M_stiff;
285 
286  //Create the indentity for F
287  matrixSmall_Type M_identity;
288 
289 };
290 
291 
292 
293 
294 
295 template <typename MeshType>
296 ExponentialMaterialNonLinear<MeshType>::ExponentialMaterialNonLinear() :
297  super ( ),
298  M_stiff ( ),
299  M_identity ( )
300 {
301 }
302 
303 
304 
305 
306 
307 template <typename MeshType>
308 ExponentialMaterialNonLinear<MeshType>::~ExponentialMaterialNonLinear()
309 {}
310 
311 
312 
313 
314 
315 template <typename MeshType>
316 void
317 ExponentialMaterialNonLinear<MeshType>::setup ( const FESpacePtr_Type& dFESpace,
318  const ETFESpacePtr_Type& dETFESpace,
319  const std::shared_ptr<const MapEpetra>& monolithicMap,
320  const UInt offset,
321  const dataPtr_Type& dataMaterial)
322 {
323 
324  this->M_dataMaterial = dataMaterial;
325  this->M_dispFESpace = dFESpace;
326  this->M_dispETFESpace = dETFESpace;
327  this->M_localMap = monolithicMap;
328  this->M_offset = offset;
329 
330  M_stiff.reset (new vector_Type (*this->M_localMap) );
331 
332  M_identity (0, 0) = 1.0;
333  M_identity (0, 1) = 0.0;
334  M_identity (0, 2) = 0.0;
335  M_identity (1, 0) = 0.0;
336  M_identity (1, 1) = 1.0;
337  M_identity (1, 2) = 0.0;
338  M_identity (2, 0) = 0.0;
339  M_identity (2, 1) = 0.0;
340  M_identity (2, 2) = 1.0;
341 
342  // The 3 is because the law uses three parameters (alpha, gamma, bulk).
343  // another way would be to set up the number of constitutive parameters of the law
344  // in the data file to get the right size. Note the comment below.
345  this->M_vectorsParameters.reset ( new vectorsParameters_Type ( 3 ) );
346 
347 
348  this->setupVectorsParameters( );
349 }
350 
351 
352 template <typename MeshType>
353 void
354 ExponentialMaterialNonLinear<MeshType>::setupVectorsParameters ( void )
355 {
356  // Paolo Tricerri: February, 20th
357  // In each class, the name of the parameters has to inserted in the law
358  // TODO: move the saving of the material parameters to more abstract objects
359  // such that in the class of the material we do not need to call explicitly
360  // the name of the parameter.
361 
362  // Number of volume on the local part of the mesh
363  UInt nbElements = this->M_dispFESpace->mesh()->numVolumes();
364 
365  // Parameter alpha
366  // 1. resize the vector in the first element of the vector.
367  (* (this->M_vectorsParameters) ) [0].resize ( nbElements );
368 
369  // Parameter gamma
370  (* (this->M_vectorsParameters) ) [1].resize ( nbElements );
371 
372  // Parameter bulk
373  (* (this->M_vectorsParameters) ) [2].resize ( nbElements );
374 
375  for (UInt i (0); i < nbElements; i++ )
376  {
377  // Extracting the marker
378  UInt markerID = this->M_dispFESpace->mesh()->element ( i ).markerID();
379 
380  Real alpha = this->M_dataMaterial->alpha ( markerID );
381  Real gamma = this->M_dataMaterial->gamma ( markerID );
382  Real bulk = this->M_dataMaterial->bulk ( markerID );
383 
384  ( (* (this->M_vectorsParameters) ) [0]) [ i ] = alpha;
385  ( (* (this->M_vectorsParameters) ) [1]) [ i ] = gamma;
386  ( (* (this->M_vectorsParameters) ) [2]) [ i ] = bulk;
387  }
388 }
389 
390 
391 template <typename MeshType>
392 void ExponentialMaterialNonLinear<MeshType>::computeLinearStiff (dataPtr_Type& /*dataMaterial*/,
393  const mapMarkerVolumesPtr_Type /*mapsMarkerVolumes*/,
394  const mapMarkerIndexesPtr_Type /*mapsMarkerIndexes*/)
395 {
396  //! Empty method for exponential material
397 }
398 
399 
400 
401 
402 
403 template <typename MeshType>
404 void ExponentialMaterialNonLinear<MeshType>::updateJacobianMatrix ( const vector_Type& disp,
405  const dataPtr_Type& dataMaterial,
406  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
407  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
408  const displayerPtr_Type& displayer )
409 {
410 
411  this->M_jacobian.reset (new matrix_Type (*this->M_localMap) );
412 
413  displayer->leaderPrint (" \n*********************************\n ");
414  updateNonLinearJacobianTerms (this->M_jacobian, disp, this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
415  displayer->leaderPrint (" \n*********************************\n ");
416 }
417 
418 template <typename MeshType>
419 void ExponentialMaterialNonLinear<MeshType>::updateNonLinearJacobianTerms ( matrixPtr_Type& jacobian,
420  const vector_Type& disp,
421  const dataPtr_Type& dataMaterial,
422  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
423  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
424  const displayerPtr_Type& displayer )
425 {
426  using namespace ExpressionAssembly;
427 
428  displayer->leaderPrint (" Non-Linear S- updating non linear terms in the Jacobian Matrix (Exponential)");
429 
430  * (jacobian) *= 0.0;
431 
432  //! Nonlinear part of jacobian
433  //! loop on volumes (i)
434 
435  // mapIterator_Type it;
436  // //mapIteratorIndex_Type itIndex;
437 
438  // vectorVolumesPtr_Type pointerListOfVolumes;
439  // vectorIndexesPtr_Type pointerListOfIndexes;
440 
441  // for( it = (*mapsMarkerVolumes).begin(); it != (*mapsMarkerVolumes).end(); it++ )
442  // {
443  //Given the marker pointed by the iterator, let's extract the material parameters
444  // UInt marker = it->first;
445 
446  // Debug
447  // UInt markerIndex = itIndex->first;
448  // ASSERT( marker == markerIndex, "The list of volumes is referring to a marker that is not the same as the marker of index!!!");
449 
450  // pointerListOfVolumes.reset( new vectorVolumes_Type(it->second) );
451  // pointerListOfIndexes.reset( new vectorIndexes_Type( (*mapsMarkerIndexes)[marker] ) );
452 
453  // Real bulk = dataMaterial->bulk(marker);
454  // Real alpha = dataMaterial->alpha(marker);
455  // Real gamma = dataMaterial->gamma(marker);
456 
457  // Definition of F
458  tensorF_Type F = ExpressionDefinitions::deformationGradient( this->M_dispETFESpace, disp, this->M_offset, this->M_identity );
459 
460  // Definition of J
461  determinantF_Type J = ExpressionDefinitions::determinantF( F );
462 
463  //Definition of C
464  tensorC_Type C = ExpressionDefinitions::tensorC( transpose(F), F );
465 
466  // Definition of F^-T
467  minusT_Type F_T = ExpressionDefinitions::minusT( F );
468 
469  // Definition I_C
470  traceTensor_Type I_C = ExpressionDefinitions::traceTensor( C );
471 
472  //Assembling Volumetric Part
473  integrate ( elements ( this->M_dispETFESpace->mesh() ),
474  this->M_dispFESpace->qr(),
475  this->M_dispETFESpace,
476  this->M_dispETFESpace,
477  value ( 1.0 / 2.0 ) * parameter ( (* (this->M_vectorsParameters) ) [2] ) * ( value (2.0) *pow (J, 2.0) - J + value (1.0) ) * dot ( F_T, grad (phi_j) ) * dot ( F_T, grad (phi_i) )
478  ) >> jacobian;
479 
480  integrate ( elements ( this->M_dispETFESpace->mesh() ),
481  this->M_dispFESpace->qr(),
482  this->M_dispETFESpace,
483  this->M_dispETFESpace,
484  value ( -1.0 / 2.0 ) * parameter ( (* (this->M_vectorsParameters) ) [2] ) * ( pow (J, 2.0) - J + log (J) ) * dot ( F_T * transpose (grad (phi_j) ) * F_T, grad (phi_i) )
485  ) >> jacobian;
486 
487  // //Assembling the Isochoric Part
488  // //! 1. Stiffness matrix : int { - 2/3 * alpha * J^(-2/3) * exp( gamma*( Ic_iso - 3) )* ( 1. + coefExp * Ic_iso )
489  // //! *( F^-T : \nabla \delta ) ( F : \nabla \v ) }
490  integrate ( elements ( this->M_dispETFESpace->mesh() ),
491  this->M_dispFESpace->qr(),
492  this->M_dispETFESpace,
493  this->M_dispETFESpace,
494  value (-2.0 / 3.0) * parameter ( (* (this->M_vectorsParameters) ) [0]) *
495  pow (J, -2.0 / 3.0) * exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow( J, -(2.0/3.0) ) * I_C - value (3.0) ) ) *
496  ( value (1.0) + parameter ( (* (this->M_vectorsParameters) ) [1]) * pow( J, -(2.0/3.0) ) * I_C ) *
497  dot ( F_T, grad (phi_j) ) * dot ( F, grad (phi_i) )
498  ) >> jacobian;
499 
500  //! 2. Stiffness matrix : int { 2 * alpha * gamma * J^(-4/3) * exp( gamma*( Ic_iso - 3) ) *
501  //! ( F : \nabla \delta ) ( F : \nabla \v ) }
502  integrate ( elements ( this->M_dispETFESpace->mesh() ),
503  this->M_dispFESpace->qr(),
504  this->M_dispETFESpace,
505  this->M_dispETFESpace,
506  value (2.0) * parameter ( (* (this->M_vectorsParameters) ) [0] ) * parameter ( (* (this->M_vectorsParameters) ) [1]) * pow (J, -4.0 / 3.0) * exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow( J, -(2.0/3.0) ) * I_C - value (3.0) ) ) *
507  dot ( F, grad (phi_j) ) * dot ( F, grad (phi_i) )
508  ) >> jacobian;
509 
510  // //! 3. Stiffness matrix : int { 2.0/9.0 * alpha * Ic_iso * exp( gamma*( Ic_iso - 3) )*
511  // //! ( 1. + gamma * Ic_iso )( F^-T : \nabla \delta ) ( F^-T : \nabla \v )}
512  integrate ( elements ( this->M_dispETFESpace->mesh() ) ,
513  this->M_dispFESpace->qr(),
514  this->M_dispETFESpace,
515  this->M_dispETFESpace,
516  value (2.0 / 9.0) * parameter ( (* (this->M_vectorsParameters) ) [0]) * ( pow( J, -(2.0/3.0) ) * I_C ) * exp ( parameter ( (* (this->M_vectorsParameters) ) [1]) * ( pow( J, (-2.0/3.0) ) * I_C - value (3.0) ) ) * ( value (1.0) + parameter ( (* (this->M_vectorsParameters) ) [1]) * pow( J, (-2.0/3.0) ) * I_C ) *
517  dot ( F_T, grad (phi_j) ) * dot ( F_T, grad (phi_i) )
518  ) >> jacobian;
519 
520  // //! 4. Stiffness matrix : int { -2.0/3.0 * alpha * J^(-2/3) * exp( gamma*( Ic_iso - 3) )
521  // //! * ( 1. + gamma * Ic_iso )( deformationGradientTensor : \nabla \delta ) ( F^-T : \nabla \v ) }
522  integrate ( elements ( this->M_dispETFESpace->mesh() ),
523  this->M_dispFESpace->qr(),
524  this->M_dispETFESpace,
525  this->M_dispETFESpace,
526  value (-2.0 / 3.0) * parameter ( (* (this->M_vectorsParameters) ) [0]) * pow (J, -2.0 / 3.0) * exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow( J, (-2.0/3.0) ) * I_C - value (3.0) ) ) *
527  ( value (1.0) + parameter ( (* (this->M_vectorsParameters) ) [1]) * pow( J, (-2.0/3.0) ) * I_C ) *
528  dot ( F, grad (phi_j) ) * dot ( F_T, grad (phi_i) )
529  ) >> jacobian;
530 
531 
532  // //! 5. Stiffness matrix : int { alpha * J^(-2/3) * exp( gamma*( Ic_iso - 3)) (\nabla \delta: \nabla \v)}
533  integrate ( elements ( this->M_dispETFESpace->mesh() ),
534  this->M_dispFESpace->qr(),
535  this->M_dispETFESpace,
536  this->M_dispETFESpace,
537  parameter ( (* (this->M_vectorsParameters) ) [0]) * pow (J, -2.0 / 3.0) *
538  exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow( J, (-2.0/3.0) ) * I_C - value (3.0) ) ) *
539  dot ( grad (phi_j), grad (phi_i) )
540  ) >> jacobian;
541 
542 
543  // //! 6. Stiffness matrix : int { 1.0/3.0 * alpha * Ic_iso * exp(gamma*( Ic_iso - 3)) *
544  // //! (F^-T [\nabla \delta]^t F^-T ) : \nabla \v }
545  integrate ( elements ( this->M_dispETFESpace->mesh() ),
546  this->M_dispFESpace->qr(),
547  this->M_dispETFESpace,
548  this->M_dispETFESpace,
549  value (1.0 / 3.0) * parameter ( (* (this->M_vectorsParameters) ) [0]) * pow(J, (-2.0/3.0) ) * I_C *
550  exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow(J, (-2.0/3.0) ) * I_C - value (3.0) ) ) *
551  dot ( F_T * transpose (grad (phi_j) ) * F_T , grad (phi_i) )
552  ) >> jacobian;
553 
554 
555  //}
556 
557  jacobian->globalAssemble();
558 }
559 
560 
561 
562 
563 
564 template <typename MeshType>
565 void ExponentialMaterialNonLinear<MeshType>::computeStiffness ( const vector_Type& disp,
566  Real /*factor*/,
567  const dataPtr_Type& dataMaterial,
568  const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
569  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,
570  const displayerPtr_Type& displayer )
571 {
572 
573  using namespace ExpressionAssembly;
574 
575  M_stiff.reset (new vector_Type (*this->M_localMap) );
576  * (M_stiff) *= 0.0;
577 
578  displayer->leaderPrint (" \n*********************************\n ");
579  displayer->leaderPrint (" Non-Linear S- Computing the Exponential nonlinear stiffness vector ");
580  displayer->leaderPrint (" \n*********************************\n ");
581 
582 
583  // This part of the code is to loop on the volumes with a certain ID. It is there to show the possibility
584  // but it is slower than the one actually implemented because we double the calls to the integrate functions
585  // mapIterator_Type it;
586  // //mapIteratorIndex_Type itIndex;
587  // vectorVolumesPtr_Type pointerListOfVolumes;
588  // vectorIndexesPtr_Type pointerListOfIndexes;
589 
590  // for( it = (*mapsMarkerVolumes).begin(); it != (*mapsMarkerVolumes).end(); it++ )
591  // {
592 
593  // //Given the marker pointed by the iterator, let's extract the material parameters
594  // UInt marker = it->first;
595 
596  // //Debug reason
597  // // UInt markerIndex = itIndex->first;
598  // // ASSERT( marker == markerIndex, "The list of volumes is referring to a marker that is not the same as the marker of index!!!");
599 
600  // pointerListOfVolumes.reset( new vectorVolumes_Type(it->second) );
601  // pointerListOfIndexes.reset( new vectorIndexes_Type( (*mapsMarkerIndexes)[marker] ) );
602 
603  // Real bulk = dataMaterial->bulk(marker);
604  // Real alpha = dataMaterial->alpha(marker);
605  // Real gamma = dataMaterial->gamma(marker);
606 
607 
608  // Here the only critical template parameter is MapEpetra. When other maps will be
609  // available in LifeV, in order to keep the same structure, we should template the
610  // StructuralConstitutiveLaw both on the MeshType and MapType.
611 
612  // Definition of F
613  tensorF_Type F = ExpressionDefinitions::deformationGradient( this->M_dispETFESpace, disp, this->M_offset, this->M_identity );
614 
615  // Definition of J
616  determinantF_Type J = ExpressionDefinitions::determinantF( F );
617 
618  //Definition of C
619  tensorC_Type C = ExpressionDefinitions::tensorC( transpose(F), F );
620 
621  // Definition of F^-T
622  minusT_Type F_T = ExpressionDefinitions::minusT( F );
623 
624  // Definition I_C
625  traceTensor_Type I_C = ExpressionDefinitions::traceTensor( C );
626 
627  //Computation of the volumetric part
628  integrate ( elements ( this->M_dispETFESpace->mesh() ) ,
629  this->M_dispFESpace->qr(),
630  this->M_dispETFESpace,
631  value (1.0 / 2.0) * parameter ( (* (this->M_vectorsParameters) ) [2] ) * ( pow ( J , 2.0) - J + log (J) ) * dot ( F_T, grad (phi_i) )
632  ) >> M_stiff;
633 
634  //Computation of the isochoric part
635  integrate ( elements ( this->M_dispETFESpace->mesh() ) ,
636  this->M_dispFESpace->qr(),
637  this->M_dispETFESpace,
638  parameter ( (* (this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) *
639  exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow (J, -2.0 / 3.0) * I_C - value (3.0) ) ) *
640  (dot ( F - value (1.0 / 3.0) * I_C * F_T, grad (phi_i) ) )
641  ) >> M_stiff;
642  //}
643 
644  this->M_stiff->globalAssemble();
645 }
646 
647 
648 template <typename MeshType>
649 void ExponentialMaterialNonLinear<MeshType>::showMe ( std::string const& fileNameStiff,
650  std::string const& fileNameJacobian )
651 {
652  this->M_stiff->spy (fileNameStiff);
654 }
655 
656 
657 template <typename MeshType>
658 void ExponentialMaterialNonLinear<MeshType>::apply ( const vector_Type& sol,
659  vector_Type& res, const mapMarkerVolumesPtr_Type mapsMarkerVolumes,
660  const mapMarkerIndexesPtr_Type mapsMarkerIndexes,const displayerPtr_Type displayer)
661 {
662  computeStiffness (sol, 0, this->M_dataMaterial, mapsMarkerVolumes, mapsMarkerIndexes, displayer);
663  res += *M_stiff;
664 }
665 
666 
667 template <typename MeshType>
668 void ExponentialMaterialNonLinear<MeshType>::computeLocalFirstPiolaKirchhoffTensor ( Epetra_SerialDenseMatrix& firstPiola,
669  const Epetra_SerialDenseMatrix& tensorF,
670  const Epetra_SerialDenseMatrix& cofactorF,
671  const std::vector<Real>& invariants,
672  const UInt marker)
673 {
674 
675  //Get the material parameters
676  Real alpha = this->M_dataMaterial->alpha (marker);
677  Real gamma = this->M_dataMaterial->gamma (marker);
678  Real bulk = this->M_dataMaterial->bulk (marker);
679 
680 
681  //Computing the first term \alphaJ^{-2/3}[F-(1/3)tr(C)F^{-T}]exp(\gamma(tr(Ciso) - 3)
682  Epetra_SerialDenseMatrix firstTerm (tensorF);
683  Epetra_SerialDenseMatrix copyCofactorF (cofactorF);
684 
685  Real scale (0.0);
686  scale = -invariants[0] / 3.0;
687  copyCofactorF.Scale ( scale );
688  firstTerm += copyCofactorF;
689 
690  //Computation trace of the isochoric C
691  Real trCiso (0.0);
692  trCiso = std::pow (invariants[3], - (2.0 / 3.0) ) * invariants[0];
693 
694  Real coef ( 0.0 );
695  coef = alpha * std::pow (invariants[3], - (2.0 / 3.0) ) * std::exp ( gamma * ( trCiso - 3 ) );
696  firstTerm.Scale ( coef );
697 
698  //Computing the second term (volumetric part) J*(bulk/2)(J-1+(1/J)*ln(J))F^{-T}
699  Epetra_SerialDenseMatrix secondTerm (cofactorF);
700  Real secCoef (0);
701  secCoef = invariants[3] * (bulk / 2.0) * (invariants[3] - 1 + (1.0 / invariants[3]) * std::log (invariants[3]) );
702 
703  secondTerm.Scale ( secCoef );
704 
705  firstPiola += firstTerm;
706  firstPiola += secondTerm;
707 
708 }
709 
710 template <typename MeshType>
711 void ExponentialMaterialNonLinear<MeshType>::computeCauchyStressTensor ( const vectorPtr_Type disp,
712  const QuadratureRule& evalQuad,
713  vectorPtr_Type sigma_1,
714  vectorPtr_Type sigma_2,
715  vectorPtr_Type sigma_3)
716 
717 {
718  using namespace ExpressionAssembly;
719 
720  // Definition that are needed
721  // Definition of F
722  tensorF_Type F = ExpressionDefinitions::deformationGradient( this->M_dispETFESpace, *disp, this->M_offset, this->M_identity );
723 
724  // Definition of J
725  determinantF_Type J = ExpressionDefinitions::determinantF( F );
726 
727  //Definition of C
728  tensorC_Type C = ExpressionDefinitions::tensorC( transpose(F), F );
729 
730  // Definition of F^-T
731  minusT_Type F_T = ExpressionDefinitions::minusT( F );
732 
733  // Definition I_C
734  traceTensor_Type I_C = ExpressionDefinitions::traceTensor( C );
735 
736  evaluateNode( elements ( this->M_dispETFESpace->mesh() ),
737  evalQuad,
738  this->M_dispETFESpace,
739  meas_K * dot ( vectorFromMatrix( ( 1 / J )*
740  ( value (1.0 / 2.0) * parameter ( (* (this->M_vectorsParameters) ) [2] ) * ( pow ( J , 2.0) - J + log (J) ) * F_T + // volumetric
741  parameter ( (* (this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) *
742  exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow (J, -2.0 / 3.0) * I_C - value (3.0) ) ) * ( F - value(1.0/3.0)* I_C* F_T ) ) *
743  transpose(F), 0 ), phi_i)
744  ) >> sigma_1;
745 
746  sigma_1->globalAssemble();
747  evaluateNode( elements ( this->M_dispETFESpace->mesh() ),
748  evalQuad,
749  this->M_dispETFESpace,
750  meas_K * dot ( vectorFromMatrix( ( 1 / J )*
751  ( value (1.0 / 2.0) * parameter ( (* (this->M_vectorsParameters) ) [2] ) * ( pow ( J , 2.0) - J + log (J) ) * F_T + // volumetric
752  parameter ( (* (this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) *
753  exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow (J, -2.0 / 3.0) * I_C - value (3.0) ) ) * ( F - value(1.0/3.0)* I_C* F_T ) ) *
754  transpose( F ) , 1 ), phi_i)
755  ) >> sigma_2;
756 
757  sigma_2->globalAssemble();
758 
759  evaluateNode( elements ( this->M_dispETFESpace->mesh() ),
760  evalQuad,
761  this->M_dispETFESpace,
762  meas_K * dot ( vectorFromMatrix( ( 1 / J )*
763  ( value (1.0 / 2.0) * parameter ( (* (this->M_vectorsParameters) ) [2] ) * ( pow ( J , 2.0) - J + log (J) ) * F_T + // volumetric
764  parameter ( (* (this->M_vectorsParameters) ) [0] ) * pow (J, -2.0 / 3.0) *
765  exp ( parameter ( (* (this->M_vectorsParameters) ) [1] ) * ( pow (J, -2.0 / 3.0) * I_C - value (3.0) ) ) * ( F - value(1.0/3.0)* I_C* F_T ) ) *
766  transpose( F ) , 2 ), phi_i)
767  ) >> sigma_3;
768 
769  sigma_3->globalAssemble();
770 
771 }
772 
773 
774 
775 template <typename MeshType>
776 inline StructuralIsotropicConstitutiveLaw<MeshType>* createExponentialMaterialNonLinear()
777 {
778  return new ExponentialMaterialNonLinear<MeshType >();
779 }
780 
781 namespace
782 {
784 }
785 
786 } //Namespace LifeV
787 
788 #endif /* __EXPONENTIALMATERIAL_H */
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)
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.
Definition: LifeV.hpp:175
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)
Definition: LifeV.hpp:191
ExpressionDphiI grad(const ExpressionPhiI &)
Simple function to be used in the construction of an expression.