3.2-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
35#include <dune/common/version.hh>
36#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
37#include <dune/common/parallel/mpicollectivecommunication.hh>
38#else
39#include <dune/common/parallel/mpicommunication.hh>
40#endif
41
42#include <dune/common/parallel/mpihelper.hh>
43#include <dune/common/std/type_traits.hh>
44#include <dune/istl/bvector.hh>
45#include <dune/istl/multitypeblockvector.hh>
46
54
55namespace Dumux {
56
67template <class Assembler, class LinearSolver>
68class LinearPDESolver : public PDESolver<Assembler, LinearSolver>
69{
71 using Scalar = typename Assembler::Scalar;
72 using JacobianMatrix = typename Assembler::JacobianMatrix;
73 using SolutionVector = typename Assembler::ResidualType;
75
76public:
80 LinearPDESolver(std::shared_ptr<Assembler> assembler,
81 std::shared_ptr<LinearSolver> linearSolver,
82 const std::string& paramGroup = "")
84 , paramGroup_(paramGroup)
85 {
86 initParams_(paramGroup);
87
88 // set the linear system (matrix & residual) in the assembler
89 this->assembler().setLinearSystem();
90 }
91
95 void solve(SolutionVector& uCurrentIter) override
96 {
97 Dune::Timer assembleTimer(false);
98 Dune::Timer solveTimer(false);
99 Dune::Timer updateTimer(false);
100
101 if (verbose_) {
102 std::cout << "Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
103 << std::flush;
104 }
105
107 // assemble
109
110 // linearize the problem at the current solution
111 assembleTimer.start();
112 this->assembler().assembleJacobianAndResidual(uCurrentIter);
113 assembleTimer.stop();
114
116 // linear solve
118
119 // Clear the current line using an ansi escape
120 // sequence. for an explanation see
121 // http://en.wikipedia.org/wiki/ANSI_escape_code
122 const char clearRemainingLine[] = { 0x1b, '[', 'K', 0 };
123
124 if (verbose_) {
125 std::cout << "\rSolve: M deltax^k = r";
126 std::cout << clearRemainingLine
127 << std::flush;
128 }
129
130 // solve the resulting linear equation system
131 solveTimer.start();
132
133 // set the delta vector to zero before solving the linear system!
134 SolutionVector deltaU(uCurrentIter);
135 deltaU = 0;
136
137 // solve by calling the appropriate implementation depending on whether the linear solver
138 // is capable of handling MultiType matrices or not
139 bool converged = solveLinearSystem_(deltaU);
140 solveTimer.stop();
141
142 if (!converged)
143 DUNE_THROW(NumericalProblem, "Linear solver didn't converge.\n");
144
146 // update
148 if (verbose_) {
149 std::cout << "\rUpdate: x^(k+1) = x^k - deltax^k";
150 std::cout << clearRemainingLine;
151 std::cout.flush();
152 }
153
154 // update the current solution and secondary variables
155 updateTimer.start();
156 uCurrentIter -= deltaU;
157 this->assembler().updateGridVariables(uCurrentIter);
158 updateTimer.stop();
159
160 if (verbose_) {
161 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
162 std::cout << "Assemble/solve/update time: "
163 << assembleTimer.elapsed() << "(" << 100*assembleTimer.elapsed()/elapsedTot << "%)/"
164 << solveTimer.elapsed() << "(" << 100*solveTimer.elapsed()/elapsedTot << "%)/"
165 << updateTimer.elapsed() << "(" << 100*updateTimer.elapsed()/elapsedTot << "%)"
166 << "\n";
167 }
168 }
169
174 void report(std::ostream& sout = std::cout) const
175 {}
176
181 Scalar suggestTimeStepSize(Scalar oldTimeStep) const
182 {
183 return oldTimeStep;
184 }
185
189 void setVerbose(bool val)
190 { verbose_ = val; }
191
195 bool verbose() const
196 { return verbose_ ; }
197
201 const std::string& paramGroup() const
202 { return paramGroup_; }
203
204private:
205
206 virtual bool solveLinearSystem_(SolutionVector& deltaU)
207 {
208 return solveLinearSystemImpl_(this->linearSolver(),
209 this->assembler().jacobian(),
210 deltaU,
211 this->assembler().residual());
212 }
213
223 template<class V = SolutionVector>
224 typename std::enable_if_t<!isMultiTypeBlockVector<V>(), bool>
225 solveLinearSystemImpl_(LinearSolver& ls,
226 JacobianMatrix& A,
227 SolutionVector& x,
228 SolutionVector& b)
229 {
236 constexpr auto blockSize = JacobianMatrix::block_type::rows;
237 using BlockType = Dune::FieldVector<Scalar, blockSize>;
238 Dune::BlockVector<BlockType> xTmp; xTmp.resize(b.size());
239 Dune::BlockVector<BlockType> bTmp(xTmp);
240 for (unsigned int i = 0; i < b.size(); ++i)
241 for (unsigned int j = 0; j < blockSize; ++j)
242 bTmp[i][j] = b[i][j];
243
244 const int converged = ls.solve(A, xTmp, bTmp);
245
246 for (unsigned int i = 0; i < x.size(); ++i)
247 for (unsigned int j = 0; j < blockSize; ++j)
248 x[i][j] = xTmp[i][j];
249
250 return converged;
251 }
252
253
263 template<class LS = LinearSolver, class V = SolutionVector>
264 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
265 isMultiTypeBlockVector<V>(), bool>
266 solveLinearSystemImpl_(LinearSolver& ls,
267 JacobianMatrix& A,
268 SolutionVector& x,
269 SolutionVector& b)
270 {
271 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
272 return ls.solve(A, x, b);
273 }
274
285 template<class LS = LinearSolver, class V = SolutionVector>
286 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
287 isMultiTypeBlockVector<V>(), bool>
288 solveLinearSystemImpl_(LinearSolver& ls,
289 JacobianMatrix& A,
290 SolutionVector& x,
291 SolutionVector& b)
292 {
293 assert(this->checkSizesOfSubMatrices(A) && "Sub-blocks of MultiTypeBlockMatrix have wrong sizes!");
294
295 // create the bcrs matrix the IterativeSolver backend can handle
297
298 // get the new matrix sizes
299 const std::size_t numRows = M.N();
300 assert(numRows == M.M());
301
302 // create the vector the IterativeSolver backend can handle
304 assert(bTmp.size() == numRows);
305
306 // create a blockvector to which the linear solver writes the solution
307 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
308 using BlockVector = typename Dune::BlockVector<VectorBlock>;
309 BlockVector y(numRows);
310
311 // solve
312 const bool converged = ls.solve(M, y, bTmp);
313
314 // copy back the result y into x
315 if(converged)
317
318 return converged;
319 }
320
322 void initParams_(const std::string& group = "")
323 {
324 verbose_ = (Dune::MPIHelper::getCollectiveCommunication().rank() == 0);
325 }
326
328 bool verbose_;
329
331 std::string paramGroup_;
332};
333
334} // end namespace Dumux
335
336#endif
Type traits to be used with vector types.
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Manages the handling of time dependent problems.
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
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: timeloop.hh:72
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:243
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:69
const std::string & paramGroup() const
Returns the parameter group.
Definition: linear/pdesolver.hh:201
void setVerbose(bool val)
Specifies if the solver ought to be chatty.
Definition: linear/pdesolver.hh:189
LinearPDESolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const std::string &paramGroup="")
The Constructor.
Definition: linear/pdesolver.hh:80
void solve(SolutionVector &uCurrentIter) override
Solve a linear PDE system.
Definition: linear/pdesolver.hh:95
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: linear/pdesolver.hh:174
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: linear/pdesolver.hh:181
bool verbose() const
Returns true if the solver ought to be chatty.
Definition: linear/pdesolver.hh:195
Defines a high-level interface for a PDESolver.