36 #ifndef STOKESFLOWSOLVER_HPP_ 37 #define STOKESFLOWSOLVER_HPP_ 39 #include "AbstractContinuumMechanicsSolver.hpp" 40 #include "LinearSystem.hpp" 41 #include "LinearBasisFunction.hpp" 42 #include "QuadraticBasisFunction.hpp" 44 #include "PetscMatTools.hpp" 45 #include "PetscVecTools.hpp" 46 #include "StokesFlowProblemDefinition.hpp" 47 #include "StokesFlowAssembler.hpp" 48 #include "StokesFlowPreconditionerAssembler.hpp" 49 #include "ContinuumMechanicsNeumannBcsAssembler.hpp" 54 template<
unsigned DIM>
57 friend class TestStokesFlowSolver;
101 std::string outputDirectory);
140 template<
unsigned DIM>
143 std::string outputDirectory)
148 assert(DIM==2 || DIM==3);
156 template<
unsigned DIM>
164 template<
unsigned DIM>
196 KSPCreate(PETSC_COMM_WORLD,&solver);
197 #if ((PETSC_VERSION_MAJOR==3) && (PETSC_VERSION_MINOR>=5)) 202 KSPSetType(solver, KSPGMRES);
206 double ksp_rel_tol = 1e-6;
207 KSPSetTolerances(solver, ksp_rel_tol, PETSC_DEFAULT, PETSC_DEFAULT, 10000 );
211 KSPSetTolerances(solver, 1e-16,
mKspAbsoluteTol, PETSC_DEFAULT, 10000 );
213 unsigned num_restarts = 100;
214 KSPGMRESSetRestart(solver,num_restarts);
217 KSPGetPC(solver, &pc);
218 PCSetType(pc, PCJACOBI);
222 KSPSetFromOptions(solver);
251 KSPConvergedReason reason;
252 KSPGetConvergedReason(solver,&reason);
259 KSPGetIterationNumber(solver, &num_iters);
260 std::cout <<
"[" <<
PetscTools::GetMyRank() <<
"]: Num iterations = " << num_iters <<
"\n" << std::flush;
267 for (
unsigned i=0; i<this->
mNumDofs; i++)
284 template<
unsigned DIM>
317 template<
unsigned DIM>
320 assert(kspAbsoluteTolerance > 0);
324 template<
unsigned DIM>
330 for (
unsigned j=0; j<DIM; j++)
339 template<
unsigned DIM>
std::vector< c_vector< double, DIM > > & rGetVelocities()
void ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem)
void WriteCurrentPressureSolution(int counterToAppend=-1)
Vec mLinearSystemRhsVector
AbstractTetrahedralMesh< DIM, DIM > & mrQuadMesh
StokesFlowAssembler< DIM > * mpStokesFlowAssembler
void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1)
static void BeginEvent(unsigned event)
std::vector< c_vector< double, DIM > > mSpatialSolution
virtual unsigned GetNumNodes() const
static void PrintAndReset(std::string message)
void AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type)
StokesFlowPreconditionerAssembler< DIM > * mpStokesFlowPreconditionerAssembler
void RemovePressureDummyValuesThroughLinearInterpolation()
virtual ~StokesFlowSolver()
std::vector< double > mCurrentSolution
StokesFlowSolver(AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, StokesFlowProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory)
ContinuumMechanicsNeumannBcsAssembler< DIM > * mpNeumannBcsAssembler
static void EndEvent(unsigned event)
void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
StokesFlowProblemDefinition< DIM > & mrProblemDefinition
std::vector< c_vector< double, DIM > > & rGetSpatialSolution()