LifeV
lifev/structure/examples/example_computingJacobian/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 #undef HAVE_HDF5
32 #ifdef TWODIM
33 #error test_structure cannot be compiled in 2D
34 #endif
35 
36 // Tell the compiler to ignore specific kind of warnings:
37 #pragma GCC diagnostic ignored "-Wunused-variable"
38 #pragma GCC diagnostic ignored "-Wunused-parameter"
39 
40 #include <Epetra_ConfigDefs.h>
41 #ifdef EPETRA_MPI
42 #include <mpi.h>
43 #include <Epetra_MpiComm.h>
44 #else
45 #include <Epetra_SerialComm.h>
46 #endif
47 
48 //Tell the compiler to restore the warning previously silented
49 #pragma GCC diagnostic warning "-Wunused-variable"
50 #pragma GCC diagnostic warning "-Wunused-parameter"
51 
52 #include <lifev/core/LifeV.hpp>
53 
54 #include <lifev/core/array/MapEpetra.hpp>
55 
56 #include <lifev/core/mesh/MeshData.hpp>
57 #include <lifev/core/mesh/MeshPartitioner.hpp>
58 
59 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp>
60 #include <lifev/structure/solver/StructuralConstitutiveLaw.hpp>
61 #include <lifev/structure/solver/StructuralOperator.hpp>
62 
63 //Not Used
64 #include <lifev/structure/solver/VenantKirchhoffMaterialLinear.hpp>
65 #include <lifev/structure/solver/VenantKirchhoffMaterialNonLinear.hpp>
66 #include <lifev/structure/solver/NeoHookeanMaterialNonLinear.hpp>
67 #include <lifev/structure/solver/ExponentialMaterialNonLinear.hpp>
68 
69 #include <lifev/core/filter/ExporterEnsight.hpp>
70 #ifdef HAVE_HDF5
71 #include <lifev/core/filter/ExporterHDF5.hpp>
72 #endif
73 #include <lifev/core/filter/ExporterEmpty.hpp>
74 
75 #include <iostream>
76 
77 
78 using namespace LifeV;
79 
80 int returnValue = EXIT_FAILURE;
81 
83 {
84  std::string stringList = list;
85  std::set<UInt> setList;
86  if ( list == "" )
87  {
88  return setList;
89  }
90  size_t commaPos = 0;
91  while ( commaPos != std::string::npos )
92  {
93  commaPos = stringList.find ( "," );
94  setList.insert ( atoi ( stringList.substr ( 0, commaPos ).c_str() ) );
95  stringList = stringList.substr ( commaPos + 1 );
96  }
97  setList.insert ( atoi ( stringList.c_str() ) );
98  return setList;
99 }
100 
101 
102 class Structure
103 {
104 public:
105 
107 
108  // Filters
109  typedef LifeV::Exporter<mesh_Type > filter_Type;
111 
112  typedef LifeV::ExporterEmpty<mesh_Type > emptyFilter_Type;
114  typedef LifeV::ExporterEnsight<mesh_Type > ensightFilter_Type;
116 
117 #ifdef HAVE_HDF5
120 #endif
121 
122 
123 
124  /** @name Constructors, destructor
125  */
126  //@{
127  Structure ( int argc,
128  char** argv,
129  std::shared_ptr<Epetra_Comm> structComm );
130 
132  {}
133  //@}
134 
135  //@{
136  void run()
137  {
138  run3d();
139  }
140  //@}
141 
142 protected:
143 
144 private:
145 
146  /**
147  * run the 2D driven cylinder simulation
148  */
149  void run2d();
150 
151  /**
152  * run the 3D driven cylinder simulation
153  */
154  void run3d();
155 
156 private:
157  struct Private;
161 };
162 
163 
164 
165 struct Structure::Private
166 {
168  rho (1),
169  poisson (1),
170  young (1),
171  bulk (1),
172  alpha (1),
173  gamma (1)
174  {}
175  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > fct_type;
176  double rho, poisson, young, bulk, alpha, gamma;
177 
179 
181 
182 };
183 
184 
185 
186 Structure::Structure ( int argc,
187  char** argv,
188  std::shared_ptr<Epetra_Comm> structComm) :
189  parameters ( new Private() )
190 {
191  GetPot command_line (argc, argv);
192  std::string data_file_name = command_line.follow ("data", 2, "-f", "--file");
193  GetPot dataFile ( data_file_name );
194  parameters->data_file_name = data_file_name;
195 
196  parameters->rho = dataFile ( "solid/physics/density", 1. );
197  parameters->young = dataFile ( "solid/physics/young", 1. );
198  parameters->poisson = dataFile ( "solid/physics/poisson", 1. );
199  parameters->bulk = dataFile ( "solid/physics/bulk", 1. );
200  parameters->alpha = dataFile ( "solid/physics/alpha", 1. );
201  parameters->gamma = dataFile ( "solid/physics/gamma", 1. );
202 
203  std::cout << "density = " << parameters->rho << std::endl
204  << "young = " << parameters->young << std::endl
205  << "poisson = " << parameters->poisson << std::endl
206  << "bulk = " << parameters->bulk << std::endl
207  << "alpha = " << parameters->alpha << std::endl
208  << "gamma = " << parameters->gamma << std::endl;
209 
210  parameters->comm = structComm;
211  int ntasks = parameters->comm->NumProc();
212 
213  if (!parameters->comm->MyPID() )
214  {
215  std::cout << "My PID = " << parameters->comm->MyPID() << " out of " << ntasks << " running." << std::endl;
216  }
217 }
218 
219 
220 
221 void
222 Structure::run2d()
223 {
224  std::cout << "2D cylinder test case is not available yet\n";
225 }
226 
227 
228 
229 void
230 Structure::run3d()
231 {
232  // General typedefs
233  typedef StructuralOperator<mesh_Type >::vector_Type vector_Type;
234  typedef std::shared_ptr<vector_Type> vectorPtr_Type;
235  typedef FESpace< mesh_Type, MapEpetra > solidFESpace_Type;
236  typedef std::shared_ptr<solidFESpace_Type> solidFESpacePtr_Type;
237 
238  typedef ETFESpace< RegionMesh<LinearTetra>, MapEpetra, 3, 3 > solidETFESpace_Type;
239  typedef std::shared_ptr<solidETFESpace_Type> solidETFESpacePtr_Type;
240 
241  bool verbose = (parameters->comm->MyPID() == 0);
242 
243  //! dataElasticStructure for parameters
244  GetPot dataFile ( parameters->data_file_name.c_str() );
245  std::shared_ptr<StructuralConstitutiveLawData> dataStructure (new StructuralConstitutiveLawData( ) );
246  dataStructure->setup (dataFile);
247 
248  //! Read and partition mesh
249  MeshData meshData;
250  meshData.setup (dataFile, "solid/space_discretization");
251 
252  std::shared_ptr<mesh_Type > fullMeshPtr (new RegionMesh<LinearTetra> ( ( parameters->comm ) ) );
253  readMesh (*fullMeshPtr, meshData);
254 
255  MeshPartitioner< mesh_Type > meshPart ( fullMeshPtr, parameters->comm );
256 
257  //! Functional spaces - needed for the computations of the gradients
258  std::string dOrder = dataFile ( "solid/space_discretization/order", "P1");
259  solidFESpacePtr_Type dFESpace ( new solidFESpace_Type (meshPart, dOrder, 3, parameters->comm) );
260  solidETFESpacePtr_Type dETFESpace ( new solidETFESpace_Type (meshPart, & (dFESpace->refFE() ), & (dFESpace->fe().geoMap() ), parameters->comm) );
261 
262 
263  if (verbose)
264  {
265  std::cout << std::endl;
266  }
267 
268  //! 1. Constructor of the class to compute the tensions
269  StructuralOperator<RegionMesh<LinearTetra> > solid;
270 
271  //! face BChandler object
272  std::shared_ptr<BCHandler> BCh ( new BCHandler() );
273 
274  //! 2. Its setup
275  solid.setup (dataStructure,
276  dFESpace,
277  dETFESpace,
278  BCh,
279  parameters->comm);
280 
281  //! 3. Creation of the importers to read the displacement field
282  std::string const filename = dataFile ( "importer/filename", "structure");
283  std::string const importerType = dataFile ( "importer/type", "hdf5");
284 
285  std::string iterationString; //useful to iterate over the strings
286 
287  if (verbose)
288  {
289  std::cout << "The filename is : " << filename << std::endl;
290  std::cout << "The importerType is: " << importerType << std::endl;
291  }
292 
293 #ifdef HAVE_HDF5
294  if (importerType.compare ("hdf5") == 0)
295  {
296  M_importer.reset ( new hdf5Filter_Type (dataFile, filename) );
297  }
298  else
299 #endif
300  {
301  if (importerType.compare ("none") == 0)
302  {
303  M_importer.reset ( new emptyFilter_Type ( dataFile, dFESpace->mesh(), "solid", dFESpace->map().comm().MyPID() ) );
304  }
305  else
306  {
307  M_importer.reset ( new ensightFilter_Type ( dataFile, filename ) );
308  }
309  }
310  M_importer->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
311 
312  // The vector where the solution will be stored
313  vectorPtr_Type solidDisp (new vector_Type (dFESpace->map(), LifeV::Unique ) );
314 
315  //! 6. Post-processing setting
316  std::shared_ptr< Exporter<RegionMesh<LinearTetra> > > exporter;
317 
318  std::string const exporterType = dataFile ( "exporter/type", "hdf5");
319  std::string const nameExporter = dataFile ( "exporter/name", "jacobian");
320 
321 #ifdef HAVE_HDF5
322  if (exporterType.compare ("hdf5") == 0)
323  {
324  M_exporter.reset ( new hdf5Filter_Type ( dataFile, nameExporter ) );
325  }
326  else
327 #endif
328  {
329  if (exporterType.compare ("none") == 0)
330  {
331  M_exporter.reset ( new emptyFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
332  }
333 
334  else
335  {
336  M_exporter.reset ( new ensightFilter_Type ( dataFile, meshPart.meshPartition(), nameExporter, parameters->comm->MyPID() ) ) ;
337  }
338  }
339 
340  M_exporter->setMeshProcId (dFESpace->mesh(), dFESpace->map().comm().MyPID() );
341 
342  vectorPtr_Type jacobianVector ( new vector_Type (solid.displacement(), LifeV::Unique ) );
343  vectorPtr_Type meshColors ( new vector_Type (solid.displacement(), LifeV::Unique ) );
344 
345  M_exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "determinantF", dFESpace, jacobianVector, UInt (0) );
346  M_exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "colors", dFESpace, meshColors, UInt (0) );
347  M_exporter->addVariable ( ExporterData<RegionMesh<LinearTetra> >::VectorField, "displacementField", dFESpace, solidDisp, UInt (0) );
348 
349  M_exporter->postProcess ( 0.0 );
350 
351  //! =================================================================================
352  //! Analysis - Istant or Interval
353  //! =================================================================================
354 
355  MPI_Barrier (MPI_COMM_WORLD);
356 
357  //! 5. For each interval, the analysis is performed
358  LifeV::Real dt = dataFile ( "solid/time_discretization/timestep", 0.0);
359  std::string const nameField = dataFile ( "importer/nameField", "displacement");
360 
361  //Get the iteration number
362  iterationString = dataFile ("importer/iteration", "00000");
363  LifeV::Real time = dataFile ("importer/time", 1.0);
364 
365  /*!Definition of the ExporterData, used to load the solution inside the previously defined vectors*/
366  LifeV::ExporterData<mesh_Type> solutionDispl (LifeV::ExporterData<mesh_Type>::VectorField, nameField + "." + iterationString, dFESpace, solidDisp, UInt (0), LifeV::ExporterData<mesh_Type>::UnsteadyRegime );
367 
368  //Read the variable
369  M_importer->readVariable (solutionDispl);
370  M_importer->closeFile();
371 
372  //Set the current solution as the displacement vector to use
373  solid.jacobianDistribution ( solidDisp, *jacobianVector);
374 
375  //color the mesh according to the marker of the volume
376  solid.colorMesh ( *meshColors );
377 
378  //Extracting the tensions
379  std::cout << std::endl;
380  std::cout << "Norm of the J = det(F) : " << jacobianVector->norm2() << std::endl;
381 
382  M_exporter->postProcess ( time );
383 
384  if (verbose )
385  {
386  std::cout << "Analysis Completed!" << std::endl;
387  }
388 
389  //Closing files
390  M_exporter->closeFile();
391 
392  if (verbose )
393  {
394  std::cout << "finished" << std::endl;
395  }
396 
397  MPI_Barrier (MPI_COMM_WORLD);
398  //!---------------------------------------------.-----------------------------------------------------
399 }
400 
401 int
402 main ( int argc, char** argv )
403 {
404 
405 #ifdef HAVE_MPI
406  MPI_Init (&argc, &argv);
407  std::shared_ptr<Epetra_MpiComm> Comm (new Epetra_MpiComm ( MPI_COMM_WORLD ) );
408  if ( Comm->MyPID() == 0 )
409  {
410  std::cout << "% using MPI" << std::endl;
411  }
412 #else
413  std::shared_ptr<Epetra_SerialComm> Comm ( new Epetra_SerialComm() );
414  std::cout << "% using serial Version" << std::endl;
415 #endif
416 
417  Structure structure ( argc, argv, Comm );
418  structure.run();
419 
420 #ifdef HAVE_MPI
421  MPI_Finalize();
422 #endif
423  return returnValue ;
424 }
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)
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
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
std::set< UInt > parseList(const std::string &list)
std::shared_ptr< LifeV::Exporter< mesh_Type > > filterPtr_Type