LifeV
lifev/core/testsuite/adr_assembler/2d/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 /*!
28  @file
29  @brief
30 
31  @author Samuel Quinodoz <samuel.quinodoz@epfl.ch>
32  @date 08-10-2010
33  */
34 
35 
36 #include <Epetra_ConfigDefs.h>
37 #ifdef EPETRA_MPI
38 #include <mpi.h>
39 #include <Epetra_MpiComm.h>
40 #else
41 #include <Epetra_SerialComm.h>
42 #endif
43 
44 
45 #include <lifev/core/LifeV.hpp>
46 
47 
48 #include <lifev/core/algorithm/PreconditionerIfpack.hpp>
49 #include <lifev/core/algorithm/PreconditionerML.hpp>
50 
51 #include <lifev/core/algorithm/SolverAztecOO.hpp>
52 #include <lifev/core/algorithm/LinearSolver.hpp>
53 
54 #include <lifev/core/array/MatrixEpetra.hpp>
55 
56 #include <lifev/core/filter/ExporterEnsight.hpp>
57 
58 #include <lifev/core/fem/FESpace.hpp>
59 #include <lifev/core/fem/BCManage.hpp>
60 
61 #include <lifev/core/mesh/MeshPartitioner.hpp>
62 #include <lifev/core/mesh/RegionMesh2DStructured.hpp>
63 #include <lifev/core/mesh/RegionMesh.hpp>
64 
65 #include <lifev/core/solver/ADRAssembler.hpp>
66 #include <lifev/core/mesh/MeshData.hpp>
67 
68 #include <Teuchos_ParameterList.hpp>
69 #include <Teuchos_XMLParameterListHelpers.hpp>
70 #include <Teuchos_RCP.hpp>
71 
72 using namespace LifeV;
73 
74 namespace
75 {
76 static bool regIF = (PRECFactory::instance().registerProduct ( "Ifpack", &createIfpack ) );
77 static bool regML = (PRECFactory::instance().registerProduct ( "ML", &createML ) );
78 }
79 
80 //#define TEST_MASS
81 //#define TEST_ADVECTION
82 #define TEST_RHS
83 
84 
85 #ifdef TEST_MASS
86 Real epsilon (1);
87 
88 Real exactSolution ( const Real& /* t */, const Real& x, const Real& /* y */, const Real& /* z */ , const ID& /* i */ )
89 {
90  Real seps (s td::sqrt (epsilon) );
91  return std::exp (seps * x) / (std::exp (seps) - std::exp (-seps) );
92 }
93 #endif
94 
95 #ifdef TEST_ADVECTION
96 Real epsilon (1);
97 
98 Real exactSolution ( const Real& /* t */, const Real& x, const Real& /* y */, const Real& /* z */, const ID& /* i */ )
99 {
100  return ( std::exp (x / epsilon) - 1 ) / ( std::exp (1 / epsilon) - 1 );
101 }
102 
103 Real betaFct ( const Real& /* t */, const Real& /* x */, const Real& /* y */, const Real& /* z */, const ID& i )
104 {
105  return (i == 0);
106 }
107 #endif
108 
109 #ifdef TEST_RHS
110 Real epsilon (1);
111 
112 Real exactSolution ( const Real& /* t */, const Real& x, const Real& y, const Real& /* z */, const ID& /* i */ )
113 {
114  return std::sin (x) + y * y / 2;
115 }
116 
117 
118 Real fRhs ( const Real& /* t */, const Real& x, const Real& /* y */, const Real& /* z */ , const ID& /* i */ )
119 {
120  return std::sin (x) - 1;
121 }
122 #endif
123 
124 
129 typedef FESpace<mesh_Type, MapEpetra> feSpace_Type;
130 
135 
137 
138 
139 int
140 main ( int argc, char** argv )
141 {
142 
143 #ifdef HAVE_MPI
144  MPI_Init (&argc, &argv);
145  std::shared_ptr<Epetra_Comm> Comm (new Epetra_MpiComm (MPI_COMM_WORLD) );
146 #else
147  std::shared_ptr<Epetra_Comm> Comm (new Epetra_SerialComm);
148 #endif
149 
150  const bool verbose (Comm->MyPID() == 0);
151 
152  // Read first the data needed
153 
154  if (verbose)
155  {
156  std::cout << " -- Reading the data ... " << std::flush;
157  }
158  GetPot dataFile ( "data" );
159  if (verbose)
160  {
161  std::cout << " done ! " << std::endl;
162  }
163 
164 
165  // Build and partition the mesh
166 
167  if (verbose)
168  {
169  std::cout << " -- Reading the mesh ... " << std::flush;
170  }
171  MeshData meshData (dataFile, "mesh");
172  std::shared_ptr< mesh_Type > fullMeshPtr ( new mesh_Type ( Comm ) );
173 
174  // Select if the mesh is structured or not
175  if ( meshData.meshType() != "structured" )
176  {
177  // Set up the mesh
178  readMesh ( *fullMeshPtr, meshData );
179  }
180  else
181  {
182  // Set up the structured mesh
183  regularMesh2D ( *fullMeshPtr, 0,
184  dataFile ( "mesh/nx", 20 ),
185  dataFile ( "mesh/ny", 20 ),
186  dataFile ( "mesh/verbose", false ),
187  dataFile ( "mesh/lx", 1. ),
188  dataFile ( "mesh/ly", 1. ) );
189  }
190 
191  if (verbose)
192  {
193  std::cout << " done ! " << std::endl;
194  }
195 
196  if (verbose)
197  {
198  std::cout << " -- Partitioning the mesh ... " << std::flush;
199  }
200  std::shared_ptr< mesh_Type > meshPtr;
201  {
202  MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, Comm);
203  meshPtr = meshPart.meshPartition();
204  }
205  if (verbose)
206  {
207  std::cout << " done ! " << std::endl;
208  }
209 
210  if (verbose)
211  {
212  std::cout << " -- Freeing the global mesh ... " << std::flush;
213  }
214  fullMeshPtr.reset();
215  if (verbose)
216  {
217  std::cout << " done ! " << std::endl;
218  }
219 
220  // Build the FESpaces
221 
222  if (verbose)
223  {
224  std::cout << " -- Building FESpaces ... " << std::flush;
225  }
226  std::string uOrder ("P1");
227  std::string bOrder ("P1");
228  feSpacePtr_Type uFESpace ( new feSpace_Type ( meshPtr, uOrder, 1, Comm ) );
229  feSpacePtr_Type betaFESpace ( new feSpace_Type ( meshPtr, bOrder, 3, Comm ) );
230  if (verbose)
231  {
232  std::cout << " done ! " << std::endl;
233  }
234  if (verbose)
235  {
236  std::cout << " ---> Dofs: " << uFESpace->dof().numTotalDof() << std::endl;
237  }
238 
239  // Build the assembler and the matrices
240 
241  if (verbose)
242  {
243  std::cout << " -- Building assembler ... " << std::flush;
244  }
245  ADRAssembler<mesh_Type, matrix_Type, vector_Type> adrAssembler;
246  if (verbose)
247  {
248  std::cout << " done! " << std::endl;
249  }
250 
251  if (verbose)
252  {
253  std::cout << " -- Setting up assembler ... " << std::flush;
254  }
255  adrAssembler.setup (uFESpace, betaFESpace);
256  if (verbose)
257  {
258  std::cout << " done! " << std::endl;
259  }
260 
261  if (verbose)
262  {
263  std::cout << " -- Defining the matrix ... " << std::flush;
264  }
265  std::shared_ptr<matrix_Type> systemMatrix (new matrix_Type ( uFESpace->map() ) );
266  *systemMatrix *= 0.0;
267  if (verbose)
268  {
269  std::cout << " done! " << std::endl;
270  }
271 
272  // Perform the assembly of the matrix
273 
274  if (verbose)
275  {
276  std::cout << " -- Adding the diffusion ... " << std::flush;
277  }
278  adrAssembler.addDiffusion (systemMatrix, epsilon);
279  if (verbose)
280  {
281  std::cout << " done! " << std::endl;
282  }
283  if (verbose)
284  {
285  std::cout << " Time needed : " << adrAssembler.diffusionAssemblyChrono().diffCumul() << std::endl;
286  }
287 
288 #ifdef TEST_ADVECTION
289  if (verbose)
290  {
291  std::cout << " -- Adding the advection ... " << std::flush;
292  }
293  vector_Type beta (betaFESpace->map(), Repeated);
294  betaFESpace->interpolate (betaFct, beta, 0.0);
295  adrAssembler.addAdvection (systemMatrix, beta);
296  if (verbose)
297  {
298  std::cout << " done! " << std::endl;
299  }
300 #endif
301 #ifdef TEST_MASS
302  if (verbose)
303  {
304  std::cout << " -- Adding the mass ... " << std::flush;
305  }
306  adrAssembler.addMass (systemMatrix, 1.0);
307  if (verbose)
308  {
309  std::cout << " done! " << std::endl;
310  }
311 #endif
312 
313  if (verbose)
314  {
315  std::cout << " -- Closing the matrix ... " << std::flush;
316  }
317  systemMatrix->globalAssemble();
318  if (verbose)
319  {
320  std::cout << " done ! " << std::endl;
321  }
322 
323 #ifdef TEST_RHS
324  Real matrixNorm (systemMatrix->norm1() );
325  if (verbose)
326  {
327  std::cout << " ---> Norm 1 : " << matrixNorm << std::endl;
328  }
329  if ( std::fabs (matrixNorm - 8 ) > 1e-3)
330  {
331  std::cout << " <!> Matrix has changed !!! <!> " << std::endl;
332  return EXIT_FAILURE;
333  }
334 #endif
335 
336  // Definition and assembly of the RHS
337 
338  if (verbose)
339  {
340  std::cout << " -- Building the RHS ... " << std::flush;
341  }
342  //vector_Type rhs(uFESpace->map(),Unique);
343  vector_Type rhs (uFESpace->map(), Repeated);
344  rhs *= 0.0;
345 
346 #ifdef TEST_RHS
347  vector_Type fInterpolated (uFESpace->map(), Repeated);
348  fInterpolated *= 0.0;
349  uFESpace->interpolate ( static_cast<feSpace_Type::function_Type> ( fRhs ), fInterpolated, 0.0 );
350  adrAssembler.addMassRhs (rhs, fInterpolated);
351  rhs.globalAssemble();
352 #endif
353 
354  if (verbose)
355  {
356  std::cout << " done ! " << std::endl;
357  }
358 
359  // Definition and application of the BCs
360 
361  if (verbose)
362  {
363  std::cout << " -- Building the BCHandler ... " << std::flush;
364  }
365  BCHandler bchandler;
366  BCFunctionBase BCu ( static_cast<feSpace_Type::function_Type> ( exactSolution ) );
367  bchandler.addBC ("Dirichlet", 1, Essential, Full, BCu, 1);
368  for (UInt i (2); i <= 4; ++i)
369  {
370  bchandler.addBC ("Dirichlet", i, Essential, Full, BCu, 1);
371  }
372  if (verbose)
373  {
374  std::cout << " done ! " << std::endl;
375  }
376 
377  if (verbose)
378  {
379  std::cout << " -- Updating the BCs ... " << std::flush;
380  }
381  bchandler.bcUpdate (*uFESpace->mesh(), uFESpace->feBd(), uFESpace->dof() );
382  if (verbose)
383  {
384  std::cout << " done ! " << std::endl;
385  }
386 
387  if (verbose)
388  {
389  std::cout << " -- Applying the BCs ... " << std::flush;
390  }
391  vectorPtr_Type rhsBC ( new vector_Type (rhs, Unique) );
392  bcManage (*systemMatrix, *rhsBC, *uFESpace->mesh(), uFESpace->dof(), bchandler, uFESpace->feBd(), 1.0, 0.0);
393  rhs = *rhsBC;
394  if (verbose)
395  {
396  std::cout << " done ! " << std::endl;
397  }
398 
399  //************* SPY ***********
400  //systemMatrix->spy("matrix");
401  //rhs.spy("vector");
402  //*****************************
403 
404  // Definition of the solver
405 
406  if (verbose)
407  {
408  std::cout << " -- Building the solver ... " << std::flush;
409  }
410  LinearSolver linearSolver;
411  if (verbose)
412  {
413  std::cout << " done ! " << std::endl;
414  }
415 
416  if (verbose)
417  {
418  std::cout << " -- Setting up the solver ... " << std::flush;
419  }
420 
421  Teuchos::RCP< Teuchos::ParameterList > belosList = Teuchos::rcp ( new Teuchos::ParameterList );
422  belosList = Teuchos::getParametersFromXmlFile ( "SolverParamList.xml" );
423 
424  linearSolver.setCommunicator ( Comm );
425  linearSolver.setParameters ( *belosList );
426 
427  prec_Type* precRawPtr;
428  basePrecPtr_Type precPtr;
429  precRawPtr = new prec_Type;
430  precRawPtr->setDataFromGetPot ( dataFile, "prec" );
431  precPtr.reset ( precRawPtr );
432  linearSolver.setPreconditioner ( precPtr );
433 
434  if (verbose)
435  {
436  std::cout << " done ! " << std::endl;
437  }
438 
439  if (verbose)
440  {
441  std::cout << " -- Setting matrix in the solver ... " << std::flush;
442  }
443  linearSolver.setOperator ( systemMatrix );
444 
445  if (verbose)
446  {
447  std::cout << " done ! " << std::endl;
448  }
449 
450  // Definition of the solution
451 
452  if (verbose)
453  {
454  std::cout << " -- Defining the solution ... " << std::flush;
455  }
456  vectorPtr_Type solution ( new vector_Type (uFESpace->map(), Unique) );
457  solution->zero();
458  if (verbose)
459  {
460  std::cout << " done ! " << std::endl;
461  }
462 
463  // Solve the solution
464 
465  if (verbose)
466  {
467  std::cout << " -- Solving the system ... " << std::flush;
468  }
469  linearSolver.setRightHandSide ( rhsBC );
470  linearSolver.solve ( solution );
471 
472  if (verbose)
473  {
474  std::cout << " done ! " << std::endl;
475  }
476 
477  //************* SPY ***********
478  //solution.spy("solution");
479  //*****************************
480 
481  // Error computation
482 
483  if (verbose)
484  {
485  std::cout << " -- Computing the error ... " << std::flush;
486  }
487  vector_Type solutionErr (*solution);
488  solutionErr *= 0.0;
489  uFESpace->interpolate ( static_cast<feSpace_Type::function_Type> ( exactSolution ), solutionErr, 0.0 );
490  solutionErr -= *solution;
491  solutionErr.abs();
492  Real l2error (uFESpace->l2Error (exactSolution, vector_Type (*solution, Repeated), 0.0) );
493  if (verbose)
494  {
495  std::cout << " -- done ! " << std::endl;
496  }
497  if (verbose)
498  {
499  std::cout << " ---> Norm L2 : " << l2error << std::endl;
500  }
501  Real linferror ( solutionErr.normInf() );
502  if (verbose)
503  {
504  std::cout << " ---> Norm Inf : " << linferror << std::endl;
505  }
506 
507 
508  if (l2error > 0.0026)
509  {
510  std::cout << " <!> Solution has changed !!! <!> " << std::endl;
511  return EXIT_FAILURE;
512  }
513  if (linferror > 0.000016)
514  {
515  std::cout << " <!> Solution has changed !!! <!> " << std::endl;
516  return EXIT_FAILURE;
517  }
518 
519  // Exporter definition and use
520 
521  if (verbose)
522  {
523  std::cout << " -- Defining the exporter ... " << std::flush;
524  }
525  ExporterEnsight<mesh_Type> exporter ( dataFile, meshPtr, "solution", Comm->MyPID() ) ;
526  if (verbose)
527  {
528  std::cout << " done ! " << std::endl;
529  }
530 
531  if (verbose)
532  {
533  std::cout << " -- Defining the exported quantities ... " << std::flush;
534  }
535 
536  std::shared_ptr<vector_Type> solutionPtr (new vector_Type (*solution, Repeated) );
537  std::shared_ptr<vector_Type> solutionErrPtr (new vector_Type (solutionErr, Repeated) );
538 
539  if (verbose)
540  {
541  std::cout << " done ! " << std::endl;
542  }
543 
544  if (verbose)
545  {
546  std::cout << " -- Updating the exporter ... " << std::flush;
547  }
548  exporter.addVariable ( ExporterData<mesh_Type>::ScalarField, "solution", uFESpace, solutionPtr, UInt (0) );
549  exporter.addVariable ( ExporterData<mesh_Type>::ScalarField, "error", uFESpace, solutionErrPtr, UInt (0) );
550  if (verbose)
551  {
552  std::cout << " done ! " << std::endl;
553  }
554 
555  if (verbose)
556  {
557  std::cout << " -- Exporting ... " << std::flush;
558  }
559  exporter.postProcess (0);
560  if (verbose)
561  {
562  std::cout << " done ! " << std::endl;
563  }
564 
565  if (verbose)
566  {
567  std::cout << "End Result: TEST PASSED" << std::endl;
568  }
569 
570 #ifdef HAVE_MPI
571  MPI_Finalize();
572 #endif
573 
574  return ( EXIT_SUCCESS );
575 }
VectorEpetra - The Epetra Vector format Wrapper.
LifeV::Preconditioner basePrec_Type
Real exactSolution(const Real &t, const Real &x, const Real &y, const Real &z, const ID &)
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
MatrixEpetra< Real > matrix_Type
Real epsilon(1)
LifeV::PreconditionerIfpack prec_Type
FESpace< mesh_Type, MapEpetra > feSpace_Type
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
uint32_type ID
IDs.
Definition: LifeV.hpp:194
double Real
Generic real data.
Definition: LifeV.hpp:175
std::shared_ptr< prec_Type > precPtr_Type
std::shared_ptr< basePrec_Type > basePrecPtr_Type
std::shared_ptr< feSpace_Type > feSpacePtr_Type
Real fRhs(const Real &, const Real &x, const Real &, const Real &, const ID &)
int main(int argc, char **argv)