LifeV
lifev/structure/examples/example_computePrincipalTensions/main.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 #include <Epetra_ConfigDefs.h>
46 #ifdef EPETRA_MPI
47 #include <mpi.h>
48 #include <Epetra_MpiComm.h>
49 #else
50 #include <Epetra_SerialComm.h>
51 #endif
52 
53 #include <lifev/core/LifeV.hpp>
54 
55 #include <lifev/core/array/MapEpetra.hpp>
56 
57 #include <lifev/core/mesh/MeshData.hpp>
58 #include <lifev/core/mesh/MeshPartitioner.hpp>
59 
60 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp>
61 #include <lifev/structure/solver/StructuralConstitutiveLaw.hpp>
62 #include <lifev/structure/solver/StructuralOperator.hpp>
63 
64 #include <lifev/structure/solver/WallTensionEstimatorData.hpp>
65 #include <lifev/structure/solver/WallTensionEstimator.hpp>
66 
67 #include <lifev/structure/fem/ExpressionDefinitions.hpp>
68 
69 #include <lifev/structure/solver/isotropic/VenantKirchhoffMaterialNonLinearPenalized.hpp>
70 #include <lifev/structure/solver/isotropic/ExponentialMaterialNonLinear.hpp>
71 #include <lifev/structure/solver/isotropic/SecondOrderExponentialMaterialNonLinear.hpp>
72 
73 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp>
74 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp>
75 #include <lifev/structure/solver/anisotropic/HolzapfelGeneralizedMaterialNonLinear.hpp>
76 
77 // Evaluation operations
78 #include <lifev/core/array/MatrixSmall.hpp>
79 #include <lifev/eta/expression/Evaluate.hpp>
80 
81 #include <lifev/core/filter/ExporterEnsight.hpp>
82 #ifdef HAVE_HDF5
83 #include <lifev/core/filter/ExporterHDF5.hpp>
84 #endif
85 #include <lifev/core/filter/ExporterEmpty.hpp>
86 
87 #include <iostream>
88 #include "ud_functions.hpp"
89 
90 using namespace LifeV;
91 
92 int returnValue = EXIT_FAILURE;
93 
95 {
96  std::string stringList = list;
97  std::set<UInt> setList;
98  if ( list == "" )
99  {
100  return setList;
101  }
102  size_t commaPos = 0;
103  while ( commaPos != std::string::npos )
104  {
105  commaPos = stringList.find ( "," );
106  setList.insert ( atoi ( stringList.substr ( 0, commaPos ).c_str() ) );
107  stringList = stringList.substr ( commaPos + 1 );
108  }
109  setList.insert ( atoi ( stringList.c_str() ) );
110  return setList;
111 }
112 
113 class Structure
114 {
115 public:
116 
118 
119  // Filters
120  typedef LifeV::Exporter<mesh_Type > filter_Type;
122 
123  typedef LifeV::ExporterEmpty<mesh_Type > emptyFilter_Type;
125  typedef LifeV::ExporterEnsight<mesh_Type > ensightFilter_Type;
127 
128 #ifdef HAVE_HDF5
131 #endif
132 
133 
134 
135  /** @name Constructors, destructor
136  */
137  //@{
138  Structure ( int argc,
139  char** argv,
140  std::shared_ptr<Epetra_Comm> structComm );
141 
143  {}
144  //@}
145 
146  //@{
147  void run()
148  {
149  run3d();
150  }
151  //@}
152 
153 protected:
154 
155 private:
156 
157  /**
158  * run the 2D driven cylinder simulation
159  */
160  void run2d();
161 
162  /**
163  * run the 3D driven cylinder simulation
164  */
165  void run3d();
166 
167 private:
168  struct Private;
171 };
172 
173 
174 
175 struct Structure::Private
176 {
178  rho (1),
179  poisson (1),
180  young (1),
181  bulk (1),
182  alpha (1),
183  gamma (1)
184  {}
185  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > fct_type;
186  double rho, poisson, young, bulk, alpha, gamma;
187 
189 
191 
192 };
193 
194 
195 
196 Structure::Structure ( int argc,
197  char** argv,
198  std::shared_ptr<Epetra_Comm> structComm) :
199  parameters ( new Private() )
200 {
201  GetPot command_line (argc, argv);
202  std::string data_file_name = command_line.follow ("data", 2, "-f", "--file");
203  GetPot dataFile ( data_file_name );
204  parameters->data_file_name = data_file_name;
205 
206  parameters->comm = structComm;
207  int ntasks = parameters->comm->NumProc();
208 
209  if (!parameters->comm->MyPID() )
210  {
211  std::cout << "My PID = " << parameters->comm->MyPID() << " out of " << ntasks << " running." << std::endl;
212  }
213 }
214 
215 
216 
217 void
218 Structure::run2d()
219 {
220  std::cout << "2D cylinder test case is not available yet\n";
221 }
222 
223 
224 
225 void
226 Structure::run3d()
227 {
228  typedef std::shared_ptr<VectorEpetra> vectorPtr_Type;
229 
230  // typedefs for fibers
231  // std function for fiber direction
232  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > fiberFunction_Type;
233  typedef std::shared_ptr<fiberFunction_Type> fiberFunctionPtr_Type;
234 
235  typedef std::vector<fiberFunctionPtr_Type> vectorFiberFunction_Type;
236  typedef std::shared_ptr<vectorFiberFunction_Type> vectorFiberFunctionPtr_Type;
237 
238  typedef std::vector<vectorPtr_Type> listOfFiberDirections_Type;
239 
240  // General typedefs
241  typedef StructuralOperator<mesh_Type >::vector_Type vector_Type;
242  typedef std::shared_ptr<vector_Type> vectorPtr_Type;
243  typedef FESpace< mesh_Type, MapEpetra > solidFESpace_Type;
244  typedef std::shared_ptr<solidFESpace_Type> solidFESpacePtr_Type;
245 
246  typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 3 > solidETFESpace_Type;
247  typedef std::shared_ptr<solidETFESpace_Type> solidETFESpacePtr_Type;
248 
249  typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 1 > scalarETFESpace_Type;
250  typedef std::shared_ptr<scalarETFESpace_Type> scalarETFESpacePtr_Type;
251 
252  bool verbose = (parameters->comm->MyPID() == 0);
253 
254  //! BChandler use to create the StructuralOperator object
255  std::shared_ptr<BCHandler> BCh ( new BCHandler() );
256 
257  //! dataElasticStructure for parameters
258  GetPot dataFile ( parameters->data_file_name.c_str() );
259  std::shared_ptr<StructuralConstitutiveLawData> dataStructure (new StructuralConstitutiveLawData( ) );
260 
261  dataStructure->setup (dataFile);
262 
263  dataStructure->showMe();
264 
265  //! Parameters for the analysis
266  std::shared_ptr<WallTensionEstimatorData> tensionData (new WallTensionEstimatorData( ) );
267  tensionData->setup (dataFile);
268 
269  //! Read and partition mesh
270  MeshData meshData;
271  meshData.setup (dataFile, "solid/space_discretization");
272 
273  std::shared_ptr<mesh_Type > fullMeshPtr (new RegionMesh<LinearTetra> ( ( parameters->comm ) ) );
274  std::shared_ptr<mesh_Type > localMeshPtr (new RegionMesh<LinearTetra> ( ( parameters->comm ) ) );
275  readMesh (*fullMeshPtr, meshData);
276 
277  MeshPartitioner< mesh_Type > meshPart ( fullMeshPtr, parameters->comm );
278  localMeshPtr = meshPart.meshPartition();
279 
280  //! Functional spaces - needed for the computations of the gradients
281  std::string dOrder = dataFile ( "solid/space_discretization/order", "P1");
282  solidFESpacePtr_Type dFESpace ( new solidFESpace_Type (localMeshPtr, dOrder, 3, parameters->comm) );
283  solidETFESpacePtr_Type dETFESpace ( new solidETFESpace_Type (meshPart, & (dFESpace->refFE() ), & (dFESpace->fe().geoMap() ), parameters->comm) );
284 
285  //! Scalar ETFEspace to evaluate scalar quantities
286  solidFESpacePtr_Type dScalarFESpace ( new solidFESpace_Type (localMeshPtr, dOrder, 1, parameters->comm) );
287  scalarETFESpacePtr_Type dScalarETFESpace ( new scalarETFESpace_Type (meshPart, & (dFESpace->refFE() ), & (dFESpace->fe().geoMap() ), parameters->comm) );
288 
289  // Setting the fibers
290  vectorFiberFunctionPtr_Type pointerToVectorOfFamilies( new vectorFiberFunction_Type( ) );
291  (*pointerToVectorOfFamilies).resize( dataStructure->numberFibersFamilies( ) );
292 
293  if( verbose )
294  std::cout << "Number of families: " << (*pointerToVectorOfFamilies).size() << std::endl;
295 
296  fibersDirectionList setOfFiberFunctions;
297  setOfFiberFunctions.setupFiberDefinitions( dataStructure->numberFibersFamilies( ) );
298 
299  // Setting the vector of fibers functions
300  for( UInt k(1); k <= pointerToVectorOfFamilies->size( ); k++ )
301  {
302  // Setting up the name of the function to define the family
303  std::string family="Family";
304  // adding the number of the family
305  std::string familyNumber;
306  std::ostringstream number;
307  number << ( k );
308  familyNumber = number.str();
309 
310  // Name of the function to create
311  std::string creationString = family + familyNumber;
312  (*pointerToVectorOfFamilies)[ k-1 ].reset( new fiberFunction_Type() );
313  (*pointerToVectorOfFamilies)[ k-1 ] = setOfFiberFunctions.fiberDefinition( creationString );
314  }
315 
316  // Interpolate fibers
317  std::vector<vectorPtr_Type> vectorInterpolated(0);
318 
319  // Interpolating fiber fields
320  vectorInterpolated.resize( (*pointerToVectorOfFamilies).size() );
321 
322  for( UInt k(0); k < (*pointerToVectorOfFamilies).size(); k++ )
323  {
324  vectorInterpolated[ k ].reset( new vector_Type( dFESpace->map() ) );
325  dFESpace->interpolate ( *( ( *(pointerToVectorOfFamilies) )[ k ] ) ,
326  * ( ( vectorInterpolated )[ k ] ),
327  0.0 );
328  }
329 
330  // Creating the StructuralOperator object
331  //! 1. Constructor of the structuralSolver
332  StructuralOperator< RegionMesh<LinearTetra> > solid;
333 
334  //! 2. Setup of the structuralSolver
335  solid.setup (dataStructure,
336  dFESpace,
337  dETFESpace,
338  BCh,
339  parameters->comm);
340 
341  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
342  {
343  //! 3.b Setting the fibers in the abstract class of Anisotropic materials
344  solid.material()->anisotropicLaw()->setupFiberDirections( pointerToVectorOfFamilies );
345  }
346 
347 
348  //! 3. Creation of the importers to read the displacement field
349  std::string const filename = dataFile ( "importer/filename", "structure");
350  std::string const importerType = dataFile ( "importer/type", "hdf5");
351 
352  // How many solution do we have to read?
353  std::string readType = dataFile ( "importer/analysis", "instant");
354  UInt numberOfSol(0);
355  UInt start(0);
356  UInt end(0);
357 
358  if( !readType.compare("instant") )
359  {
360  numberOfSol = dataFile.vector_variable_size ( "importer/iteration" );
361  ASSERT( numberOfSol, "You did not set any solution to read!! ");
362  }
363  else
364  {
365  start = dataFile ( "importer/iterationStart" , 0 );
366  end = dataFile ( "importer/iterationEnd" , 0 );
367  numberOfSol = end - start;
368  }
369  ASSERT( numberOfSol, "You did not set any solution to read!! ");
370 
371 
372  if (verbose)
373  {
374  std::cout << "The filename is : " << filename << std::endl;
375  std::cout << "The importerType is: " << importerType << std::endl;
376  }
377 
378 #ifdef HAVE_HDF5
379  if (importerType.compare ("hdf5") == 0)
380  {
381  M_importer.reset ( new hdf5Filter_Type (dataFile, filename) );
382  }
383  else
384 #endif
385  {
386  if (importerType.compare ("none") == 0)
387  {
388  M_importer.reset ( new emptyFilter_Type ( dataFile, dFESpace->mesh(), "solid", dFESpace->map().comm().MyPID() ) );
389  }
390  else
391  {
392  M_importer.reset ( new ensightFilter_Type ( dataFile, filename ) );
393  }
394  }
395  M_importer->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
396 
397  //! 6. Post-processing setting
398  std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporter;
399 
400  std::string const exporterType = dataFile ( "exporter/type", "hdf5");
401  std::string const nameExporter = dataFile ( "exporter/name", "jacobian");
402 
403 #ifdef HAVE_HDF5
404  if (exporterType.compare ("hdf5") == 0)
405  {
406  exporter.reset ( new hdf5Filter_Type ( dataFile, nameExporter) );
407  }
408  else
409 #endif
410  {
411  if (exporterType.compare ("none") == 0)
412  {
413  exporter.reset ( new emptyFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
414  }
415 
416  else
417  {
418  exporter.reset ( new ensightFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
419  }
420  }
421  exporter->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
422 
423  // Patch area vector for the reconstruction
424  vectorPtr_Type patchAreaVectorScalar;
425  vectorPtr_Type patchAreaVector;
426 
427  /*
428  Vector to pass to the StructuralOperator and to export
429  what is read.
430  */
431  vectorPtr_Type disp;
432  vectorPtr_Type dispRead;
433 
434  // Vectors to export Cauchy stress tensors
435  vectorPtr_Type sigma_1;
436  vectorPtr_Type sigma_2;
437  vectorPtr_Type sigma_3;
438 
439  // Vector for the eigenvalues
440  vectorPtr_Type vectorEigenvalues;
441 
442  // Vector to save the stretch
443  std::vector< vectorPtr_Type > stretch(0);
444  stretch.resize( dataStructure->numberFibersFamilies( ) );
445 
446  patchAreaVector.reset ( new vector_Type ( dETFESpace->map() ) );
447  patchAreaVectorScalar.reset ( new vector_Type ( dScalarETFESpace->map() ) );
448 
449  disp.reset( new vector_Type( dFESpace->map() ) );
450  dispRead.reset( new vector_Type( dFESpace->map() ) );
451  sigma_1.reset( new vector_Type( dFESpace->map() ) );
452  sigma_2.reset( new vector_Type( dFESpace->map() ) );
453  sigma_3.reset( new vector_Type( dFESpace->map() ) );
454 
455  vectorEigenvalues.reset( new vector_Type( dFESpace->map() ) );
456 
457  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "displacement",
458  dFESpace, dispRead, UInt (0) );
459 
460  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "sigma_1",
461  dFESpace, sigma_1, UInt (0) );
462 
463  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "sigma_2",
464  dFESpace, sigma_2, UInt (0) );
465 
466  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "sigma_3",
467  dFESpace, sigma_3, UInt (0) );
468 
469  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "vectorEigenvalues",
470  dFESpace, vectorEigenvalues, UInt (0) );
471 
472 
473  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
474  {
475  for( UInt i(0); i < dataStructure->numberFibersFamilies( ); i++ )
476  {
477  std::string familyNumber;
478  std::ostringstream number;
479  number << ( i + 1 );
480  familyNumber = number.str();
481 
482  stretch[ i ].reset( new vector_Type( dScalarFESpace->map() ) );
483 
484  std::string stringNameS("stretch");
485 
486  stringNameS += "-"; stringNameS += familyNumber;
487 
488  exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::ScalarField, stringNameS,
489  dScalarFESpace, stretch[ i ], UInt (0) );
490  }
491  }
492 
493  exporter->postProcess ( 0.0 );
494 
495  //! =================================================================================
496  //! Analysis - Istant or Interval
497  //! =================================================================================
498 
499  MPI_Barrier (MPI_COMM_WORLD);
500 
501  QuadratureRule fakeQuadratureRule;
502 
503  Real refElemArea (0); //area of reference element
504  //compute the area of reference element
505  for (UInt iq = 0; iq < dFESpace->qr().nbQuadPt(); iq++)
506  {
507  refElemArea += dFESpace->qr().weight (iq);
508  }
509 
510  Real wQuad (refElemArea / dFESpace->refFE().nbDof() );
511 
512  //Setting the quadrature Points = DOFs of the element and weight = 1
513  std::vector<GeoVector> coords = dFESpace->refFE().refCoor();
514  std::vector<Real> weights (dFESpace->fe().nbFEDof(), wQuad);
515  fakeQuadratureRule.setDimensionShape ( shapeDimension (dFESpace->refFE().shape() ), dFESpace->refFE().shape() );
516  fakeQuadratureRule.setPoints (coords, weights);
517 
518  //fakeQuadratureRule.showMe();
519  using namespace ExpressionAssembly;
520 
521  // Trying to compute the Jacobian using ET
522  MatrixSmall<3,3> identity;
523  identity (0, 0) = 1.0;
524  identity (0, 1) = 0.0;
525  identity (0, 2) = 0.0;
526  identity (1, 0) = 0.0;
527  identity (1, 1) = 1.0;
528  identity (1, 2) = 0.0;
529  identity (2, 0) = 0.0;
530  identity (2, 1) = 0.0;
531  identity (2, 2) = 1.0;
532 
533  // Computing areas
534  evaluateNode( elements ( dScalarETFESpace->mesh() ),
535  fakeQuadratureRule,
536  dScalarETFESpace,
537  meas_K * phi_i
538  ) >> patchAreaVectorScalar;
539  patchAreaVectorScalar->globalAssemble();
540 
542  evaluateNode( elements ( dETFESpace->mesh() ),
543  fakeQuadratureRule,
544  dETFESpace,
545  dot( vMeas , phi_i )
546  ) >> patchAreaVector;
547  patchAreaVector->globalAssemble();
548 
549  std::string const nameField = dataFile ( "importer/nameField", "displacement");
550 
551  UInt i,k;
552 
553  for( i=start,k =0; i < numberOfSol; i++, k++ )
554  {
555  *sigma_1 *=0.0;
556  *sigma_2 *=0.0;
557  *sigma_3 *=0.0;
558  *vectorEigenvalues *=0.0;
559 
560  // Reading the solution
561  // resetting the element of the vector
562  *disp *= 0.0;
563  *dispRead *= 0.0;
564 
565  UInt current(0);
566  if( !readType.compare("interval") )
567  {
568  current = i ;
569  }
570  else
571  {
572  current = dataFile ( "importer/iteration" , 100000, i );
573  }
574  // Transform current ingot a string
575  std::string iterationString;
576  std::ostringstream number;
577  number.fill ( '0' );
578  number << std::setw (5) << ( current );
579  iterationString = number.str();
580 
581  std::cout << "Current reading: " << iterationString << std::endl;
582 
583  /*!Definition of the ExporterData, used to load the solution inside the previously defined vectors*/
584  LifeV::ExporterData<mesh_Type> solutionDispl (LifeV::ExporterData<mesh_Type>::VectorField,nameField + "." + iterationString,
585  dFESpace, disp, UInt (0), LifeV::ExporterData<mesh_Type>::UnsteadyRegime );
586 
587  //Read the variable
588  M_importer->readVariable (solutionDispl);
589 
590  *dispRead = *disp;
591 
592  solid.computeCauchyStressTensor( disp, fakeQuadratureRule, sigma_1, sigma_2, sigma_3 );
593 
594  // Concluding reconstruction
595  *sigma_1 = *sigma_1 / *patchAreaVector;
596  *sigma_2 = *sigma_2 / *patchAreaVector;
597  *sigma_3 = *sigma_3 / *patchAreaVector;
598 
599  // Computing eigenvalues
600  solid.computePrincipalTensions( sigma_1, sigma_2, sigma_3, vectorEigenvalues );
601 
602  if( !dataStructure->constitutiveLaw().compare("anisotropic") )
603  {
604  // Defining expressions
606  ExpressionDefinitions::deformationGradient ( dETFESpace, *disp, 0, identity );
607 
608  // Definition of C = F^T F
610  ExpressionDefinitions::tensorC( transpose(F), F );
611 
612  // Definition of J
614  ExpressionDefinitions::determinantF( F );
615 
616  // Definition of J^-(2/3) = det( C ) using the isochoric/volumetric decomposition
618  ExpressionDefinitions::powerExpression( J , (-2.0/3.0) );
619 
620  // Evaluating quantities that are related to fibers
621  for( UInt j(0); j < dataStructure->numberFibersFamilies( ); j++ )
622  {
623  // Fibers definitions
624  ExpressionDefinitions::interpolatedValue_Type fiberIth =
625  ExpressionDefinitions::interpolateFiber( dETFESpace, *(vectorInterpolated[ j ] ) );
626 
627  // f /otimes f
628  ExpressionDefinitions::outerProduct_Type Mith = ExpressionDefinitions::fiberTensor( fiberIth );
629 
630  // Definition of the fourth invariant : I_4^i = C:Mith
631  ExpressionDefinitions::stretch_Type IVith = ExpressionDefinitions::fiberStretch( C, Mith );
632 
633  // Definition of the fouth isochoric invariant : J^(-2.0/3.0) * I_4^i
634  ExpressionDefinitions::isochoricStretch_Type IVithBar =
635  ExpressionDefinitions::isochoricFourthInvariant( Jel, IVith );
636 
637  *stretch[ j ] *= 0.0;
638 
639  evaluateNode( elements ( dScalarETFESpace->mesh() ),
640  fakeQuadratureRule,
641  dScalarETFESpace,
642  meas_K * IVithBar * phi_i
643  ) >> stretch[ j ];
644  stretch[ j ]->globalAssemble();
645  *(stretch[ j ]) = *(stretch[ j ]) / *patchAreaVectorScalar;
646 
647  }
648  }
649  exporter->postProcess( dataStructure->dataTime()->initialTime() + k * dataStructure->dataTime()->timeStep() );
650  }
651 
652 
653  if (verbose )
654  {
655  std::cout << "Analysis Completed!" << std::endl;
656  }
657 
658  //Closing files
659  exporter->closeFile( );
660  M_importer->closeFile( );
661 
662  MPI_Barrier (MPI_COMM_WORLD);
663  //!---------------------------------------------.-----------------------------------------------------
664 }
665 
666 int
667 main ( int argc, char** argv )
668 {
669 
670 #ifdef HAVE_MPI
671  MPI_Init (&argc, &argv);
672  std::shared_ptr<Epetra_MpiComm> Comm (new Epetra_MpiComm ( MPI_COMM_WORLD ) );
673  if ( Comm->MyPID() == 0 )
674  {
675  std::cout << "% using MPI" << std::endl;
676  }
677 #else
678  std::shared_ptr<Epetra_SerialComm> Comm ( new Epetra_SerialComm() );
679  std::cout << "% using serial Version" << std::endl;
680 #endif
681 
682  Structure structure ( argc, argv, Comm );
683  structure.run();
684 
685 #ifdef HAVE_MPI
686  MPI_Finalize();
687 #endif
688  return returnValue ;
689 }
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)
class ExpressionEmult Class for representing the transpose operation of an expression ...
std::shared_ptr< emptyFilter_Type > emptyFilterPtr_Type
FESpace - Short description here please!
Definition: FESpace.hpp:78
int main(int argc, char **argv)
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
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
LifeV::RegionMesh< LinearTetra > mesh_Type
ExpressionProduct< ExpressionTranspose< ExpressionAddition< ExpressionInterpolateGradient< MeshType, MapEpetra, 3, 3 >, ExpressionMatrix< 3, 3 > > >, ExpressionAddition< ExpressionInterpolateGradient< MeshType, MapEpetra, 3, 3 >, ExpressionMatrix< 3, 3 > > > rightCauchyGreenTensor_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
ExpressionAddition< ExpressionInterpolateGradient< MeshType, MapEpetra, 3, 3 >, ExpressionMatrix< 3, 3 > > deformationGradient_Type
std::set< UInt > parseList(const std::string &list)
int operator()(const char *VarName, int Default, unsigned Idx) const
Definition: GetPot.hpp:2137
int operator()(const char *VarName, int Default) const
Definition: GetPot.hpp:2021
std::shared_ptr< LifeV::Exporter< mesh_Type > > filterPtr_Type
ExpressionDeterminant< ExpressionAddition< ExpressionInterpolateGradient< MeshType, MapEpetra, 3, 3 >, ExpressionMatrix< 3, 3 > > > determinantTensorF_Type
QuadratureRule - The basis class for storing and accessing quadrature rules.
ExpressionPower< ExpressionDeterminant< ExpressionAddition< ExpressionInterpolateGradient< MeshType, MapEpetra, 3, 3 >, ExpressionMatrix< 3, 3 > > > > powerExpression_Type
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