LifeV
2d/darcy.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 darcy.cpp
28  @author A. Fumagalli <alessio.fumagalli@mail.polimi.it>
29  @date 2012-06-13
30  */
31 
32 // ===================================================
33 //! Includes
34 // ===================================================
35 
36 #include "darcy.hpp"
37 #include "user_fun.hpp"
38 
39 //#if __cplusplus < 201103L
40 //#error "ERRORE C++"
41 //#endif
42 
43 // ===================================================
44 //! Namespaces & define
45 // ===================================================
46 
47 using namespace LifeV;
48 
49 // ===================================================
50 //! Standard functions
51 // ===================================================
52 
53 Real UOne ( const Real& /* t */,
54  const Real& /* x */,
55  const Real& /* y */,
56  const Real& /* z */,
57  const ID& /* icomp */)
58 {
59  return 1.;
60 }
61 
62 // ===================================================
63 //! Private Members
64 // ===================================================
65 
67 {
68  Private() {}
69 
70  // Policy for scalar functions
71  typedef std::function < Real ( const Real&, const Real&,
72  const Real&, const Real&, const ID& ) >
74 
78 
80 
81  // Function Types
82 
84  {
85  fct_type f;
86  f = std::bind ( &UOne, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, std::placeholders::_4, std::placeholders::_5 );
87  return f;
88  }
89 
91  {
92  fct_type f;
93  f = std::bind ( &dataProblem::analyticalSolution, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, std::placeholders::_4, std::placeholders::_5 );
94  return f;
95  }
96 
98  {
99  fct_type f;
100  f = std::bind ( &dataProblem::analyticalFlux, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, std::placeholders::_4, std::placeholders::_5 );
101  return f;
102  }
103 
104 };
105 
106 // ===================================================
107 //! Constructors
108 // ===================================================
109 
110 darcy_nonlinear::darcy_nonlinear ( int argc, char** argv )
111  : Members ( new Private )
112 {
113  GetPot command_line (argc, argv);
114  Members->data_file_name = command_line.follow ("data", 2, "-f", "--file");
115  Members->xml_file_name = command_line.follow ("parameterList.xml", "--xml");
116 
117  GetPot dataFile ( Members->data_file_name );
118 
119  Members->discretization_section = "darcy";
120 
121 #ifdef EPETRA_MPI
122  Members->comm.reset ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
123 #else
124  Members->comm.reset ( new Epetra_SerialComm() );
125 #endif
126 
127 }
128 
129 // ===================================================
130 //! Methods
131 // ===================================================
132 
133 Real
135 {
136  using namespace dataProblem;
137 
138  // Life chonos
139  LifeChrono chronoTotal;
140  LifeChrono chronoReadAndPartitionMesh;
141  LifeChrono chronoBoundaryCondition;
142  LifeChrono chronoFiniteElementSpace;
143  LifeChrono chronoProblem;
144  LifeChrono chronoProcess;
145  LifeChrono chronoError;
146 
147  // Displayer to print on screen
148  std::shared_ptr < Displayer > displayer ( new Displayer ( Members->comm ) );
149 
150  // Start chronoTotal for measure the total time for the computation
151  chronoTotal.start();
152 
153  // Reading from data file
154  GetPot dataFile ( Members->data_file_name.c_str() );
155 
156  //
157  // The Darcy Solver
158  //
159 
160  displayer->leaderPrint ( "The Darcy solver\n" );
161 
162  // Start chronoReadAndPartitionMesh for measure the total time for the creation of the local meshes
163  chronoReadAndPartitionMesh.start();
164 
165  // Create the data file
166  darcyDataPtr_Type darcyData ( new darcyData_Type );
167 
168  // Set up the data
169  darcyData->setup ( dataFile );
170 
171  // Create a Teuchos parameter list for the linear algebra.
172  darcyData_Type::paramListPtr_Type linearAlgebraList = Teuchos::rcp ( new darcyData_Type::paramList_Type );
173  linearAlgebraList = Teuchos::getParametersFromXmlFile ( Members->xml_file_name );
174 
175  // Set the parameter list into the data darcy.
176  darcyData->setLinearAlgebraList ( linearAlgebraList );
177 
178  // Create the mesh file handler
179  MeshData meshData;
180 
181  // Set up the mesh file
182  meshData.setup ( dataFile, Members->discretization_section + "/space_discretization");
183 
184  // Create the the mesh
185  regionMeshPtr_Type fullMeshPtr ( new regionMesh_Type ( Members->comm ) );
186 
187  // Select if the mesh is structured or not
188  if ( meshData.meshType() != "structured" )
189  {
190  // Set up the mesh
191  readMesh ( *fullMeshPtr, meshData );
192  }
193  else
194  {
195  // Section of the structured mesh
196  const std::string structuredSection = Members->discretization_section + "/space_discretization/";
197 
198  // Set up the structured mesh
199  regularMesh2D ( *fullMeshPtr, 0,
200  dataFile ( ( structuredSection + "nx" ).data(), 4 ),
201  dataFile ( ( structuredSection + "ny" ).data(), 4 ),
202  dataFile ( ( structuredSection + "verbose" ).data(), false ),
203  dataFile ( ( structuredSection + "lx" ).data(), 1. ),
204  dataFile ( ( structuredSection + "ly" ).data(), 1. ) );
205  }
206 
207  // Create the local mesh ( local scope used to delete the meshPart object )
208  regionMeshPtr_Type meshPtr;
209  {
210  // Create the partitioner
211  MeshPartitioner< regionMesh_Type > meshPart;
212 
213  // Partition the mesh using ParMetis
214  meshPart.doPartition ( fullMeshPtr, Members->comm );
215 
216  // Get the local mesh
217  meshPtr = meshPart.meshPartition();
218  }
219 
220  // Stop chronoReadAndPartitionMesh
221  chronoReadAndPartitionMesh.stop();
222 
223  // The leader process print chronoReadAndPartitionMesh
224  displayer->leaderPrint ( "Time for read and partition the mesh ",
225  chronoReadAndPartitionMesh.diff(), "\n" );
226 
227  // Create the boundary conditions
228 
229  // Start chronoBoundaryCondition for measure the total time for create the boundary conditions
230  chronoBoundaryCondition.start();
231 
232  bcHandlerPtr_Type bcDarcy ( new bcHandler_Type );
233 
234  setBoundaryConditions ( bcDarcy );
235 
236  // Stop chronoBoundaryCondition
237  chronoBoundaryCondition.stop();
238 
239  // The leader process print chronoBoundaryCondition
240  displayer->leaderPrint ( "Time for create the boundary conditions handler ",
241  chronoBoundaryCondition.diff(), "\n" );
242 
243  // Create the solution spaces
244 
245  // Start chronoFiniteElementSpace for measure the total time for create the finite element spaces
246  chronoFiniteElementSpace.start();
247 
248  // Primal solution parameters
249  const ReferenceFE* refFE_primal ( static_cast<ReferenceFE*> (NULL) );
250  const QuadratureRule* qR_primal ( static_cast<QuadratureRule*> (NULL) );
251  const QuadratureRule* bdQr_primal ( static_cast<QuadratureRule*> (NULL) );
252 
253  refFE_primal = &feTriaP0;
254  qR_primal = &quadRuleTria4pt;
255  bdQr_primal = &quadRuleSeg1pt;
256 
257  // Dual solution parameters
258  const ReferenceFE* refFE_dual ( static_cast<ReferenceFE*> (NULL) );
259  const QuadratureRule* qR_dual ( static_cast<QuadratureRule*> (NULL) );
260  const QuadratureRule* bdQr_dual ( static_cast<QuadratureRule*> (NULL) );
261 
262  refFE_dual = &feTriaRT0;
263  qR_dual = &quadRuleTria4pt;
264  bdQr_dual = &quadRuleSeg1pt;
265 
266  // Interpolate of dual solution parameters
267  const ReferenceFE* refFE_dualInterpolate ( static_cast<ReferenceFE*> (NULL) );
268  const QuadratureRule* qR_dualInterpolate ( static_cast<QuadratureRule*> (NULL) );
269  const QuadratureRule* bdQr_dualInterpolate ( static_cast<QuadratureRule*> (NULL) );
270 
271  refFE_dualInterpolate = &feTriaP0;
272  qR_dualInterpolate = &quadRuleTria4pt;
273  bdQr_dualInterpolate = &quadRuleSeg1pt;
274 
275  // Hybrid solution parameters
276  // hybrid
277  const ReferenceFE* refFE_hybrid ( static_cast<ReferenceFE*> (NULL) );
278  const QuadratureRule* qR_hybrid ( static_cast<QuadratureRule*> (NULL) );
279  const QuadratureRule* bdQr_hybrid ( static_cast<QuadratureRule*> (NULL) );
280 
281  refFE_hybrid = &feTriaRT0Hyb;
282  qR_hybrid = &quadRuleTria4pt;
283  bdQr_hybrid = &quadRuleSeg1pt;
284 
285  // Finite element space of the primal variable
286  fESpacePtr_Type p_FESpacePtr ( new fESpace_Type ( meshPtr, *refFE_primal, *qR_primal,
287  *bdQr_primal, 1, Members->comm ) );
288 
289  // Finite element space of the dual variable
290  fESpacePtr_Type u_FESpacePtr ( new fESpace_Type ( meshPtr, *refFE_dual, *qR_dual,
291  *bdQr_dual, 1, Members->comm ) );
292 
293  // Finite element space of the hybrid variable
294  fESpacePtr_Type hybrid_FESpacePtr ( new fESpace_Type ( meshPtr, *refFE_hybrid, *qR_hybrid,
295  *bdQr_hybrid, 1, Members->comm ) );
296 
297  // Finite element space of the interpolation of dual variable
298  fESpacePtr_Type uInterpolate_FESpacePtr ( new fESpace_Type ( meshPtr, *refFE_dualInterpolate, *qR_dualInterpolate,
299  *bdQr_dualInterpolate, 2, Members->comm ) );
300 
301  // Stop chronoFiniteElementSpace
302  chronoFiniteElementSpace.stop();
303 
304  // The leader process print chronoFiniteElementSpace
305  displayer->leaderPrint ( "Time for create the finite element spaces ",
306  chronoFiniteElementSpace.diff(), "\n" );
307 
308  // Start chronoProblem for measure the total time for create the problem
309  chronoProblem.start();
310 
311  // Instantiation of the DarcySolver class
312  darcySolver_Type darcySolver;
313 
314  // Stop chronoProblem
315  chronoProblem.stop();
316 
317  // The leader process print chronoProblem
318  displayer->leaderPrint ( "Time for create the problem ", chronoProblem.diff(), "\n" );
319 
320  // Process the problem
321 
322  // Start chronoProcess for measure the total time for the simulation
323  chronoProcess.start ();
324 
325  // Set the data for the solver.
326  darcySolver.setData ( darcyData );
327 
328  // Set the displayer.
329  darcySolver.setDisplayer ( displayer );
330 
331  // Setup phase for the linear solver
332  darcySolver.setup ();
333 
334  // Create the fields for the solver
335 
336  // Scalar field for primal variable
337  scalarFieldPtr_Type primalField ( new scalarField_Type ( p_FESpacePtr ) );
338 
339  // Scalar field for dual variable
340  scalarFieldPtr_Type dualField ( new scalarField_Type ( u_FESpacePtr ) );
341 
342  // Scalar field for hybrid variable
343  scalarFieldPtr_Type hybridField ( new scalarField_Type ( hybrid_FESpacePtr ) );
344 
345  // Vector for the interpolated dual solution
346  vectorFieldPtr_Type dualInterpolated ( new vectorField_Type ( uInterpolate_FESpacePtr, Repeated ) );
347 
348  // Set the field for the solver
349 
350  // Set the h
351  darcySolver.setHybridField ( hybridField );
352 
353  // Set the primal field
354  darcySolver.setPrimalField ( primalField );
355 
356  // Set the dual field
357  darcySolver.setDualField ( dualField );
358 
359  // Set the source term
360  scalarFctPtr_Type scalarSourceFct ( new scalarSource );
361  darcySolver.setScalarSource ( scalarSourceFct );
362 
363  // Set the vector source term
364  vectorFctPtr_Type vectorSourceFct ( new vectorSource );
365  darcySolver.setVectorSource ( vectorSourceFct );
366 
367  // Set the reaction term
368  scalarFctPtr_Type reactionTermFct ( new reactionTerm );
369  darcySolver.setReactionTerm ( reactionTermFct );
370 
371  // Set the inverse of the permeability
372  matrixFctPtr_Type inversePermeabilityFct ( new inversePermeability );
373  darcySolver.setInversePermeability ( inversePermeabilityFct );
374 
375  // Set the initial primal variable
376  scalarFctPtr_Type initialPrimalFct ( new initialCondition );
377  darcySolver.setInitialPrimal ( initialPrimalFct );
378 
379  // Set the mass function
380  scalarFctPtr_Type massFct ( new massFunction );
381  darcySolver.setMass ( massFct );
382 
383  // Set the boudary conditions
384  darcySolver.setBoundaryConditions ( bcDarcy );
385 
386  // Set the exporter for the solution
387  std::shared_ptr< Exporter< regionMesh_Type > > exporter;
388 
389  // Shared pointer used in the exporter for the primal solution
390  vectorPtr_Type primalExporter;
391 
392  // Shared pointer used in the exporter for the dual solution
393  vectorPtr_Type dualExporter;
394 
395  // Type of the exporter
396  std::string const exporterType = dataFile ( "exporter/type", "none" );
397 
398  // The name of the file
399  const std::string exporterFileName = dataFile ( "exporter/file_name", "PressureVelocity" );
400 
401  // Choose the exporter
402 #ifdef HAVE_HDF5
403  if ( exporterType.compare ("hdf5") == 0 )
404  {
405  exporter.reset ( new ExporterHDF5 < regionMesh_Type > ( dataFile, exporterFileName ) );
406  }
407  else
408 #endif
409  {
410  exporter.reset ( new ExporterEmpty < regionMesh_Type > ( dataFile, exporterFileName ) );
411  }
412 
413  // Set directory where to save the solution
414  exporter->setPostDir ( dataFile ( "exporter/folder", "./" ) );
415 
416  exporter->setMeshProcId ( meshPtr, Members->comm->MyPID() );
417 
418  // Set the exporter primal pointer
419  primalExporter.reset ( new vector_Type ( primalField->getVector(), exporter->mapType() ) );
420 
421  // Add the primal variable to the exporter
422  exporter->addVariable ( ExporterData < regionMesh_Type >::ScalarField,
423  dataFile ( "exporter/name_primal", "Pressure" ),
424  p_FESpacePtr,
425  primalExporter,
426  static_cast<UInt> ( 0 ),
427  ExporterData < regionMesh_Type >::UnsteadyRegime,
428  ExporterData < regionMesh_Type >::Cell );
429 
430  // Set the exporter dual pointer
431  dualExporter.reset ( new vector_Type ( dualInterpolated->getVector(), exporter->mapType() ) );
432 
433  // Add the variable to the exporter
434  exporter->addVariable ( ExporterData < regionMesh_Type >::VectorField,
435  dataFile ( "exporter/name_dual", "Velocity" ),
436  uInterpolate_FESpacePtr,
437  dualExporter,
438  static_cast<UInt> ( 0 ),
439  ExporterData < regionMesh_Type >::UnsteadyRegime,
440  ExporterData < regionMesh_Type >::Cell );
441 
442  // Display the total number of unknowns
443  displayer->leaderPrint ( "Number of unknowns : ",
444  hybrid_FESpacePtr->map().map ( Unique )->NumGlobalElements(), "\n" );
445 
446  // Export the partitioning
447  exporter->exportPID ( meshPtr, Members->comm );
448 
449  // Copy the initial primal to the exporter
450  *primalExporter = primalField->getVector ();
451 
452  // Save the initial primal solution into the exporter
453  exporter->postProcess ( darcyData->dataTimePtr()->initialTime() );
454 
455  // Define if the current time is the last time step.
456  bool isLastTimeStep ( false );
457 
458  // Define the end time
459  const Real endTime = darcyData->dataTimePtr()->endTime();
460 
461  // Define the tolerance for the elapsed time
462  const Real tolTime = 1e-10;
463 
464  // Take the current time.
465  Real currentTime = darcyData->dataTimePtr()->time();
466 
467  // A loop for the simulation
468  while ( darcyData->dataTimePtr()->time() < endTime && !isLastTimeStep )
469  {
470  // Take the left time.
471  const Real leftTime = endTime - currentTime;
472 
473  // Check if the time step is consistent
474  if ( darcyData->dataTimePtr()->timeStep() > leftTime + tolTime )
475  {
476  // Compute the last time step.
477  darcyData->dataTimePtr()->setTimeStep ( leftTime );
478 
479  // This is the last time step in the simulation
480  isLastTimeStep = true;
481  }
482 
483  // Advance the current time of \Delta t.
484  darcyData->dataTimePtr()->updateTime();
485 
486  // Take the current time.
487  currentTime = darcyData->dataTimePtr()->time();
488 
489  // The leader process prints the temporal data.
490  if ( displayer->isLeader() )
491  {
492  darcyData->dataTimePtr()->showMe();
493  }
494 
495  // Solve the problem.
496  darcySolver.solve ();
497 
498  // Save the solutions.
499 
500  // Copy the primal solution to the exporter.
501  *primalExporter = primalField->getVector ();
502 
503  // Interpolate the dual vector field spammed as Raviart-Thomas into a P0 vector field.
504  dualInterpolated->getVector() = uInterpolate_FESpacePtr->feToFEInterpolate (
505  *u_FESpacePtr, dualField->getVector() );
506 
507  // Copy the dual interpolated solution to the exporter.
508  *dualExporter = dualInterpolated->getVector();
509 
510  // Save the solution into the exporter.
511  exporter->postProcess ( currentTime );
512 
513  }
514 
515  // Stop chronoProcess
516  chronoProcess.stop();
517 
518  // The leader process print chronoProcess
519  displayer->leaderPrint ( "Time for process ", chronoProcess.diff(), "\n" );
520 
521  // Compute the errors
522  // For non time dependences problem the time where the errors are computed is useless,
523  // but thanks to common interface we handle both cases.
524 
525  // Start chronoError for measure the total time for computing the errors.
526  chronoError.start();
527 
528  // Compute the error L2 norms
529  Real primalL2Norm (0), exactPrimalL2Norm (0), primalL2Error (0), primalL2RelativeError (0);
530  Real dualL2Norm (0), exactDualL2Norm (0), dualL2Error (0), dualL2RelativeError (0);
531 
532  // Norms and errors for the pressure
533  displayer->leaderPrint ( "\nPRESSURE ERROR\n" );
534 
535  // Compute the L2 norm for the primal solution
536  primalL2Norm = p_FESpacePtr->l2Norm ( primalField->getVector () );
537 
538  // Display the L2 norm for the primal solution
539  displayer->leaderPrint ( " L2 norm of primal unknown: ", primalL2Norm, "\n" );
540 
541  // Compute the L2 norm for the analytical primal
542  exactPrimalL2Norm = p_FESpacePtr->l2NormFunction ( Members->getAnalyticalSolution(),
543  darcyData->dataTimePtr()->endTime() );
544 
545  // Display the L2 norm for the analytical primal
546  displayer->leaderPrint ( " L2 norm of primal exact: ", exactPrimalL2Norm, "\n" );
547 
548  // Compute the L2 error for the primal solution
549  primalL2Error = p_FESpacePtr->l2ErrorWeighted ( Members->getAnalyticalSolution(),
550  primalField->getVector(),
551  Members->getUOne(),
552  darcyData->dataTimePtr()->endTime() );
553 
554  // Display the L2 error for the primal solution
555  displayer->leaderPrint ( " L2 error of primal unknown: ", primalL2Error, "\n" );
556 
557  // Compute the L2 realative error for the primal solution
558  primalL2RelativeError = primalL2Error / exactPrimalL2Norm;
559 
560  // Display the L2 relative error for the primal solution
561  displayer->leaderPrint ( " L2 relative error of primal unknown: ", primalL2RelativeError, "\n" );
562 
563  // Norms and errors for the interpolated dual
564  displayer->leaderPrint ( "\nINTERPOLATED DARCY VELOCITY ERROR\n" );
565 
566  // Compute the L2 norm for the interpolated dual solution
567  dualL2Norm = uInterpolate_FESpacePtr->l2Norm ( dualInterpolated->getVector() );
568 
569  // Display the L2 norm for the interpolated dual solution
570  displayer->leaderPrint ( " L2 norm of dual unknown: ", dualL2Norm, "\n" );
571 
572  // Compute the L2 norm for the analytical dual
573  exactDualL2Norm = uInterpolate_FESpacePtr->l2NormFunction ( Members->getAnalyticalFlux(),
574  darcyData->dataTimePtr()->endTime() );
575 
576  // Display the L2 norm for the analytical dual
577  displayer->leaderPrint ( " L2 norm of dual exact: ", exactDualL2Norm, "\n" );
578 
579  // Compute the L2 error for the dual solution
580  dualL2Error = uInterpolate_FESpacePtr->l2Error ( Members->getAnalyticalFlux(),
581  dualInterpolated->getVector(),
582  darcyData->dataTimePtr()->endTime(),
583  NULL );
584 
585  // Display the L2 error for the dual solution
586  displayer->leaderPrint ( " L2 error of dual unknown: ", dualL2Error, "\n" );
587 
588  // Compute the L2 relative error for the dual solution
589  dualL2RelativeError = dualL2Error / exactDualL2Norm;
590 
591  // Display the L2 relative error for the dual solution
592  displayer->leaderPrint ( " L2 relative error of Dual unknown: ", dualL2RelativeError, "\n" );
593 
594  // Stop chronoError
595  chronoError.stop();
596 
597  // The leader process print chronoError
598  displayer->leaderPrint ( "Time for compute errors ", chronoError.diff(), "\n" );
599 
600  // Stop chronoTotal
601  chronoTotal.stop();
602 
603  // The leader process print chronoTotal
604  displayer->leaderPrint ( "Total time for the computation ", chronoTotal.diff(), "\n" );
605 
606  // Return the error, needed for the succes/failure of the test
607  return primalL2Error;
608 
609 } // run
Real UOne(const Real &, const Real &, const Real &, const Real &, const ID &)
Standard functions.
Definition: 3d/darcy.cpp:49
Private Members.
Definition: 2d/darcy.cpp:66
darcySolver_Type::dataPtr_Type darcyDataPtr_Type
fct_type getAnalyticalSolution()
Definition: 2d/darcy.cpp:90
std::shared_ptr< Epetra_Comm > comm
Definition: 2d/darcy.cpp:79
LifeV::Real run()
To lunch the simulation.
Definition: 2d/darcy.cpp:134
darcy_nonlinear(int argc, char **argv)
Constructor.
Definition: 2d/darcy.cpp:110
std::function< Real(const Real &, const Real &, const Real &, const Real &, const ID &) > fct_type
Definition: 2d/darcy.cpp:73
uint32_type ID
IDs.
Definition: LifeV.hpp:194
fct_type getAnalyticalFlux()
Definition: 2d/darcy.cpp:97
double Real
Generic real data.
Definition: LifeV.hpp:175
std::string discretization_section
Definition: 2d/darcy.cpp:77
Includes.
Definition: 2d/darcy.hpp:143
std::string data_file_name
Definition: 2d/darcy.cpp:75
std::string xml_file_name
Definition: 2d/darcy.cpp:76