37 #ifndef ASSEMBLYPOLICYNAVIERSTOKESNEWTON_HPP 38 #define ASSEMBLYPOLICYNAVIERSTOKESNEWTON_HPP 42 #include <boost/shared_ptr.hpp> 46 #include <Epetra_MpiComm.h> 48 #include <Epetra_SerialComm.h> 51 #include <Teuchos_ParameterList.hpp> 52 #include <Teuchos_XMLParameterListHelpers.hpp> 53 #include <Teuchos_RCP.hpp> 56 #include <lifev/core/LifeV.hpp> 57 #include <lifev/core/array/MatrixEpetra.hpp> 58 #include <lifev/core/array/VectorEpetra.hpp> 59 #include <lifev/core/util/Displayer.hpp> 60 #include <lifev/core/mesh/RegionMesh.hpp> 61 #include <lifev/core/fem/FESpace.hpp> 62 #include <lifev/core/fem/TimeAdvanceBDF.hpp> 63 #include <lifev/core/algorithm/Preconditioner.hpp> 64 #include <lifev/core/util/LifeChrono.hpp> 66 #include <lifev/navier_stokes/algorithm/PreconditionerPCD.hpp> 67 #include <lifev/navier_stokes/solver/NavierStokesSolver/AssemblyPolicyStokes.hpp> 68 #include <lifev/navier_stokes/solver/NavierStokesSolver/NavierStokesProblem.hpp> 74 template<
class mesh_Type >
111 template<
class mesh_Type >
115 AssemblyPolicyStokes< mesh_Type >::initAssembly ( list );
120 displayer().leaderPrint (
"Creating the mass matrix... " );
122 M_massMatrix.reset (
new matrix_Type ( solutionMap ) );
123 M_massMatrix->zero();
124 AssemblyPolicyStokes< mesh_Type >::M_assembler->addMass ( *M_massMatrix, 1.0 );
125 M_massMatrix->globalAssemble();
129 displayer().leaderPrintMax (
"Matrices assembly time: ", assemblyChrono
.diff(),
" s.\n");
132 template<
class mesh_Type >
139 AssemblyPolicyStokes< mesh_Type >::assembleSystem ( systemMatrix, rhs, solution, preconditioner );
142 *rhs += *M_massMatrix * bdf()->rhsContributionFirstDerivative();
145 *systemMatrix += *M_massMatrix * alpha;
147 vector_Type beta ( systemMatrix->map(), Repeated );
149 AssemblyPolicyStokes< mesh_Type >::M_assembler->addConvection ( *systemMatrix, 1.0, beta );
150 AssemblyPolicyStokes< mesh_Type >::M_assembler->addNewtonConvection ( *systemMatrix, beta );
151 AssemblyPolicyStokes< mesh_Type >::M_assembler->addConvectionRhs ( *rhs, 1.0, beta );
153 if ( preconditioner->preconditionerType() ==
"PCD" )
155 PreconditionerPCD* pcdPtr =
dynamic_cast<PreconditionerPCD*> ( preconditioner.get() );
virtual Real timestep() const =0
VectorEpetra - The Epetra Vector format Wrapper.
std::shared_ptr< NavierStokesProblem< mesh_Type > > NSProblemPtr_Type
MeshPartitioner< mesh_Type > meshPartitioner_Type
void initAssembly(Teuchos::ParameterList &list)
PreconditionerPCD(std::shared_ptr< Epetra_Comm > comm=std::shared_ptr< Epetra_Comm >(new Epetra_MpiComm(MPI_COMM_WORLD)))
default constructor
std::shared_ptr< VectorEpetra > vectorPtr_Type
std::shared_ptr< map_Type > mapPtr_Type
void start()
Start the timer.
TimeAdvanceBDF< vector_Type > bdf_Type
matrixPtr_Type M_massMatrix
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
void updateInverseJacobian(const UInt &iQuadPt)
virtual fespacePtr_Type pFESpace() const =0
Epetra_Import const & importer()
Getter for the Epetra_Import.
std::shared_ptr< matrix_Type > matrixPtr_Type
Preconditioner preconditioner_Type
Real diff()
Compute the difference in time between start and stop.
FESpace< mesh_Type, map_Type > fespace_Type
std::shared_ptr< bdf_Type > bdfPtr_Type
std::shared_ptr< fespace_Type > fespacePtr_Type
double Real
Generic real data.
Preconditioner - Abstract preconditioner class.
void stop()
Stop the timer.
virtual bdfPtr_Type bdf() const =0
virtual fespacePtr_Type uFESpace() const =0
virtual Displayer displayer()=0
std::shared_ptr< preconditioner_Type > preconditionerPtr_Type
void updateBeta(const vector_Type &beta)
Update the vector beta of the convective term in Fp.
virtual NSProblemPtr_Type problem() const =0
Displayer - This class is used to display messages in parallel simulations.
MatrixEpetra< Real > matrix_Type
void assembleSystem(matrixPtr_Type systemMatrix, vectorPtr_Type rhs, vectorPtr_Type solution, preconditionerPtr_Type preconditioner)