24#ifndef DUMUX_LINEAR_PDE_SOLVER_HH
25#define DUMUX_LINEAR_PDE_SOLVER_HH
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>
60template <
class Assembler,
class LinearSolver>
64 using Scalar =
typename Assembler::Scalar;
65 using JacobianMatrix =
typename Assembler::JacobianMatrix;
66 using SolutionVector =
typename Assembler::ResidualType;
89 void solve(SolutionVector& uCurrentIter)
override
91 Dune::Timer assembleTimer(
false);
92 Dune::Timer solveTimer(
false);
93 Dune::Timer updateTimer(
false);
96 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
105 assembleTimer.start();
107 this->
assembler().assembleResidual(uCurrentIter);
109 this->
assembler().assembleJacobianAndResidual(uCurrentIter);
110 assembleTimer.stop();
119 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
122 std::cout <<
"\rSolve: M deltax^k = r";
123 std::cout << clearRemainingLine
131 SolutionVector deltaU(uCurrentIter);
136 bool converged = solveLinearSystem_(deltaU);
146 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k";
147 std::cout << clearRemainingLine;
153 uCurrentIter -= deltaU;
154 this->
assembler().updateGridVariables(uCurrentIter);
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 <<
"%)"
171 void report(std::ostream& sout = std::cout)
const
193 {
return verbose_ ; }
199 {
return paramGroup_; }
208 { reuseMatrix_ = reuse; }
212 virtual bool solveLinearSystem_(SolutionVector& deltaU)
229 template<
class V = SolutionVector>
230 typename std::enable_if_t<!isMultiTypeBlockVector<V>(),
bool>
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];
250 const int converged = ls.solve(A, xTmp, bTmp);
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];
269 template<
class LS = LinearSolver,
class V = SolutionVector>
270 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
271 isMultiTypeBlockVector<V>(),
bool>
278 return ls.solve(A, x, b);
291 template<
class LS = LinearSolver,
class V = SolutionVector>
292 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
293 isMultiTypeBlockVector<V>(),
bool>
305 const std::size_t numRows = M.N();
306 assert(numRows == M.M());
310 assert(bTmp.size() == numRows);
313 using VectorBlock =
typename Dune::FieldVector<Scalar, 1>;
314 using BlockVector =
typename Dune::BlockVector<VectorBlock>;
315 BlockVector y(numRows);
318 const bool converged = ls.solve(M, y, bTmp);
328 void initParams_(
const std::string& group =
"")
330 verbose_ = (Dune::MPIHelper::getCommunication().rank() == 0);
337 std::string paramGroup_;
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.
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: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: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 ¶mGroup="")
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.