3.3.0
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 {
79 initParams_(paramGroup);
80
81 // set the linear system (matrix & residual) in the assembler
82 this->assembler().setLinearSystem();
83 }
84
88 void solve(SolutionVector& uCurrentIter) override
89 {
90 Dune::Timer assembleTimer(false);
91 Dune::Timer solveTimer(false);
92 Dune::Timer updateTimer(false);
93
94 if (verbose_) {
95 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
96 << std::flush;
97 }
98
100 // assemble
102
103 // linearize the problem at the current solution
104 assembleTimer.start();
105 this->assembler().assembleJacobianAndResidual(uCurrentIter);
106 assembleTimer.stop();
107
109 // linear solve
111
112 // Clear the current line using an ansi escape
113 // sequence. for an explanation see
114 // http://en.wikipedia.org/wiki/ANSI_escape_code
115 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
116
117 if (verbose_) {
118 std::cout << "\rSolve: M deltax^k = r";
119 std::cout << clearRemainingLine
120 << std::flush;
121 }
122
123 // solve the resulting linear equation system
124 solveTimer.start();
125
126 // set the delta vector to zero before solving the linear system!
127 SolutionVector deltaU(uCurrentIter);
128 deltaU = 0;
129
130 // solve by calling the appropriate implementation depending on whether the linear solver
131 // is capable of handling MultiType matrices or not
132 bool converged = solveLinearSystem_(deltaU);
133 solveTimer.stop();
134
135 if (!converged)
136 DUNE_THROW(NumericalProblem, "Linear solver didn't converge.\n");
137
139 // update
141 if (verbose_) {
142 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k";
143 std::cout << clearRemainingLine;
144 std::cout.flush();
145 }
146
147 // update the current solution and secondary variables
148 updateTimer.start();
149 uCurrentIter -= deltaU;
150 this->assembler().updateGridVariables(uCurrentIter);
151 updateTimer.stop();
152
153 if (verbose_) {
154 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
155 std::cout << "Assemble/solve/update time: "
156 << assembleTimer.elapsed() << "(" << 100*assembleTimer.elapsed()/elapsedTot << "%)/"
157 << solveTimer.elapsed() << "(" << 100*solveTimer.elapsed()/elapsedTot << "%)/"
158 << updateTimer.elapsed() << "(" << 100*updateTimer.elapsed()/elapsedTot << "%)"
159 << "\n";
160 }
161 }
162
167 void report(std::ostream& sout = std::cout) const
168 {}
169
174 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
175 {
176 return oldTimeStep;
177 }
178
182 void setVerbose(bool val)
183 { verbose_ = val; }
184
188 bool verbose() const
189 { return verbose_ ; }
190
194 const std::string& paramGroup() const
195 { return paramGroup_; }
196
197private:
198
199 virtual bool solveLinearSystem_(SolutionVector& deltaU)
200 {
201 return solveLinearSystemImpl_(this->linearSolver(),
202 this->assembler().jacobian(),
203 deltaU,
204 this->assembler().residual());
205 }
206
216 template<class V = SolutionVector>
217 typename std::enable_if_t<!isMultiTypeBlockVector<V>(), bool>
218 solveLinearSystemImpl_(LinearSolver& ls,
219 JacobianMatrix& A,
220 SolutionVector& x,
221 SolutionVector& b)
222 {
229 constexpr auto blockSize = JacobianMatrix::block_type::rows;
230 using BlockType = Dune::FieldVector<Scalar, blockSize>;
231 Dune::BlockVector<BlockType> xTmp; xTmp.resize(b.size());
232 Dune::BlockVector<BlockType> bTmp(xTmp);
233 for (unsigned int i = 0; i < b.size(); ++i)
234 for (unsigned int j = 0; j < blockSize; ++j)
235 bTmp[i][j] = b[i][j];
236
237 const int converged = ls.solve(A, xTmp, bTmp);
238
239 for (unsigned int i = 0; i < x.size(); ++i)
240 for (unsigned int j = 0; j < blockSize; ++j)
241 x[i][j] = xTmp[i][j];
242
243 return converged;
244 }
245
246
256 template<class LS = LinearSolver, class V = SolutionVector>
257 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
258 isMultiTypeBlockVector<V>(), bool>
259 solveLinearSystemImpl_(LinearSolver& ls,
260 JacobianMatrix& A,
261 SolutionVector& x,
262 SolutionVector& b)
263 {
264 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
265 return ls.solve(A, x, b);
266 }
267
278 template<class LS = LinearSolver, class V = SolutionVector>
279 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
280 isMultiTypeBlockVector<V>(), bool>
281 solveLinearSystemImpl_(LinearSolver& ls,
282 JacobianMatrix& A,
283 SolutionVector& x,
284 SolutionVector& b)
285 {
286 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
287
288 // create the bcrs matrix the IterativeSolver backend can handle
290
291 // get the new matrix sizes
292 const std::size_t numRows = M.N();
293 assert(numRows == M.M());
294
295 // create the vector the IterativeSolver backend can handle
297 assert(bTmp.size() == numRows);
298
299 // create a blockvector to which the linear solver writes the solution
300 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
301 using BlockVector = typename Dune::BlockVector<VectorBlock>;
302 BlockVector y(numRows);
303
304 // solve
305 const bool converged = ls.solve(M, y, bTmp);
306
307 // copy back the result y into x
308 if(converged)
310
311 return converged;
312 }
313
315 void initParams_(const std::string& group = "")
316 {
317 verbose_ = (Dune::MPIHelper::getCollectiveCommunication().rank() == 0);
318 }
319
321 bool verbose_;
322
324 std::string paramGroup_;
325};
326
327} // end namespace Dumux
328
329#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 classe that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Definition: adapt.hh:29
constexpr std::size_t blockSize()
Definition: nonlinear/newtonsolver.hh:156
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:55
const Assembler & assembler() const
Access the assembler.
Definition: common/pdesolver.hh:92
const LinearSolver & linearSolver() const
Access the linear solver.
Definition: common/pdesolver.hh:104
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:117
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
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)
Copys 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:194
void setVerbose(bool val)
Specifies if the solver ought to be chatty.
Definition: linear/pdesolver.hh:182
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:88
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: linear/pdesolver.hh:167
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: linear/pdesolver.hh:174
bool verbose() const
Returns true if the solver ought to be chatty.
Definition: linear/pdesolver.hh:188
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.