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>
35#include <dune/common/version.hh>
36#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
37#include <dune/common/parallel/mpicollectivecommunication.hh>
39#include <dune/common/parallel/mpicommunication.hh>
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>
67template <
class Assembler,
class LinearSolver>
71 using Scalar =
typename Assembler::Scalar;
72 using JacobianMatrix =
typename Assembler::JacobianMatrix;
73 using SolutionVector =
typename Assembler::ResidualType;
95 void solve(SolutionVector& uCurrentIter)
override
97 Dune::Timer assembleTimer(
false);
98 Dune::Timer solveTimer(
false);
99 Dune::Timer updateTimer(
false);
102 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
111 assembleTimer.start();
112 this->
assembler().assembleJacobianAndResidual(uCurrentIter);
113 assembleTimer.stop();
122 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
125 std::cout <<
"\rSolve: M deltax^k = r";
126 std::cout << clearRemainingLine
134 SolutionVector deltaU(uCurrentIter);
139 bool converged = solveLinearSystem_(deltaU);
149 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k";
150 std::cout << clearRemainingLine;
156 uCurrentIter -= deltaU;
157 this->
assembler().updateGridVariables(uCurrentIter);
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 <<
"%)"
174 void report(std::ostream& sout = std::cout)
const
196 {
return verbose_ ; }
202 {
return paramGroup_; }
206 virtual bool solveLinearSystem_(SolutionVector& deltaU)
223 template<
class V = SolutionVector>
224 typename std::enable_if_t<!isMultiTypeBlockVector<V>(),
bool>
225 solveLinearSystemImpl_(LinearSolver& ls,
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];
244 const int converged = ls.solve(A, xTmp, bTmp);
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];
263 template<
class LS = LinearSolver,
class V = SolutionVector>
264 typename std::enable_if_t<linearSolverAcceptsMultiTypeMatrix<LS>() &&
265 isMultiTypeBlockVector<V>(),
bool>
266 solveLinearSystemImpl_(LinearSolver& ls,
272 return ls.solve(A, x, b);
285 template<
class LS = LinearSolver,
class V = SolutionVector>
286 typename std::enable_if_t<!linearSolverAcceptsMultiTypeMatrix<LS>() &&
287 isMultiTypeBlockVector<V>(),
bool>
288 solveLinearSystemImpl_(LinearSolver& ls,
299 const std::size_t numRows = M.N();
300 assert(numRows == M.M());
304 assert(bTmp.size() == numRows);
307 using VectorBlock =
typename Dune::FieldVector<Scalar, 1>;
308 using BlockVector =
typename Dune::BlockVector<VectorBlock>;
309 BlockVector y(numRows);
312 const bool converged = ls.solve(M, y, bTmp);
322 void initParams_(
const std::string& group =
"")
324 verbose_ = (Dune::MPIHelper::getCollectiveCommunication().rank() == 0);
331 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.
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.
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 ¶mGroup="")
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.