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)
145 mrProblemDefinition(rProblemDefinition),
148 assert(DIM==2 || DIM==3);
156 template<
unsigned DIM>
159 delete mpStokesFlowAssembler;
160 delete mpStokesFlowPreconditionerAssembler;
161 delete mpNeumannBcsAssembler;
164 template<
unsigned DIM>
167 mrProblemDefinition.Validate();
192 VecDuplicate(this->mLinearSystemRhsVector,&solution);
196 KSPCreate(PETSC_COMM_WORLD,&solver);
197 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
198 KSPSetOperators(solver, this->mSystemLhsMatrix, this->mPreconditionMatrix);
200 KSPSetOperators(solver, this->mSystemLhsMatrix, this->mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN );
202 KSPSetType(solver, KSPGMRES);
204 if (mKspAbsoluteTol < 0)
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);
249 KSPSolve(solver,this->mLinearSystemRhsVector,solution);
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++)
269 this->mCurrentSolution[i] = solution_repl[i];
275 this->RemovePressureDummyValuesThroughLinearInterpolation();
280 this->WriteCurrentSpatialSolution(
"flow_solution",
"nodes");
281 this->WriteCurrentPressureSolution();
286 template<
unsigned DIM>
290 mpStokesFlowAssembler->SetMatrixToAssemble(this->mSystemLhsMatrix,
true);
291 mpStokesFlowAssembler->SetVectorToAssemble(this->mLinearSystemRhsVector,
true);
292 mpStokesFlowAssembler->Assemble();
294 mpStokesFlowPreconditionerAssembler->SetMatrixToAssemble(this->mPreconditionMatrix,
true);
295 mpStokesFlowPreconditionerAssembler->AssembleMatrix();
297 mpNeumannBcsAssembler->SetVectorToAssemble(this->mLinearSystemRhsVector,
false );
298 mpNeumannBcsAssembler->AssembleVector();
311 this->AddIdentityBlockForDummyPressureVariables(LINEAR_PROBLEM);
312 this->ApplyDirichletBoundaryConditions(LINEAR_PROBLEM,
true);
319 template<
unsigned DIM>
322 assert(kspAbsoluteTolerance > 0);
323 mKspAbsoluteTol = kspAbsoluteTolerance;
326 template<
unsigned DIM>
329 this->mSpatialSolution.resize(this->mrQuadMesh.GetNumNodes(), zero_vector<double>(DIM));
330 for (
unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
332 for (
unsigned j=0; j<DIM; j++)
335 this->mSpatialSolution[i](j) = this->mCurrentSolution[(DIM+1)*i+j];
338 return this->mSpatialSolution;
341 template<
unsigned DIM>
344 return rGetSpatialSolution();
std::vector< c_vector< double, DIM > > & rGetVelocities()
AbstractTetrahedralMesh< DIM, DIM > & mrQuadMesh
StokesFlowAssembler< DIM > * mpStokesFlowAssembler
static void BeginEvent(unsigned event)
static void PrintAndReset(std::string message)
StokesFlowPreconditionerAssembler< DIM > * mpStokesFlowPreconditionerAssembler
virtual ~StokesFlowSolver()
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()