36 #include <Epetra_ConfigDefs.h> 39 #include <Epetra_MpiComm.h> 41 #include <Epetra_SerialComm.h> 45 #include <lifev/core/LifeV.hpp> 48 #include <lifev/core/algorithm/PreconditionerIfpack.hpp> 49 #include <lifev/core/algorithm/PreconditionerML.hpp> 51 #include <lifev/core/algorithm/SolverAztecOO.hpp> 52 #include <lifev/core/algorithm/LinearSolver.hpp> 54 #include <lifev/core/array/MatrixEpetra.hpp> 56 #include <lifev/core/filter/ExporterEnsight.hpp> 58 #include <lifev/core/fem/FESpace.hpp> 59 #include <lifev/core/fem/BCManage.hpp> 61 #include <lifev/core/mesh/MeshPartitioner.hpp> 62 #include <lifev/core/mesh/RegionMesh2DStructured.hpp> 63 #include <lifev/core/mesh/RegionMesh.hpp> 65 #include <lifev/core/solver/ADRAssembler.hpp> 66 #include <lifev/core/mesh/MeshData.hpp> 68 #include <Teuchos_ParameterList.hpp> 69 #include <Teuchos_XMLParameterListHelpers.hpp> 70 #include <Teuchos_RCP.hpp> 72 using namespace LifeV;
88 Real exactSolution (
const Real& ,
const Real& x,
const Real& ,
const Real& ,
const ID& )
90 Real seps (s td::sqrt (epsilon) );
91 return std::exp (seps * x) / (std::exp (seps) - std::exp (-seps) );
98 Real exactSolution (
const Real& ,
const Real& x,
const Real& ,
const Real& ,
const ID& )
100 return ( std::exp (x / epsilon) - 1 ) / ( std::exp (1 / epsilon) - 1 );
103 Real betaFct (
const Real& ,
const Real& ,
const Real& ,
const Real& ,
const ID& i )
114 return std::sin (x) + y * y / 2;
120 return std::sin (x) - 1;
144 MPI_Init (&argc, &argv);
145 std::shared_ptr<Epetra_Comm> Comm (
new Epetra_MpiComm (MPI_COMM_WORLD) );
147 std::shared_ptr<Epetra_Comm> Comm (
new Epetra_SerialComm);
150 const bool verbose (Comm->MyPID() == 0);
156 std::cout <<
" -- Reading the data ... " << std::flush;
158 GetPot dataFile (
"data" );
161 std::cout <<
" done ! " << std::endl;
169 std::cout <<
" -- Reading the mesh ... " << std::flush;
171 MeshData meshData (dataFile,
"mesh");
172 std::shared_ptr< mesh_Type > fullMeshPtr (
new mesh_Type ( Comm ) );
175 if ( meshData.meshType() !=
"structured" )
178 readMesh ( *fullMeshPtr, meshData );
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. ) );
193 std::cout <<
" done ! " << std::endl;
198 std::cout <<
" -- Partitioning the mesh ... " << std::flush;
200 std::shared_ptr< mesh_Type > meshPtr;
202 MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, Comm);
203 meshPtr = meshPart.meshPartition();
207 std::cout <<
" done ! " << std::endl;
212 std::cout <<
" -- Freeing the global mesh ... " << std::flush;
217 std::cout <<
" done ! " << std::endl;
224 std::cout <<
" -- Building FESpaces ... " << std::flush;
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 ) );
232 std::cout <<
" done ! " << std::endl;
236 std::cout <<
" ---> Dofs: " << uFESpace->dof().numTotalDof() << std::endl;
243 std::cout <<
" -- Building assembler ... " << std::flush;
245 ADRAssembler<mesh_Type, matrix_Type, vector_Type> adrAssembler;
248 std::cout <<
" done! " << std::endl;
253 std::cout <<
" -- Setting up assembler ... " << std::flush;
255 adrAssembler.setup (uFESpace, betaFESpace);
258 std::cout <<
" done! " << std::endl;
263 std::cout <<
" -- Defining the matrix ... " << std::flush;
265 std::shared_ptr<matrix_Type> systemMatrix (
new matrix_Type ( uFESpace->map() ) );
266 *systemMatrix *= 0.0;
269 std::cout <<
" done! " << std::endl;
276 std::cout <<
" -- Adding the diffusion ... " << std::flush;
278 adrAssembler.addDiffusion (systemMatrix, epsilon);
281 std::cout <<
" done! " << std::endl;
285 std::cout <<
" Time needed : " << adrAssembler.diffusionAssemblyChrono().diffCumul() << std::endl;
288 #ifdef TEST_ADVECTION 291 std::cout <<
" -- Adding the advection ... " << std::flush;
293 vector_Type beta (betaFESpace->map(), Repeated);
294 betaFESpace->interpolate (betaFct, beta, 0.0);
295 adrAssembler.addAdvection (systemMatrix, beta);
298 std::cout <<
" done! " << std::endl;
304 std::cout <<
" -- Adding the mass ... " << std::flush;
306 adrAssembler.addMass (systemMatrix, 1.0);
309 std::cout <<
" done! " << std::endl;
315 std::cout <<
" -- Closing the matrix ... " << std::flush;
317 systemMatrix->globalAssemble();
320 std::cout <<
" done ! " << std::endl;
324 Real matrixNorm (systemMatrix->norm1() );
327 std::cout <<
" ---> Norm 1 : " << matrixNorm << std::endl;
329 if ( std::fabs (matrixNorm - 8 ) > 1e-3)
331 std::cout <<
" <!> Matrix has changed !!! <!> " << std::endl;
340 std::cout <<
" -- Building the RHS ... " << std::flush;
343 vector_Type rhs (uFESpace->map(), Repeated);
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();
356 std::cout <<
" done ! " << std::endl;
363 std::cout <<
" -- Building the BCHandler ... " << std::flush;
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)
370 bchandler.addBC (
"Dirichlet", i, Essential, Full, BCu, 1);
374 std::cout <<
" done ! " << std::endl;
379 std::cout <<
" -- Updating the BCs ... " << std::flush;
381 bchandler.bcUpdate (*uFESpace->mesh(), uFESpace->feBd(), uFESpace->dof() );
384 std::cout <<
" done ! " << std::endl;
389 std::cout <<
" -- Applying the BCs ... " << std::flush;
391 vectorPtr_Type rhsBC (
new vector_Type (rhs, Unique) );
392 bcManage (*systemMatrix, *rhsBC, *uFESpace->mesh(), uFESpace->dof(), bchandler, uFESpace->feBd(), 1.0, 0.0);
396 std::cout <<
" done ! " << std::endl;
408 std::cout <<
" -- Building the solver ... " << std::flush;
410 LinearSolver linearSolver;
413 std::cout <<
" done ! " << std::endl;
418 std::cout <<
" -- Setting up the solver ... " << std::flush;
421 Teuchos::RCP< Teuchos::ParameterList > belosList = Teuchos::rcp (
new Teuchos::ParameterList );
422 belosList = Teuchos::getParametersFromXmlFile (
"SolverParamList.xml" );
424 linearSolver.setCommunicator ( Comm );
425 linearSolver.setParameters ( *belosList );
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 );
436 std::cout <<
" done ! " << std::endl;
441 std::cout <<
" -- Setting matrix in the solver ... " << std::flush;
443 linearSolver.setOperator ( systemMatrix );
447 std::cout <<
" done ! " << std::endl;
454 std::cout <<
" -- Defining the solution ... " << std::flush;
456 vectorPtr_Type solution (
new vector_Type (uFESpace->map(), Unique) );
460 std::cout <<
" done ! " << std::endl;
467 std::cout <<
" -- Solving the system ... " << std::flush;
469 linearSolver.setRightHandSide ( rhsBC );
470 linearSolver.solve ( solution );
474 std::cout <<
" done ! " << std::endl;
485 std::cout <<
" -- Computing the error ... " << std::flush;
487 vector_Type solutionErr (*solution);
489 uFESpace->interpolate (
static_cast<feSpace_Type::function_Type> ( exactSolution ), solutionErr, 0.0 );
490 solutionErr -= *solution;
492 Real l2error (uFESpace->l2Error (exactSolution, vector_Type (*solution, Repeated), 0.0) );
495 std::cout <<
" -- done ! " << std::endl;
499 std::cout <<
" ---> Norm L2 : " << l2error << std::endl;
501 Real linferror ( solutionErr.normInf() );
504 std::cout <<
" ---> Norm Inf : " << linferror << std::endl;
508 if (l2error > 0.0026)
510 std::cout <<
" <!> Solution has changed !!! <!> " << std::endl;
513 if (linferror > 0.000016)
515 std::cout <<
" <!> Solution has changed !!! <!> " << std::endl;
523 std::cout <<
" -- Defining the exporter ... " << std::flush;
525 ExporterEnsight<mesh_Type> exporter ( dataFile, meshPtr,
"solution", Comm->MyPID() ) ;
528 std::cout <<
" done ! " << std::endl;
533 std::cout <<
" -- Defining the exported quantities ... " << std::flush;
536 std::shared_ptr<vector_Type> solutionPtr (
new vector_Type (*solution, Repeated) );
537 std::shared_ptr<vector_Type> solutionErrPtr (
new vector_Type (solutionErr, Repeated) );
541 std::cout <<
" done ! " << std::endl;
546 std::cout <<
" -- Updating the exporter ... " << std::flush;
548 exporter.addVariable ( ExporterData<mesh_Type>::ScalarField,
"solution", uFESpace, solutionPtr, UInt (0) );
549 exporter.addVariable ( ExporterData<mesh_Type>::ScalarField,
"error", uFESpace, solutionErrPtr, UInt (0) );
552 std::cout <<
" done ! " << std::endl;
557 std::cout <<
" -- Exporting ... " << std::flush;
559 exporter.postProcess (0);
562 std::cout <<
" done ! " << std::endl;
567 std::cout <<
"End Result: TEST PASSED" << std::endl;
574 return ( EXIT_SUCCESS );
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.
std::shared_ptr< vector_Type > vectorPtr_Type
MatrixEpetra< Real > matrix_Type
static const LifeV::UInt elm_nodes_num[]
LifeV::PreconditionerIfpack prec_Type
FESpace< mesh_Type, MapEpetra > feSpace_Type
Epetra_Import const & importer()
Getter for the Epetra_Import.
RegionMesh< LinearTetra > mesh_Type
double Real
Generic real data.
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)