3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
linear/pdesolver.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_LINEAR_PDE_SOLVER_HH
25#define DUMUX_LINEAR_PDE_SOLVER_HH
26
27#include <cmath>
28#include <memory>
29#include <iostream>
30#include <type_traits>
31
32#include <dune/common/timer.hh>
33#include <dune/common/exceptions.hh>
34#include <dune/common/parallel/mpicommunication.hh>
35#include <dune/common/parallel/mpihelper.hh>
36#include <dune/common/std/type_traits.hh>
37#include <dune/istl/bvector.hh>
38#include <dune/istl/multitypeblockvector.hh>
39
47
48namespace Dumux {
49
60template <class Assembler, class LinearSolver>
61class LinearPDESolver : public PDESolver<Assembler, LinearSolver>
62{
64 using Scalar = typename Assembler::Scalar;
65 using JacobianMatrix = typename Assembler::JacobianMatrix;
66 using SolutionVector = typename Assembler::ResidualType;
68
69public:
73 LinearPDESolver(std::shared_ptr<Assembler> assembler,
74 std::shared_ptr<LinearSolver> linearSolver,
75 const std::string& paramGroup = "")
77 , paramGroup_(paramGroup)
78 , reuseMatrix_(false)
79 {
80 initParams_(paramGroup);
81
82 // set the linear system (matrix & residual) in the assembler
83 this->assembler().setLinearSystem();
84 }
85
89 void solve(SolutionVector& uCurrentIter) override
90 {
91 Dune::Timer assembleTimer(false);
92 Dune::Timer solveTimer(false);
93 Dune::Timer updateTimer(false);
94
95 if (verbose_) {
96 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
97 << std::flush;
98 }
99
101 // assemble
103
104 // linearize the problem at the current solution
105 assembleTimer.start();
106 if (reuseMatrix_)
107 this->assembler().assembleResidual(uCurrentIter);
108 else
109 this->assembler().assembleJacobianAndResidual(uCurrentIter);
110 assembleTimer.stop();
111
113 // linear solve
115
116 // Clear the current line using an ansi escape
117 // sequence. for an explanation see
118 // http://en.wikipedia.org/wiki/ANSI_escape_code
119 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
120
121 if (verbose_) {
122 std::cout << "\rSolve: M deltax^k = r";
123 std::cout << clearRemainingLine
124 << std::flush;
125 }
126
127 // solve the resulting linear equation system
128 solveTimer.start();
129
130 // set the delta vector to zero before solving the linear system!
131 SolutionVector deltaU(uCurrentIter);
132 deltaU = 0;
133
134 // solve by calling the appropriate implementation depending on whether the linear solver
135 // is capable of handling MultiType matrices or not
136 bool converged = solveLinearSystem_(deltaU);
137 solveTimer.stop();
138
139 if (!converged)
140 DUNE_THROW(NumericalProblem, "Linear solver didn't converge.\n");
141
143 // update
145 if (verbose_) {
146 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k";
147 std::cout << clearRemainingLine;
148 std::cout.flush();
149 }
150
151 // update the current solution and secondary variables
152 updateTimer.start();
153 uCurrentIter -= deltaU;
154 this->assembler().updateGridVariables(uCurrentIter);
155 updateTimer.stop();
156
157 if (verbose_) {
158 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
159 std::cout << "Assemble/solve/update time: "
160 << assembleTimer.elapsed() << "(" << 100*assembleTimer.elapsed()/elapsedTot << "%)/"
161 << solveTimer.elapsed() << "(" << 100*solveTimer.elapsed()/elapsedTot << "%)/"
162 << updateTimer.elapsed() << "(" << 100*updateTimer.elapsed()/elapsedTot << "%)"
163 << "\n";
164 }
165 }
166
171 void report(std::ostream& sout = std::cout) const
172 {}
173
178 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
179 {
180 return oldTimeStep;
181 }
182
186 void setVerbose(bool val)
187 { verbose_ = val; }
188
192 bool verbose() const
193 { return verbose_ ; }
194
198 const std::string& paramGroup() const
199 { return paramGroup_; }
200
207 void reuseMatrix(bool reuse = true)
208 { reuseMatrix_ = reuse; }
209
210private:
211
212 virtual bool solveLinearSystem_(SolutionVector& deltaU)
213 {
214 return solveLinearSystemImpl_(this->linearSolver(),
215 this->assembler().jacobian(),
216 deltaU,
217 this->assembler().residual());
218 }
219
229 template<class V = SolutionVector>
230 typename std::enable_if_t<!isMultiTypeBlockVector<V>(), bool>
231 solveLinearSystemImpl_(LinearSolver& ls,
232 JacobianMatrix& A,
233 SolutionVector& x,
234 SolutionVector& b)
235 {
242 constexpr auto blockSize = JacobianMatrix::block_type::rows;
243 using BlockType = Dune::FieldVector<Scalar, blockSize>;
244 Dune::BlockVector<BlockType> xTmp; xTmp.resize(b.size());
245 Dune::BlockVector<BlockType> bTmp(xTmp);
246 for (unsigned int i = 0; i < b.size(); ++i)
247 for (unsigned int j = 0; j < blockSize; ++j)
248 bTmp[i][j] = b[i][j];
249
250 const int converged = ls.solve(A, xTmp, bTmp);
251
252 for (unsigned int i = 0; i < x.size(); ++i)
253 for (unsigned int j = 0; j < blockSize; ++j)
254 x[i][j] = xTmp[i][j];
255
256 return converged;
257 }
258
259
269 template<class LS = LinearSolver, class V = SolutionVector>
270 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
271 isMultiTypeBlockVector<V>(), bool>
272 solveLinearSystemImpl_(LinearSolver& ls,
273 JacobianMatrix& A,
274 SolutionVector& x,
275 SolutionVector& b)
276 {
277 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
278 return ls.solve(A, x, b);
279 }
280
291 template<class LS = LinearSolver, class V = SolutionVector>
292 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
293 isMultiTypeBlockVector<V>(), bool>
294 solveLinearSystemImpl_(LinearSolver& ls,
295 JacobianMatrix& A,
296 SolutionVector& x,
297 SolutionVector& b)
298 {
299 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
300
301 // create the bcrs matrix the IterativeSolver backend can handle
303
304 // get the new matrix sizes
305 const std::size_t numRows = M.N();
306 assert(numRows == M.M());
307
308 // create the vector the IterativeSolver backend can handle
310 assert(bTmp.size() == numRows);
311
312 // create a blockvector to which the linear solver writes the solution
313 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
314 using BlockVector = typename Dune::BlockVector<VectorBlock>;
315 BlockVector y(numRows);
316
317 // solve
318 const bool converged = ls.solve(M, y, bTmp);
319
320 // copy back the result y into x
321 if(converged)
323
324 return converged;
325 }
326
328 void initParams_(const std::string& group = "")
329 {
330 verbose_ = (Dune::MPIHelper::getCommunication().rank() == 0);
331 }
332
334 bool verbose_;
335
337 std::string paramGroup_;
338
340 bool reuseMatrix_;
341};
342
343} // end namespace Dumux
344
345#endif
Type traits to be used with vector types.
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Trait checking if linear solvers accept Dune::MultiTypeBlockMatrix or we need to convert the matrix.
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:185
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:72
LinearSolver LinearSolver
Definition: common/pdesolver.hh:79
bool checkSizesOfSubMatrices(const Dune::MultiTypeBlockMatrix< FirstRow, Args... > &matrix) const
Helper function to assure the MultiTypeBlockMatrix's sub-blocks have the correct sizes.
Definition: common/pdesolver.hh:148
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:133
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:121
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:68
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:58
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copies the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:241
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:215
An implementation of a linear PDE solver.
Definition: linear/pdesolver.hh:62
const std::string & paramGroup() const
Returns the parameter group.
Definition: linear/pdesolver.hh:198
void setVerbose(bool val)
Specifies if the solver ought to be chatty.
Definition: linear/pdesolver.hh:186
LinearPDESolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const std::string &paramGroup="")
The Constructor.
Definition: linear/pdesolver.hh:73
void solve(SolutionVector &uCurrentIter) override
Solve a linear PDE system.
Definition: linear/pdesolver.hh:89
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: linear/pdesolver.hh:171
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: linear/pdesolver.hh:178
bool verbose() const
Returns true if the solver ought to be chatty.
Definition: linear/pdesolver.hh:192
void reuseMatrix(bool reuse=true)
Set whether the matrix should be reused.
Definition: linear/pdesolver.hh:207
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.