LifeV
mainAnisotropic.cpp
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 by
12  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  \file main.cpp
28  \author Paolo Tricerri <paolo.tricerri@epfl.ch>
29  \date 2012-05-01
30 
31 
32  The general procedure to use this example is the following:
33  - read the set of solutions from the hdf5 file
34  - compute the sigma tensor using the reconstruction procedure
35  shown in testsuite/evaluateNodalETA/
36  - use the method StructuralOperator::computeEigenvalue to,
37  given a certain dof to compute the set of eigenvalues.
38  - post process the vector of eigenvalues
39  */
40 #undef HAVE_HDF5
41 #ifdef TWODIM
42 #error test_structure cannot be compiled in 2D
43 #endif
44 
45 // Tell the compiler to ignore specific kind of warnings:
46 #pragma GCC diagnostic ignored "-Wunused-variable"
47 #pragma GCC diagnostic ignored "-Wunused-parameter"
48 
49 #include <Epetra_ConfigDefs.h>
50 #ifdef EPETRA_MPI
51 #include <mpi.h>
52 #include <Epetra_MpiComm.h>
53 #else
54 #include <Epetra_SerialComm.h>
55 #endif
56 
57 //Tell the compiler to restore the warning previously silented
58 #pragma GCC diagnostic warning "-Wunused-variable"
59 #pragma GCC diagnostic warning "-Wunused-parameter"
60 
61 #include <lifev/core/LifeV.hpp>
62 
63 #include <lifev/core/array/MapEpetra.hpp>
64 
65 #include <lifev/core/mesh/MeshData.hpp>
66 #include <lifev/core/mesh/MeshPartitioner.hpp>
67 
68 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp>
69 #include <lifev/structure/solver/StructuralConstitutiveLaw.hpp>
70 #include <lifev/structure/solver/StructuralOperator.hpp>
71 
72 #include <lifev/structure/fem/ExpressionDefinitions.hpp>
73 
74 #include <lifev/structure/solver/isotropic/ExponentialMaterialNonLinear.hpp>
75 #include <lifev/structure/solver/isotropic/NeoHookeanMaterialNonLinear.hpp>
76 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp>
77 
78 // Evaluation operations
79 #include <lifev/core/array/MatrixSmall.hpp>
80 #include <lifev/eta/expression/Evaluate.hpp>
81 
82 #include <lifev/core/filter/ExporterEnsight.hpp>
83 #ifdef HAVE_HDF5
84 #include <lifev/core/filter/ExporterHDF5.hpp>
85 #endif
86 #include <lifev/core/filter/ExporterEmpty.hpp>
87 
88 #include <iostream>
89 #include "ud_functions.hpp"
90 
91 #include <functional>
92 
93 using namespace LifeV;
94 
95 int returnValue = EXIT_FAILURE;
96 
98 {
99  std::string stringList = list;
100  std::set<UInt> setList;
101  if ( list == "" )
102  {
103  return setList;
104  }
105  size_t commaPos = 0;
106  while ( commaPos != std::string::npos )
107  {
108  commaPos = stringList.find ( "," );
109  setList.insert ( atoi ( stringList.substr ( 0, commaPos ).c_str() ) );
110  stringList = stringList.substr ( commaPos + 1 );
111  }
112  setList.insert ( atoi ( stringList.c_str() ) );
113  return setList;
114 }
115 
116 class Structure
117 {
118 public:
119 
121 
122  // Filters
123  typedef LifeV::Exporter<mesh_Type > filter_Type;
125 
126  typedef LifeV::ExporterEmpty<mesh_Type > emptyFilter_Type;
128  typedef LifeV::ExporterEnsight<mesh_Type > ensightFilter_Type;
130 
131 #ifdef HAVE_HDF5
134 #endif
135 
136 
137 
138  /** @name Constructors, destructor
139  */
140  //@{
141  Structure ( int argc,
142  char** argv,
143  std::shared_ptr<Epetra_Comm> structComm );
144 
146  {}
147  //@}
148 
149  //@{
150  void run()
151  {
152  run3d();
153  }
154  //@}
155 
156 protected:
157 
158 private:
159 
160  /**
161  * run the 2D driven cylinder simulation
162  */
163  void run2d();
164 
165  /**
166  * run the 3D driven cylinder simulation
167  */
168  void run3d();
169  void checkResultHolzapfelModel (const Real tensNorm);
170  void resultCorrect ( void );
171 private:
172  struct Private;
174 };
175 
176 
177 
178 struct Structure::Private
179 {
181  rho (1),
182  poisson (1),
183  young (1),
184  bulk (1),
185  alpha (1),
186  gamma (1)
187  {}
188  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > fct_type;
189  double rho, poisson, young, bulk, alpha, gamma;
190 
192 
194 
195  static Real displacementHolzapfelModel (const Real& /*t*/, const Real& x, const Real& y, const Real& z, const ID& i)
196  {
197 
198  switch (i)
199  {
200  case 0:
201  return -0.025697538387539 * ( x - 0.5 );
202  break;
203  case 1:
204  return 0.083508698014296 / 2.0 * ( y );
205  break;
206  case 2:
207  return -0.013480042559669 * ( z + 0.5 );
208  break;
209  default:
210  ERROR_MSG ("This entry is not allowed: ud_functions.hpp");
211  return 0.;
212  break;
213  }
214 
215  }
216 
217 
218 };
219 
220 
221 
222 Structure::Structure ( int argc,
223  char** argv,
224  std::shared_ptr<Epetra_Comm> structComm) :
225  parameters ( new Private() )
226 {
227  GetPot command_line (argc, argv);
228  std::string data_file_name = command_line.follow ("data", 2, "-f", "--file");
229  GetPot dataFile ( data_file_name );
230  parameters->data_file_name = data_file_name;
231 
232  parameters->comm = structComm;
233  int ntasks = parameters->comm->NumProc();
234 
235  if (!parameters->comm->MyPID() )
236  {
237  std::cout << "My PID = " << parameters->comm->MyPID() << " out of " << ntasks << " running." << std::endl;
238  }
239 }
240 
241 
242 
243 void
244 Structure::run2d()
245 {
246  std::cout << "2D cylinder test case is not available yet\n";
247 }
248 
249 
250 
251 void
252 Structure::run3d()
253 {
254  typedef std::shared_ptr<VectorEpetra> vectorPtr_Type;
255 
256  // typedefs for fibers
257  // std function for fiber direction
258  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > fiberFunction_Type;
259  typedef std::shared_ptr<fiberFunction_Type> fiberFunctionPtr_Type;
260 
261  typedef std::vector<fiberFunctionPtr_Type> vectorFiberFunction_Type;
262  typedef std::shared_ptr<vectorFiberFunction_Type> vectorFiberFunctionPtr_Type;
263 
264  // General typedefs
265  typedef StructuralOperator<mesh_Type >::vector_Type vector_Type;
266  typedef std::shared_ptr<vector_Type> vectorPtr_Type;
267  typedef FESpace< mesh_Type, MapEpetra > solidFESpace_Type;
268  typedef std::shared_ptr<solidFESpace_Type> solidFESpacePtr_Type;
269 
270  typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 3 > solidETFESpace_Type;
271  typedef std::shared_ptr<solidETFESpace_Type> solidETFESpacePtr_Type;
272 
273  typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 1 > scalarETFESpace_Type;
274  typedef std::shared_ptr<scalarETFESpace_Type> scalarETFESpacePtr_Type;
275 
276  bool verbose = (parameters->comm->MyPID() == 0);
277 
278  //! BChandler use to create the StructuralOperator object
279  std::shared_ptr<BCHandler> BCh ( new BCHandler() );
280 
281  //! dataElasticStructure for parameters
282  GetPot dataFile ( parameters->data_file_name.c_str() );
283  std::shared_ptr<StructuralConstitutiveLawData> dataStructure (new StructuralConstitutiveLawData( ) );
284 
285  dataStructure->setup (dataFile);
286 
287  dataStructure->showMe();
288 
289  //! Read and partition mesh
290  MeshData meshData;
291  meshData.setup (dataFile, "solid/space_discretization");
292 
293  std::shared_ptr<mesh_Type > fullMeshPtr (new RegionMesh<LinearTetra> ( ( parameters->comm ) ) );
294  std::shared_ptr<mesh_Type > localMeshPtr (new RegionMesh<LinearTetra> ( ( parameters->comm ) ) );
295  readMesh (*fullMeshPtr, meshData);
296 
297  MeshPartitioner< mesh_Type > meshPart ( fullMeshPtr, parameters->comm );
298  localMeshPtr = meshPart.meshPartition();
299 
300  //! Functional spaces - needed for the computations of the gradients
301  std::string dOrder = dataFile ( "solid/space_discretization/order", "P1");
302  solidFESpacePtr_Type dFESpace ( new solidFESpace_Type (localMeshPtr, dOrder, 3, parameters->comm) );
303  solidETFESpacePtr_Type dETFESpace ( new solidETFESpace_Type (meshPart, & (dFESpace->refFE() ), & (dFESpace->fe().geoMap() ), parameters->comm) );
304 
305  //! Scalar ETFEspace to evaluate scalar quantities
306  solidFESpacePtr_Type dScalarFESpace ( new solidFESpace_Type (localMeshPtr, dOrder, 1, parameters->comm) );
307  scalarETFESpacePtr_Type dScalarETFESpace ( new scalarETFESpace_Type (meshPart, & (dFESpace->refFE() ), & (dFESpace->fe().geoMap() ), parameters->comm) );
308 
309  // Setting the fibers
310  vectorFiberFunctionPtr_Type pointerToVectorOfFamilies( new vectorFiberFunction_Type( ) );
311  (*pointerToVectorOfFamilies).resize( dataStructure->numberFibersFamilies( ) );
312 
313  if( verbose )
314  std::cout << "Number of families: " << (*pointerToVectorOfFamilies).size() << std::endl;
315 
316  fibersDirectionList setOfFiberFunctions;
317  setOfFiberFunctions.setupFiberDefinitions( dataStructure->numberFibersFamilies( ) );
318 
319  // Setting the vector of fibers functions
320  for( UInt k(1); k <= pointerToVectorOfFamilies->size( ); k++ )
321  {
322  // Setting up the name of the function to define the family
323  std::string family="Family";
324  // adding the number of the family
325  std::string familyNumber;
326  std::ostringstream number;
327  number << ( k );
328  familyNumber = number.str();
329 
330  // Name of the function to create
331  std::string creationString = family + familyNumber;
332  (*pointerToVectorOfFamilies)[ k-1 ].reset( new fiberFunction_Type() );
333  (*pointerToVectorOfFamilies)[ k-1 ] = setOfFiberFunctions.fiberDefinition( creationString );
334  }
335 
336  // Interpolate fibers
337  std::vector<vectorPtr_Type> vectorInterpolated(0);
338 
339  // Interpolating fiber fields
340  vectorInterpolated.resize( (*pointerToVectorOfFamilies).size() );
341 
342  for( UInt k(0); k < (*pointerToVectorOfFamilies).size(); k++ )
343  {
344  vectorInterpolated[ k ].reset( new vector_Type( dFESpace->map() ) );
345  dFESpace->interpolate ( *( ( *(pointerToVectorOfFamilies) )[ k ] ) ,
346  * ( ( vectorInterpolated )[ k ] ),
347  0.0 );
348  }
349 
350  // Creating the StructuralOperator object
351  //! 1. Constructor of the structuralSolver
352  StructuralOperator< RegionMesh<LinearTetra> > solid;
353 
354  //! 2. Setup of the structuralSolver
355  solid.setup (dataStructure,
356  dFESpace,
357  dETFESpace,
358  BCh,
359  parameters->comm);
360 
361  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
362  {
363  //! 3.b Setting the fibers in the abstract class of Anisotropic materials
364  solid.material()->anisotropicLaw()->setupFiberDirections( pointerToVectorOfFamilies );
365  }
366 
367 
368  // How many solution do we have to read?
369  UInt numberOfSol(0);
370  numberOfSol = dataFile.vector_variable_size ( "importer/iteration" );
371 
372  ASSERT( numberOfSol, "You did not set any solution to read!! ");
373 
374  //! 6. Post-processing setting
375  std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporter;
376 
377  std::string const exporterType = dataFile ( "exporter/type", "hdf5");
378  std::string const nameExporter = dataFile ( "exporter/name", "jacobian");
379 
380 #ifdef HAVE_HDF5
381  if (exporterType.compare ("hdf5") == 0)
382  {
383  exporter.reset ( new hdf5Filter_Type ( dataFile, nameExporter) );
384  }
385  else
386 #endif
387  {
388  if (exporterType.compare ("none") == 0)
389  {
390  exporter.reset ( new emptyFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
391  }
392 
393  else
394  {
395  exporter.reset ( new ensightFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
396  }
397  }
398  exporter->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
399 
400  vectorPtr_Type patchAreaVector;
401 
402  /*
403  Vector to pass to the StructuralOperator and to export
404  what is read.
405  */
406  vectorPtr_Type disp;
407  vectorPtr_Type dispRead;
408 
409  // Vectors to export Cauchy stress tensors
410  vectorPtr_Type sigma_1;
411  vectorPtr_Type sigma_2;
412  vectorPtr_Type sigma_3;
413 
414  // Vector for the eigenvalues
415  vectorPtr_Type vectorEigenvalues;
416 
417  patchAreaVector.reset ( new vector_Type ( dETFESpace->map() ) );
418 
419  disp.reset( new vector_Type( dFESpace->map() ) );
420  dispRead.reset( new vector_Type( dFESpace->map() ) );
421 
422  sigma_1.reset( new vector_Type( dFESpace->map() ) );
423  sigma_2.reset( new vector_Type( dFESpace->map() ) );
424  sigma_3.reset( new vector_Type( dFESpace->map() ) );
425 
426  vectorEigenvalues.reset( new vector_Type( dFESpace->map() ) );
427 
428  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "displacement",
429  dFESpace, dispRead, UInt (0) );
430 
431  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "sigma_1",
432  dFESpace, sigma_1, UInt (0) );
433 
434  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "sigma_2",
435  dFESpace, sigma_2, UInt (0) );
436 
437  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "sigma_3",
438  dFESpace, sigma_3, UInt (0) );
439 
440  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "vectorEigenvalues",
441  dFESpace, vectorEigenvalues, UInt (0) );
442 
443  exporter->postProcess ( 0.0 );
444 
445  //! =================================================================================
446  //! Analysis - Istant or Interval
447  //! =================================================================================
448 
449  MPI_Barrier (MPI_COMM_WORLD);
450 
451  QuadratureRule fakeQuadratureRule;
452 
453  Real refElemArea (0); //area of reference element
454  //compute the area of reference element
455  for (UInt iq = 0; iq < dFESpace->qr().nbQuadPt(); iq++)
456  {
457  refElemArea += dFESpace->qr().weight (iq);
458  }
459 
460  Real wQuad (refElemArea / dFESpace->refFE().nbDof() );
461 
462  //Setting the quadrature Points = DOFs of the element and weight = 1
463  std::vector<GeoVector> coords = dFESpace->refFE().refCoor();
464  std::vector<Real> weights (dFESpace->fe().nbFEDof(), wQuad);
465  fakeQuadratureRule.setDimensionShape ( shapeDimension (dFESpace->refFE().shape() ), dFESpace->refFE().shape() );
466  fakeQuadratureRule.setPoints (coords, weights);
467 
468 
469  //fakeQuadratureRule.showMe();
470  using namespace ExpressionAssembly;
471 
472  // Trying to compute the Jacobian using ET
473  MatrixSmall<3,3> identity;
474  identity (0, 0) = 1.0;
475  identity (0, 1) = 0.0;
476  identity (0, 2) = 0.0;
477  identity (1, 0) = 0.0;
478  identity (1, 1) = 1.0;
479  identity (1, 2) = 0.0;
480  identity (2, 0) = 0.0;
481  identity (2, 1) = 0.0;
482  identity (2, 2) = 1.0;
483 
485  evaluateNode( elements ( dETFESpace->mesh() ),
486  fakeQuadratureRule,
487  dETFESpace,
488  dot( vMeas , phi_i )
489  ) >> patchAreaVector;
490  patchAreaVector->globalAssemble();
491 
492  std::string const nameField = dataFile ( "importer/nameField", "displacement");
493 
494  *sigma_1 *=0.0;
495  *sigma_2 *=0.0;
496  *sigma_3 *=0.0;
497  *vectorEigenvalues *=0.0;
498 
499  // Reading the solution
500  // resetting the element of the vector
501  *disp *= 0.0;
502  *dispRead *= 0.0;
503 
504  // Interpolating the displacement
505  dFESpace->interpolate( static_cast<solidFESpace_Type::function_Type> ( Private::displacementHolzapfelModel ),
506  *disp, 0.0 );
507 
508  *dispRead = *disp;
509 
510  solid.computeCauchyStressTensor( disp, fakeQuadratureRule, sigma_1, sigma_2, sigma_3 );
511 
512  // Concluding reconstruction
513  *sigma_1 = *sigma_1 / *patchAreaVector;
514  *sigma_2 = *sigma_2 / *patchAreaVector;
515  *sigma_3 = *sigma_3 / *patchAreaVector;
516 
517  // Computing eigenvalues
518  solid.computePrincipalTensions( sigma_1, sigma_2, sigma_3, vectorEigenvalues );
519 
520  exporter->postProcess( 1.0 );
521 
522  // Check the result
523  checkResultHolzapfelModel ( vectorEigenvalues->norm2() );
524  // End of check
525 
526  std::cout.precision(16);
527  std::cout << "Norm 2 of the vector: " << vectorEigenvalues->norm2() << std::endl;
528 
529  if (verbose )
530  {
531  std::cout << "Analysis Completed!" << std::endl;
532  }
533 
534  //Closing files
535  exporter->closeFile( );
536 
537  MPI_Barrier (MPI_COMM_WORLD);
538  //!---------------------------------------------.-----------------------------------------------------
539 }
540 
542 {
543  if ( ( ( std::fabs (testETA - 7802999.6608) / 7802999.6608) <= 1e-5 ) )
544  {
545  this->resultCorrect( );
546  }
547 }
548 
549 void Structure::resultCorrect ( void )
550 {
551  std::cout << "Correct Result after the Analysis" << std::endl;
552  returnValue = EXIT_SUCCESS;
553 }
554 
555 int
556 main ( int argc, char** argv )
557 {
558 
559 #ifdef HAVE_MPI
560  MPI_Init (&argc, &argv);
561  std::shared_ptr<Epetra_MpiComm> Comm (new Epetra_MpiComm ( MPI_COMM_WORLD ) );
562  if ( Comm->MyPID() == 0 )
563  {
564  std::cout << "% using MPI" << std::endl;
565  }
566 #else
567  std::shared_ptr<Epetra_SerialComm> Comm ( new Epetra_SerialComm() );
568  std::cout << "% using serial Version" << std::endl;
569 #endif
570 
571  Structure structure ( argc, argv, Comm );
572  structure.run();
573 
574 #ifdef HAVE_MPI
575  MPI_Finalize();
576 #endif
577  return returnValue ;
578 }
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
Structure(int argc, char **argv, std::shared_ptr< Epetra_Comm > structComm)
void resultCorrect(void)
class ExpressionEmult Class for representing the transpose operation of an expression ...
static Real displacementHolzapfelModel(const Real &, const Real &x, const Real &y, const Real &z, const ID &i)
std::shared_ptr< emptyFilter_Type > emptyFilterPtr_Type
FESpace - Short description here please!
Definition: FESpace.hpp:78
int main(int argc, char **argv)
void checkResultHolzapfelModel(const Real tensNorm)
GetPot(const int argc_, char **argv_, const char *FieldSeparator=0x0)
Definition: GetPot.hpp:507
std::function< Real(Real const &, Real const &, Real const &, Real const &, ID const &) > fct_type
LifeV::ExporterEnsight< mesh_Type > ensightFilter_Type
LifeV::Exporter< mesh_Type > filter_Type
std::shared_ptr< hdf5Filter_Type > hdf5FilterPtr_Type
#define ERROR_MSG(A)
Definition: LifeAssert.hpp:69
LifeV::ExporterEmpty< mesh_Type > emptyFilter_Type
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
void run3d()
run the 3D driven cylinder simulation
#define ASSERT(X, A)
Definition: LifeAssert.hpp:90
uint32_type ID
IDs.
Definition: LifeV.hpp:194
LifeV::RegionMesh< LinearTetra > mesh_Type
std::shared_ptr< std::vector< Int > > M_isOnProc
void run2d()
run the 2D driven cylinder simulation
std::shared_ptr< ensightFilter_Type > ensightFilterPtr_Type
MeshData - class for handling spatial discretization.
Definition: MeshData.hpp:72
ExpressionMeas - Expression for the measure of the element.
double Real
Generic real data.
Definition: LifeV.hpp:175
std::set< UInt > parseList(const std::string &list)
std::shared_ptr< LifeV::Exporter< mesh_Type > > filterPtr_Type
QuadratureRule - The basis class for storing and accessing quadrature rules.
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
unsigned vector_variable_size(const char *VarName) const
Definition: GetPot.hpp:2291