42 #error test_structure cannot be compiled in 2D
45 #include <Epetra_ConfigDefs.h> 48 #include <Epetra_MpiComm.h> 50 #include <Epetra_SerialComm.h> 53 #include <lifev/core/LifeV.hpp> 55 #include <lifev/core/array/MapEpetra.hpp> 57 #include <lifev/core/mesh/MeshData.hpp> 58 #include <lifev/core/mesh/MeshPartitioner.hpp> 60 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp> 61 #include <lifev/structure/solver/StructuralConstitutiveLaw.hpp> 62 #include <lifev/structure/solver/StructuralOperator.hpp> 64 #include <lifev/structure/solver/WallTensionEstimatorData.hpp> 65 #include <lifev/structure/solver/WallTensionEstimator.hpp> 67 #include <lifev/structure/fem/ExpressionDefinitions.hpp> 69 #include <lifev/structure/solver/isotropic/VenantKirchhoffMaterialNonLinearPenalized.hpp> 70 #include <lifev/structure/solver/isotropic/ExponentialMaterialNonLinear.hpp> 71 #include <lifev/structure/solver/isotropic/SecondOrderExponentialMaterialNonLinear.hpp> 73 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp> 74 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp> 75 #include <lifev/structure/solver/anisotropic/HolzapfelGeneralizedMaterialNonLinear.hpp> 78 #include <lifev/core/array/MatrixSmall.hpp> 79 #include <lifev/eta/expression/Evaluate.hpp> 81 #include <lifev/core/filter/ExporterEnsight.hpp> 83 #include <lifev/core/filter/ExporterHDF5.hpp> 85 #include <lifev/core/filter/ExporterEmpty.hpp> 90 using namespace LifeV;
96 std::string stringList = list;
97 std::set<UInt> setList;
103 while ( commaPos != std::string::npos )
105 commaPos = stringList.find (
"," );
106 setList.insert ( atoi ( stringList.substr ( 0, commaPos ).c_str() ) );
107 stringList = stringList.substr ( commaPos + 1 );
109 setList.insert ( atoi ( stringList.c_str() ) );
140 std::shared_ptr<Epetra_Comm> structComm );
198 std::shared_ptr<Epetra_Comm> structComm) :
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;
206 parameters->comm = structComm;
207 int ntasks = parameters->comm->NumProc();
209 if (!parameters->comm->MyPID() )
211 std::cout <<
"My PID = " << parameters->comm->MyPID() <<
" out of " << ntasks <<
" running." << std::endl;
220 std::cout <<
"2D cylinder test case is not available yet\n";
228 typedef std::shared_ptr<VectorEpetra> vectorPtr_Type;
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;
235 typedef std::vector<fiberFunctionPtr_Type> vectorFiberFunction_Type;
236 typedef std::shared_ptr<vectorFiberFunction_Type> vectorFiberFunctionPtr_Type;
238 typedef std::vector<vectorPtr_Type> listOfFiberDirections_Type;
241 typedef StructuralOperator<mesh_Type >::vector_Type vector_Type;
242 typedef std::shared_ptr<vector_Type> vectorPtr_Type;
244 typedef std::shared_ptr<solidFESpace_Type> solidFESpacePtr_Type;
246 typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 3 > solidETFESpace_Type;
247 typedef std::shared_ptr<solidETFESpace_Type> solidETFESpacePtr_Type;
249 typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 1 > scalarETFESpace_Type;
250 typedef std::shared_ptr<scalarETFESpace_Type> scalarETFESpacePtr_Type;
252 bool verbose = (parameters->comm->MyPID() == 0);
255 std::shared_ptr<BCHandler> BCh (
new BCHandler() );
258 GetPot dataFile ( parameters->data_file_name.c_str() );
259 std::shared_ptr<StructuralConstitutiveLawData> dataStructure (
new StructuralConstitutiveLawData( ) );
261 dataStructure->setup (dataFile);
263 dataStructure->showMe();
266 std::shared_ptr<WallTensionEstimatorData> tensionData (
new WallTensionEstimatorData( ) );
267 tensionData->setup (dataFile);
271 meshData.setup (dataFile,
"solid/space_discretization");
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);
278 localMeshPtr = meshPart.meshPartition();
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) );
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) );
290 vectorFiberFunctionPtr_Type pointerToVectorOfFamilies(
new vectorFiberFunction_Type( ) );
291 (*pointerToVectorOfFamilies).resize( dataStructure->numberFibersFamilies( ) );
294 std::cout <<
"Number of families: " << (*pointerToVectorOfFamilies).size() << std::endl;
297 setOfFiberFunctions.setupFiberDefinitions( dataStructure->numberFibersFamilies( ) );
300 for( UInt k(1); k <= pointerToVectorOfFamilies->size( ); k++ )
303 std::string family=
"Family";
305 std::string familyNumber;
306 std::ostringstream number;
308 familyNumber = number.str();
311 std::string creationString = family + familyNumber;
312 (*pointerToVectorOfFamilies)[ k-1 ].reset(
new fiberFunction_Type() );
313 (*pointerToVectorOfFamilies)[ k-1 ] = setOfFiberFunctions.fiberDefinition( creationString );
317 std::vector<vectorPtr_Type> vectorInterpolated(0);
320 vectorInterpolated.resize( (*pointerToVectorOfFamilies).size() );
322 for( UInt k(0); k < (*pointerToVectorOfFamilies).size(); k++ )
324 vectorInterpolated[ k ].reset(
new vector_Type( dFESpace->map() ) );
325 dFESpace->interpolate ( *( ( *(pointerToVectorOfFamilies) )[ k ] ) ,
326 * ( ( vectorInterpolated )[ k ] ),
332 StructuralOperator< RegionMesh<LinearTetra> > solid;
335 solid.setup (dataStructure,
341 if( !dataStructure->constitutiveLaw().compare(
"anisotropic") )
344 solid.material()->anisotropicLaw()->setupFiberDirections( pointerToVectorOfFamilies );
349 std::string
const filename = dataFile (
"importer/filename",
"structure");
350 std::string
const importerType = dataFile (
"importer/type",
"hdf5");
353 std::string readType = dataFile (
"importer/analysis",
"instant");
358 if( !readType.compare(
"instant") )
361 ASSERT( numberOfSol,
"You did not set any solution to read!! ");
365 start = dataFile
( "importer/iterationStart" , 0
);
366 end = dataFile
( "importer/iterationEnd" , 0
);
367 numberOfSol = end - start;
369 ASSERT( numberOfSol,
"You did not set any solution to read!! ");
374 std::cout <<
"The filename is : " << filename << std::endl;
375 std::cout <<
"The importerType is: " << importerType << std::endl;
379 if (importerType.compare (
"hdf5") == 0)
381 M_importer.reset (
new hdf5Filter_Type (dataFile, filename) );
386 if (importerType.compare (
"none") == 0)
388 M_importer.reset (
new emptyFilter_Type ( dataFile, dFESpace->mesh(),
"solid", dFESpace->map().comm().MyPID() ) );
392 M_importer.reset (
new ensightFilter_Type ( dataFile, filename ) );
395 M_importer->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
398 std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporter;
400 std::string
const exporterType = dataFile (
"exporter/type",
"hdf5");
401 std::string
const nameExporter = dataFile (
"exporter/name",
"jacobian");
404 if (exporterType.compare (
"hdf5") == 0)
406 exporter.reset (
new hdf5Filter_Type ( dataFile, nameExporter) );
411 if (exporterType.compare (
"none") == 0)
413 exporter.reset (
new emptyFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
418 exporter.reset (
new ensightFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
421 exporter->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
424 vectorPtr_Type patchAreaVectorScalar;
425 vectorPtr_Type patchAreaVector;
432 vectorPtr_Type dispRead;
435 vectorPtr_Type sigma_1;
436 vectorPtr_Type sigma_2;
437 vectorPtr_Type sigma_3;
440 vectorPtr_Type vectorEigenvalues;
443 std::vector< vectorPtr_Type > stretch(0);
444 stretch.resize( dataStructure->numberFibersFamilies( ) );
446 patchAreaVector.reset (
new vector_Type ( dETFESpace->map() ) );
447 patchAreaVectorScalar.reset (
new vector_Type ( dScalarETFESpace->map() ) );
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() ) );
455 vectorEigenvalues.reset(
new vector_Type( dFESpace->map() ) );
457 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"displacement",
458 dFESpace, dispRead, UInt (0) );
460 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"sigma_1",
461 dFESpace, sigma_1, UInt (0) );
463 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"sigma_2",
464 dFESpace, sigma_2, UInt (0) );
466 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"sigma_3",
467 dFESpace, sigma_3, UInt (0) );
469 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"vectorEigenvalues",
470 dFESpace, vectorEigenvalues, UInt (0) );
473 if( !dataStructure->constitutiveLaw().compare(
"anisotropic") )
475 for( UInt i(0); i < dataStructure->numberFibersFamilies( ); i++ )
477 std::string familyNumber;
478 std::ostringstream number;
480 familyNumber = number.str();
482 stretch[ i ].reset(
new vector_Type( dScalarFESpace->map() ) );
484 std::string stringNameS(
"stretch");
486 stringNameS +=
"-"; stringNameS += familyNumber;
488 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::ScalarField, stringNameS,
489 dScalarFESpace, stretch[ i ], UInt (0) );
493 exporter->postProcess ( 0.0 );
499 MPI_Barrier (MPI_COMM_WORLD);
503 Real refElemArea (0);
505 for (UInt iq = 0; iq < dFESpace->qr().nbQuadPt(); iq++)
507 refElemArea += dFESpace->qr().weight (iq);
510 Real wQuad (refElemArea / dFESpace->refFE().nbDof() );
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);
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;
534 evaluateNode( elements ( dScalarETFESpace->mesh() ),
538 ) >> patchAreaVectorScalar;
539 patchAreaVectorScalar->globalAssemble();
542 evaluateNode( elements ( dETFESpace->mesh() ),
546 ) >> patchAreaVector;
547 patchAreaVector->globalAssemble();
549 std::string
const nameField = dataFile (
"importer/nameField",
"displacement");
553 for( i=start,k =0; i < numberOfSol; i++, k++ )
558 *vectorEigenvalues *=0.0;
566 if( !readType.compare(
"interval") )
572 current = dataFile
( "importer/iteration" , 100000
, i );
575 std::string iterationString;
576 std::ostringstream number;
578 number << std::setw (5) << ( current );
579 iterationString = number.str();
581 std::cout <<
"Current reading: " << iterationString << std::endl;
584 LifeV::ExporterData<
mesh_Type> solutionDispl (LifeV::ExporterData<mesh_Type>::VectorField,nameField +
"." + iterationString,
585 dFESpace, disp, UInt (0), LifeV::ExporterData<mesh_Type>::UnsteadyRegime );
588 M_importer->readVariable (solutionDispl);
592 solid.computeCauchyStressTensor( disp, fakeQuadratureRule, sigma_1, sigma_2, sigma_3 );
595 *sigma_1 = *sigma_1 / *patchAreaVector;
596 *sigma_2 = *sigma_2 / *patchAreaVector;
597 *sigma_3 = *sigma_3 / *patchAreaVector;
600 solid.computePrincipalTensions( sigma_1, sigma_2, sigma_3, vectorEigenvalues );
602 if( !dataStructure->constitutiveLaw().compare(
"anisotropic") )
606 ExpressionDefinitions::deformationGradient ( dETFESpace, *disp, 0, identity );
610 ExpressionDefinitions::tensorC( transpose(F), F );
614 ExpressionDefinitions::determinantF( F );
618 ExpressionDefinitions::powerExpression( J , (-2.0/3.0) );
621 for( UInt j(0); j < dataStructure->numberFibersFamilies( ); j++ )
624 ExpressionDefinitions::interpolatedValue_Type fiberIth =
625 ExpressionDefinitions::interpolateFiber( dETFESpace, *(vectorInterpolated[ j ] ) );
628 ExpressionDefinitions::outerProduct_Type Mith = ExpressionDefinitions::fiberTensor( fiberIth );
631 ExpressionDefinitions::stretch_Type IVith = ExpressionDefinitions::fiberStretch( C, Mith );
634 ExpressionDefinitions::isochoricStretch_Type IVithBar =
635 ExpressionDefinitions::isochoricFourthInvariant( Jel, IVith );
637 *stretch[ j ] *= 0.0;
639 evaluateNode( elements ( dScalarETFESpace->mesh() ),
642 meas_K * IVithBar * phi_i
644 stretch[ j ]->globalAssemble();
645 *(stretch[ j ]) = *(stretch[ j ]) / *patchAreaVectorScalar;
649 exporter->postProcess( dataStructure->dataTime()->initialTime() + k * dataStructure->dataTime()->timeStep() );
655 std::cout <<
"Analysis Completed!" << std::endl;
659 exporter->closeFile( );
660 M_importer->closeFile( );
662 MPI_Barrier (MPI_COMM_WORLD);
671 MPI_Init (&argc, &argv);
672 std::shared_ptr<Epetra_MpiComm> Comm (
new Epetra_MpiComm ( MPI_COMM_WORLD ) );
673 if ( Comm->MyPID() == 0 )
675 std::cout <<
"% using MPI" << std::endl;
678 std::shared_ptr<Epetra_SerialComm> Comm (
new Epetra_SerialComm() );
679 std::cout <<
"% using serial Version" << std::endl;
682 Structure structure ( argc, argv, Comm );
std::shared_ptr< Epetra_Comm > comm
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!
int main(int argc, char **argv)
GetPot(const int argc_, char **argv_, const char *FieldSeparator=0x0)
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
filterPtr_Type M_importer
LifeV::ExporterEmpty< mesh_Type > emptyFilter_Type
Epetra_Import const & importer()
Getter for the Epetra_Import.
void run3d()
run the 3D driven cylinder simulation
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::string data_file_name
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.
std::shared_ptr< Private > parameters
ExpressionMeas - Expression for the measure of the element.
double Real
Generic real data.
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
int operator()(const char *VarName, int Default) const
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)
unsigned vector_variable_size(const char *VarName) const