42 #error test_structure cannot be compiled in 2D
46 #pragma GCC diagnostic ignored "-Wunused-variable" 47 #pragma GCC diagnostic ignored "-Wunused-parameter" 49 #include <Epetra_ConfigDefs.h> 52 #include <Epetra_MpiComm.h> 54 #include <Epetra_SerialComm.h> 58 #pragma GCC diagnostic warning "-Wunused-variable" 59 #pragma GCC diagnostic warning "-Wunused-parameter" 61 #include <lifev/core/LifeV.hpp> 63 #include <lifev/core/array/MapEpetra.hpp> 65 #include <lifev/core/mesh/MeshData.hpp> 66 #include <lifev/core/mesh/MeshPartitioner.hpp> 68 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp> 69 #include <lifev/structure/solver/StructuralConstitutiveLaw.hpp> 70 #include <lifev/structure/solver/StructuralOperator.hpp> 72 #include <lifev/structure/fem/ExpressionDefinitions.hpp> 74 #include <lifev/structure/solver/isotropic/ExponentialMaterialNonLinear.hpp> 75 #include <lifev/structure/solver/isotropic/NeoHookeanMaterialNonLinear.hpp> 76 #include <lifev/structure/solver/anisotropic/HolzapfelMaterialNonLinear.hpp> 79 #include <lifev/core/array/MatrixSmall.hpp> 80 #include <lifev/eta/expression/Evaluate.hpp> 82 #include <lifev/core/filter/ExporterEnsight.hpp> 84 #include <lifev/core/filter/ExporterHDF5.hpp> 86 #include <lifev/core/filter/ExporterEmpty.hpp> 93 using namespace LifeV;
99 std::string stringList = list;
100 std::set<UInt> setList;
106 while ( commaPos != std::string::npos )
108 commaPos = stringList.find (
"," );
109 setList.insert ( atoi ( stringList.substr ( 0, commaPos ).c_str() ) );
110 stringList = stringList.substr ( commaPos + 1 );
112 setList.insert ( atoi ( stringList.c_str() ) );
143 std::shared_ptr<Epetra_Comm> structComm );
201 return -0.025697538387539 * ( x - 0.5 );
204 return 0.083508698014296 / 2.0 * ( y );
207 return -0.013480042559669 * ( z + 0.5 );
210 ERROR_MSG (
"This entry is not allowed: ud_functions.hpp");
224 std::shared_ptr<Epetra_Comm> structComm) :
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;
232 parameters->comm = structComm;
233 int ntasks = parameters->comm->NumProc();
235 if (!parameters->comm->MyPID() )
237 std::cout <<
"My PID = " << parameters->comm->MyPID() <<
" out of " << ntasks <<
" running." << std::endl;
246 std::cout <<
"2D cylinder test case is not available yet\n";
254 typedef std::shared_ptr<VectorEpetra> vectorPtr_Type;
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;
261 typedef std::vector<fiberFunctionPtr_Type> vectorFiberFunction_Type;
262 typedef std::shared_ptr<vectorFiberFunction_Type> vectorFiberFunctionPtr_Type;
265 typedef StructuralOperator<mesh_Type >::vector_Type vector_Type;
266 typedef std::shared_ptr<vector_Type> vectorPtr_Type;
268 typedef std::shared_ptr<solidFESpace_Type> solidFESpacePtr_Type;
270 typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 3 > solidETFESpace_Type;
271 typedef std::shared_ptr<solidETFESpace_Type> solidETFESpacePtr_Type;
273 typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 1 > scalarETFESpace_Type;
274 typedef std::shared_ptr<scalarETFESpace_Type> scalarETFESpacePtr_Type;
276 bool verbose = (parameters->comm->MyPID() == 0);
279 std::shared_ptr<BCHandler> BCh (
new BCHandler() );
282 GetPot dataFile ( parameters->data_file_name.c_str() );
283 std::shared_ptr<StructuralConstitutiveLawData> dataStructure (
new StructuralConstitutiveLawData( ) );
285 dataStructure->setup (dataFile);
287 dataStructure->showMe();
291 meshData.setup (dataFile,
"solid/space_discretization");
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);
298 localMeshPtr = meshPart.meshPartition();
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) );
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) );
310 vectorFiberFunctionPtr_Type pointerToVectorOfFamilies(
new vectorFiberFunction_Type( ) );
311 (*pointerToVectorOfFamilies).resize( dataStructure->numberFibersFamilies( ) );
314 std::cout <<
"Number of families: " << (*pointerToVectorOfFamilies).size() << std::endl;
317 setOfFiberFunctions.setupFiberDefinitions( dataStructure->numberFibersFamilies( ) );
320 for( UInt k(1); k <= pointerToVectorOfFamilies->size( ); k++ )
323 std::string family=
"Family";
325 std::string familyNumber;
326 std::ostringstream number;
328 familyNumber = number.str();
331 std::string creationString = family + familyNumber;
332 (*pointerToVectorOfFamilies)[ k-1 ].reset(
new fiberFunction_Type() );
333 (*pointerToVectorOfFamilies)[ k-1 ] = setOfFiberFunctions.fiberDefinition( creationString );
337 std::vector<vectorPtr_Type> vectorInterpolated(0);
340 vectorInterpolated.resize( (*pointerToVectorOfFamilies).size() );
342 for( UInt k(0); k < (*pointerToVectorOfFamilies).size(); k++ )
344 vectorInterpolated[ k ].reset(
new vector_Type( dFESpace->map() ) );
345 dFESpace->interpolate ( *( ( *(pointerToVectorOfFamilies) )[ k ] ) ,
346 * ( ( vectorInterpolated )[ k ] ),
352 StructuralOperator< RegionMesh<LinearTetra> > solid;
355 solid.setup (dataStructure,
361 if( !dataStructure->constitutiveLaw().compare(
"anisotropic") )
364 solid.material()->anisotropicLaw()->setupFiberDirections( pointerToVectorOfFamilies );
372 ASSERT( numberOfSol,
"You did not set any solution to read!! ");
375 std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporter;
377 std::string
const exporterType = dataFile (
"exporter/type",
"hdf5");
378 std::string
const nameExporter = dataFile (
"exporter/name",
"jacobian");
381 if (exporterType.compare (
"hdf5") == 0)
383 exporter.reset (
new hdf5Filter_Type ( dataFile, nameExporter) );
388 if (exporterType.compare (
"none") == 0)
390 exporter.reset (
new emptyFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
395 exporter.reset (
new ensightFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
398 exporter->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
400 vectorPtr_Type patchAreaVector;
407 vectorPtr_Type dispRead;
410 vectorPtr_Type sigma_1;
411 vectorPtr_Type sigma_2;
412 vectorPtr_Type sigma_3;
415 vectorPtr_Type vectorEigenvalues;
417 patchAreaVector.reset (
new vector_Type ( dETFESpace->map() ) );
419 disp.reset(
new vector_Type( dFESpace->map() ) );
420 dispRead.reset(
new vector_Type( dFESpace->map() ) );
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() ) );
426 vectorEigenvalues.reset(
new vector_Type( dFESpace->map() ) );
428 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"displacement",
429 dFESpace, dispRead, UInt (0) );
431 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"sigma_1",
432 dFESpace, sigma_1, UInt (0) );
434 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"sigma_2",
435 dFESpace, sigma_2, UInt (0) );
437 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"sigma_3",
438 dFESpace, sigma_3, UInt (0) );
440 exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField,
"vectorEigenvalues",
441 dFESpace, vectorEigenvalues, UInt (0) );
443 exporter->postProcess ( 0.0 );
449 MPI_Barrier (MPI_COMM_WORLD);
453 Real refElemArea (0);
455 for (UInt iq = 0; iq < dFESpace->qr().nbQuadPt(); iq++)
457 refElemArea += dFESpace->qr().weight (iq);
460 Real wQuad (refElemArea / dFESpace->refFE().nbDof() );
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);
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;
485 evaluateNode( elements ( dETFESpace->mesh() ),
489 ) >> patchAreaVector;
490 patchAreaVector->globalAssemble();
492 std::string
const nameField = dataFile (
"importer/nameField",
"displacement");
497 *vectorEigenvalues *=0.0;
505 dFESpace->interpolate(
static_cast<solidFESpace_Type::function_Type> ( Private::displacementHolzapfelModel ),
510 solid.computeCauchyStressTensor( disp, fakeQuadratureRule, sigma_1, sigma_2, sigma_3 );
513 *sigma_1 = *sigma_1 / *patchAreaVector;
514 *sigma_2 = *sigma_2 / *patchAreaVector;
515 *sigma_3 = *sigma_3 / *patchAreaVector;
518 solid.computePrincipalTensions( sigma_1, sigma_2, sigma_3, vectorEigenvalues );
520 exporter->postProcess( 1.0 );
523 checkResultHolzapfelModel ( vectorEigenvalues->norm2() );
526 std::cout.precision(16);
527 std::cout <<
"Norm 2 of the vector: " << vectorEigenvalues->norm2() << std::endl;
531 std::cout <<
"Analysis Completed!" << std::endl;
535 exporter->closeFile( );
537 MPI_Barrier (MPI_COMM_WORLD);
543 if ( ( ( std::fabs (testETA - 7802999.6608) / 7802999.6608) <= 1e-5 ) )
551 std::cout <<
"Correct Result after the Analysis" << std::endl;
560 MPI_Init (&argc, &argv);
561 std::shared_ptr<Epetra_MpiComm> Comm (
new Epetra_MpiComm ( MPI_COMM_WORLD ) );
562 if ( Comm->MyPID() == 0 )
564 std::cout <<
"% using MPI" << std::endl;
567 std::shared_ptr<Epetra_SerialComm> Comm (
new Epetra_SerialComm() );
568 std::cout <<
"% using serial Version" << std::endl;
571 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 ...
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!
int main(int argc, char **argv)
void checkResultHolzapfelModel(const Real tensNorm)
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
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
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.
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)
unsigned vector_variable_size(const char *VarName) const