12#ifndef DUMUX_LINEAR_PDE_SOLVER_HH
13#define DUMUX_LINEAR_PDE_SOLVER_HH
20#include <dune/common/timer.hh>
21#include <dune/common/exceptions.hh>
22#include <dune/common/parallel/mpicommunication.hh>
23#include <dune/common/parallel/mpihelper.hh>
24#include <dune/common/std/type_traits.hh>
25#include <dune/istl/bvector.hh>
26#include <dune/istl/multitypeblockvector.hh>
39template <
class Solver,
class Matrix>
40using SetMatrixDetector =
decltype(std::declval<Solver>().setMatrix(std::declval<std::shared_ptr<Matrix>>()));
42template<
class Solver,
class Matrix>
44{
return Dune::Std::is_detected<SetMatrixDetector, Solver, Matrix>::value; }
61 class Comm = Dune::Communication<Dune::MPIHelper::MPICommunicator>>
65 using Scalar =
typename Assembler::Scalar;
66 using JacobianMatrix =
typename Assembler::JacobianMatrix;
67 using SolutionVector =
typename Assembler::SolutionVector;
68 using ResidualVector =
typename Assembler::ResidualType;
72 static constexpr bool assemblerExportsVariables = Detail::PDESolver::assemblerExportsVariables<Assembler>;
83 const Communication& comm = Dune::MPIHelper::getCommunication(),
110 Dune::Timer assembleTimer(
false);
111 Dune::Timer solveTimer(
false);
112 Dune::Timer updateTimer(
false);
114 if (verbosity_ >= 1 && enableDynamicOutput_)
115 std::cout <<
"Assemble: r(x^k) = dS/dt + div F - q; M = grad r"
123 assembleTimer.start();
125 this->
assembler().assembleResidual(vars);
127 this->
assembler().assembleJacobianAndResidual(vars);
128 assembleTimer.stop();
137 const char clearRemainingLine[] = { 0x1b,
'[',
'K', 0 };
139 if (verbosity_ >= 1 && enableDynamicOutput_)
140 std::cout <<
"\rSolve: M deltax^k = r"
141 << clearRemainingLine << std::flush;
147 ResidualVector deltaU = LinearAlgebraNativeBackend::zeros(Backend::size(Backend::dofs(vars)));
151 bool converged = solveLinearSystem_(deltaU);
160 if (verbosity_ >= 1 && enableDynamicOutput_)
161 std::cout <<
"\rUpdate: x^(k+1) = x^k - deltax^k"
162 << clearRemainingLine << std::flush;
166 auto uCurrent = Backend::dofs(vars);
167 Backend::axpy(-1.0, deltaU, uCurrent);
168 Backend::update(vars, uCurrent);
169 if constexpr (!assemblerExportsVariables)
170 this->
assembler().updateGridVariables(Backend::dofs(vars));
175 const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed();
176 if (enableDynamicOutput_)
178 std::cout <<
"Assemble/solve/update time: "
179 << assembleTimer.elapsed() <<
"(" << 100*assembleTimer.elapsed()/elapsedTot <<
"%)/"
180 << solveTimer.elapsed() <<
"(" << 100*solveTimer.elapsed()/elapsedTot <<
"%)/"
181 << updateTimer.elapsed() <<
"(" << 100*updateTimer.elapsed()/elapsedTot <<
"%)"
193 bool converged =
apply(vars);
202 void report(std::ostream& sout = std::cout)
const
218 { verbosity_ = val; }
224 {
return verbosity_ ; }
230 {
return paramGroup_; }
240 reuseMatrix_ = reuse;
242 if constexpr (Detail::LinearPDESolver::linearSolverHasSetMatrix<LinearSolver, JacobianMatrix>())
249 virtual bool solveLinearSystem_(ResidualVector& deltaU)
251 if constexpr (Detail::LinearPDESolver::linearSolverHasSetMatrix<LinearSolver, JacobianMatrix>())
258 if (reuseMatrix_ && comm_.size() > 1)
259 DUNE_THROW(Dune::NotImplemented,
260 "Reuse matrix for parallel runs with a solver that doesn't support the setMatrix interface"
274 void initParams_(
const std::string& group =
"")
276 verbosity_ = comm_.rank() == 0 ? getParamFromGroup<int>(group,
"LinearPDESolver.Verbosity", 2) : 0;
277 enableDynamicOutput_ = getParamFromGroup<bool>(group,
"LinearPDESolver.EnableDynamicOutput",
true);
284 bool enableDynamicOutput_;
287 std::string paramGroup_;
Definition: variablesbackend.hh:187
An implementation of a linear PDE solver.
Definition: linear/pdesolver.hh:63
LinearPDESolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const std::string ¶mGroup)
The Constructor.
Definition: linear/pdesolver.hh:99
void reuseMatrix(bool reuse=true)
Set whether the matrix should be reused.
Definition: linear/pdesolver.hh:238
void solve(Variables &vars) override
Solve a linear PDE system.
Definition: linear/pdesolver.hh:191
LinearPDESolver(std::shared_ptr< Assembler > assembler, std::shared_ptr< LinearSolver > linearSolver, const Communication &comm=Dune::MPIHelper::getCommunication(), const std::string ¶mGroup="")
The Constructor.
Definition: linear/pdesolver.hh:81
int verbosity() const
Returns true if the solver ought to be chatty.
Definition: linear/pdesolver.hh:223
void setVerbosity(int val)
Specifies if the solver ought to be chatty.
Definition: linear/pdesolver.hh:217
Comm Communication
Definition: linear/pdesolver.hh:76
const std::string & paramGroup() const
Returns the parameter group.
Definition: linear/pdesolver.hh:229
Scalar suggestTimeStepSize(Scalar oldTimeStep) const
Suggest a new time-step size based on the old time-step size.
Definition: linear/pdesolver.hh:209
void report(std::ostream &sout=std::cout) const
output statistics / report
Definition: linear/pdesolver.hh:202
bool apply(Variables &vars) override
Solve a linear PDE system.
Definition: linear/pdesolver.hh:108
Base class for linear solvers.
Definition: solver.hh:27
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
A high-level interface for a PDESolver.
Definition: common/pdesolver.hh:61
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
Detail::PDESolver::AssemblerVariables< Assembler > Variables
export the type of variables that represent a numerical solution
Definition: common/pdesolver.hh:71
Defines a high-level interface for a PDESolver.
Manages the handling of time dependent problems.
Some exceptions thrown in DuMux
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
Definition: linear/pdesolver.hh:37
decltype(std::declval< Solver >().setMatrix(std::declval< std::shared_ptr< Matrix > >())) SetMatrixDetector
Definition: linear/pdesolver.hh:40
static constexpr bool linearSolverHasSetMatrix()
Definition: linear/pdesolver.hh:43
Definition: common/pdesolver.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Backends for operations on different solution vector types or more generic variable classes to be use...
Type traits to be used with vector types.