Chaste Commit::1fd4e48e3990e67db148bc1bc4cf6991a0049d0c
StokesFlowSolver.hpp
1/*
2
3Copyright (c) 2005-2024, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, are permitted provided that the following conditions are met:
14 * Redistributions of source code must retain the above copyright notice,
15 this list of conditions and the following disclaimer.
16 * Redistributions in binary form must reproduce the above copyright notice,
17 this list of conditions and the following disclaimer in the documentation
18 and/or other materials provided with the distribution.
19 * Neither the name of the University of Oxford nor the names of its
20 contributors may be used to endorse or promote products derived from this
21 software without specific prior written permission.
22
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#ifndef STOKESFLOWSOLVER_HPP_
37#define STOKESFLOWSOLVER_HPP_
38
39#include "AbstractContinuumMechanicsSolver.hpp"
40#include "LinearSystem.hpp"
41#include "LinearBasisFunction.hpp"
42#include "QuadraticBasisFunction.hpp"
43#include "Timer.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"
50
54template<unsigned DIM>
56{
57 friend class TestStokesFlowSolver;
58
59private:
62
65
71
77
84
88 void AssembleSystem();
89
90public:
91
100 StokesFlowProblemDefinition<DIM>& rProblemDefinition,
101 std::string outputDirectory);
102
106 virtual ~StokesFlowSolver();
107
111 void Solve();
112
119 void SetKspAbsoluteTolerance(double kspAbsoluteTolerance);
120
121
126 std::vector<c_vector<double,DIM> >& rGetSpatialSolution();
127
132 std::vector<c_vector<double,DIM> >& rGetVelocities();
133};
134
136// Implementation
138
139
140template<unsigned DIM>
142 StokesFlowProblemDefinition<DIM>& rProblemDefinition,
143 std::string outputDirectory)
144 : AbstractContinuumMechanicsSolver<DIM>(rQuadMesh, rProblemDefinition, outputDirectory, INCOMPRESSIBLE),
145 mrProblemDefinition(rProblemDefinition),
146 mKspAbsoluteTol(-1)
147{
148 assert(DIM==2 || DIM==3);
149 assert(!mrProblemDefinition.rGetDirichletNodes().empty());
150
154}
155
156template<unsigned DIM>
158{
159 delete mpStokesFlowAssembler;
160 delete mpStokesFlowPreconditionerAssembler;
161 delete mpNeumannBcsAssembler;
162}
163
164template<unsigned DIM>
166{
167 mrProblemDefinition.Validate();
168
169 if (this->mVerbose)
170 {
171 Timer::Reset();
172 }
173
174 // Assemble Jacobian (and preconditioner)
175 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
176 AssembleSystem();
177 MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
178
179 if (this->mVerbose)
180 {
181 Timer::PrintAndReset("AssembleSystem");
182 }
183
184 /*
185 * Solve the linear system using PETSc GMRES and an LU factorisation
186 * of the preconditioner. Note we don't call Solve on the linear_system
187 * as we want to set PETSc options.
188 */
189 MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
190
191 Vec solution;
192 VecDuplicate(this->mLinearSystemRhsVector,&solution);
193 PetscVecTools::Zero(solution);
194
195 KSP solver;
196 KSPCreate(PETSC_COMM_WORLD,&solver);
197#if ((PETSC_VERSION_MAJOR==3) && (PETSC_VERSION_MINOR>=5))
198 KSPSetOperators(solver, this->mSystemLhsMatrix, this->mPreconditionMatrix);
199#else
200 KSPSetOperators(solver, this->mSystemLhsMatrix, this->mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
201#endif
202 KSPSetType(solver, KSPGMRES);
203
204 if (mKspAbsoluteTol < 0)
205 {
206 double ksp_rel_tol = 1e-6;
207 KSPSetTolerances(solver, ksp_rel_tol, PETSC_DEFAULT, PETSC_DEFAULT, 10000 /*max iter*/); //hopefully with the preconditioner this max is way too high
208 }
209 else
210 {
211 KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, 10000 /*max iter*/); //hopefully with the preconditioner this max is way too high
212 }
213 unsigned num_restarts = 100;
214 KSPGMRESSetRestart(solver,num_restarts); // gmres num restarts
215
216 PC pc;
217 KSPGetPC(solver, &pc);
218 PCSetType(pc, PCJACOBI);
219
220 KSPSetUp(solver);
221
222 KSPSetFromOptions(solver);
223
224// ///// For printing matrix when debugging
225// OutputFileHandler handler("TEMP");
226// out_stream p_file = handler.OpenOutputFile("matrix.txt");
227// for (unsigned i=0; i<this->mNumDofs; i++)
228// {
229// for (unsigned j=0; j<this->mNumDofs; j++)
230// {
231// *p_file << PetscMatTools::GetElement(this->mSystemLhsMatrix, i, j) << " ";
232// }
233// *p_file << "\n";
234// }
235// p_file->close();
236//
237// out_stream p_file2 = handler.OpenOutputFile("rhs.txt");
238// for (unsigned i=0; i<this->mNumDofs; i++)
239// {
240// *p_file2 << PetscVecTools::GetElement(this->mLinearSystemRhsVector, i) << "\n";
241// }
242// p_file2->close();
243
244 if (this->mVerbose)
245 {
246 Timer::PrintAndReset("KSP Setup");
247 }
248
249 KSPSolve(solver,this->mLinearSystemRhsVector,solution);
250
251 KSPConvergedReason reason;
252 KSPGetConvergedReason(solver,&reason);
253 KSPEXCEPT(reason);
254
255 if (this->mVerbose)
256 {
257 Timer::PrintAndReset("KSP Solve");
258 int num_iters;
259 KSPGetIterationNumber(solver, &num_iters);
260 std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
261 }
262
263 MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
264
265 // Copy solution into the std::vector
266 ReplicatableVector solution_repl(solution);
267 for (unsigned i=0; i<this->mNumDofs; i++)
268 {
269 this->mCurrentSolution[i] = solution_repl[i];
270 }
271
272 // Remove pressure dummy values (P=0 at internal nodes, which should have been
273 // been the result of the solve above), by linear interpolating from vertices of
274 // edges to the internal node
275 this->RemovePressureDummyValuesThroughLinearInterpolation();
276
277 PetscTools::Destroy(solution);
278 KSPDestroy(PETSC_DESTROY_PARAM(solver));
279
280 this->WriteCurrentSpatialSolution("flow_solution", "nodes");
281 this->WriteCurrentPressureSolution();
282}
283
284template<unsigned DIM>
286{
287 // Use assembler to assemble volume integral part....
288 mpStokesFlowAssembler->SetMatrixToAssemble(this->mSystemLhsMatrix, true);
289 mpStokesFlowAssembler->SetVectorToAssemble(this->mLinearSystemRhsVector, true);
290 mpStokesFlowAssembler->Assemble();
291
292 mpStokesFlowPreconditionerAssembler->SetMatrixToAssemble(this->mPreconditionMatrix, true);
293 mpStokesFlowPreconditionerAssembler->AssembleMatrix();
294
295 mpNeumannBcsAssembler->SetVectorToAssemble(this->mLinearSystemRhsVector, false /*don't zero!*/);
296 mpNeumannBcsAssembler->AssembleVector();
297
298 PetscVecTools::Finalise(this->mLinearSystemRhsVector);
299 PetscMatTools::SwitchWriteMode(this->mSystemLhsMatrix);
300 PetscMatTools::SwitchWriteMode(this->mPreconditionMatrix);
301
302 // Note: maintaining symmetry for Dirichlet BCs is possible at the moment (the second parameter)
303 // but doing so for the identity block is not yet implemented (the second parameter must be false)
304 // Not sure if maintaining symmetry is worth it - may allow CG to work, but matrix is indefinite
305 // so GC not guaranteed to work..
306 //
307 // Note: the identity block needs to be added before the BCs - see comments in
308 // PetscMatTools::ZeroRowsWithValueOnDiagonal()
309 this->AddIdentityBlockForDummyPressureVariables(LINEAR_PROBLEM);
310 this->ApplyDirichletBoundaryConditions(LINEAR_PROBLEM, true);
311
312 PetscVecTools::Finalise(this->mLinearSystemRhsVector);
313 PetscMatTools::Finalise(this->mSystemLhsMatrix);
314 PetscMatTools::Finalise(this->mPreconditionMatrix);
315}
316
317template<unsigned DIM>
318void StokesFlowSolver<DIM>::SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
319{
320 assert(kspAbsoluteTolerance > 0);
321 mKspAbsoluteTol = kspAbsoluteTolerance;
322}
323
324template<unsigned DIM>
325std::vector<c_vector<double,DIM> >& StokesFlowSolver<DIM>::rGetSpatialSolution()
326{
327 this->mSpatialSolution.resize(this->mrQuadMesh.GetNumNodes(), zero_vector<double>(DIM));
328 for (unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
329 {
330 for (unsigned j=0; j<DIM; j++)
331 {
332 // DIM+1 is the problem dimension
333 this->mSpatialSolution[i](j) = this->mCurrentSolution[(DIM+1)*i+j];
334 }
335 }
336 return this->mSpatialSolution;
337}
338
339template<unsigned DIM>
340std::vector<c_vector<double,DIM> >& StokesFlowSolver<DIM>::rGetVelocities()
341{
342 return rGetSpatialSolution();
343}
344
345#endif /* STOKESFLOWSOLVER_HPP_ */
#define PETSC_DESTROY_PARAM(x)
AbstractTetrahedralMesh< DIM, DIM > & mrQuadMesh
static void SwitchWriteMode(Mat matrix)
static void Finalise(Mat matrix)
static void Destroy(Vec &rVec)
static unsigned GetMyRank()
static void Finalise(Vec vector)
static void Zero(Vec vector)
StokesFlowSolver(AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, StokesFlowProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory)
StokesFlowProblemDefinition< DIM > & mrProblemDefinition
ContinuumMechanicsNeumannBcsAssembler< DIM > * mpNeumannBcsAssembler
StokesFlowAssembler< DIM > * mpStokesFlowAssembler
StokesFlowPreconditionerAssembler< DIM > * mpStokesFlowPreconditionerAssembler
std::vector< c_vector< double, DIM > > & rGetSpatialSolution()
std::vector< c_vector< double, DIM > > & rGetVelocities()
void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
static void PrintAndReset(std::string message)
Definition Timer.cpp:70
static void Reset()
Definition Timer.cpp:44