24#ifndef DUMUX_FV_ASSEMBLER_HH
25#define DUMUX_FV_ASSEMBLER_HH
29#include <dune/istl/matrixindexset.hh>
50template<
class TypeTag, DiffMethod diffMethod,
bool isImplicit = true>
55 using Element =
typename GridView::template Codim<0>::Entity;
63 using LocalAssembler = std::conditional_t<isBox, BoxLocalAssembler<TypeTag, ThisType, diffMethod, isImplicit>,
88 , isStationaryProblem_(true)
90 static_assert(isImplicit,
"Explicit assembler for stationary problem doesn't make sense!");
97 [[deprecated(
"Please use constructor taking the previous solution instead. Will be removed after release 3.2!")]]
101 std::shared_ptr<const TimeLoop> timeLoop)
105 , timeLoop_(timeLoop)
106 , isStationaryProblem_(!timeLoop)
117 std::shared_ptr<const TimeLoop> timeLoop,
122 , timeLoop_(timeLoop)
124 , isStationaryProblem_(!timeLoop)
131 template<
class PartialReassembler = DefaultPartialReassembler>
134 checkAssemblerState_();
135 resetJacobian_(partialReassembler);
138 assemble_([&](
const Element& element)
140 LocalAssembler localAssembler(*
this, element, curSol);
141 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
144 enforcePeriodicConstraints_(*jacobian_, *residual_, *gridGeometry_);
152 checkAssemblerState_();
155 assemble_([&](
const Element& element)
157 LocalAssembler localAssembler(*
this, element, curSol);
158 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
172 checkAssemblerState_();
175 gridVariables_->update(curSol);
177 assemble_([&](
const Element& element)
179 LocalAssembler localAssembler(*
this, element, curSol);
180 localAssembler.assembleResidual(r);
191 if (isBox &&
gridView().comm().size() > 1)
193 using VertexMapper =
typename GridGeometry::VertexMapper;
195 sumResidualHandle(
residual, gridGeometry_->vertexMapper());
196 gridView().communicate(sumResidualHandle,
197 Dune::InteriorBorder_InteriorBorder_Interface,
198 Dune::ForwardCommunication);
204 result2 =
gridView().comm().sum(result2);
207 return sqrt(result2);
216 std::shared_ptr<SolutionVector> r)
222 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
223 jacobian_->setBuildMode(JacobianMatrix::random);
224 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
225 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
237 jacobian_ = std::make_shared<JacobianMatrix>();
238 jacobian_->setBuildMode(JacobianMatrix::random);
239 residual_ = std::make_shared<SolutionVector>();
255 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
258 occupationPattern.exportIdx(*jacobian_);
263 { residual_->resize(
numDofs()); }
267 {
return gridGeometry_->numDofs(); }
271 {
return *problem_; }
274 [[deprecated(
"Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
280 {
return *gridGeometry_; }
288 {
return *gridVariables_; }
292 {
return *gridVariables_; }
296 {
return *jacobian_; }
300 {
return *residual_; }
304 {
return *prevSol_; }
311 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
324 {
return isStationaryProblem_; }
330 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
350 void resetResidual_()
354 residual_ = std::make_shared<SolutionVector>();
362 template <
class PartialReassembler = DefaultPartialReassembler>
363 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
367 jacobian_ = std::make_shared<JacobianMatrix>();
368 jacobian_->setBuildMode(JacobianMatrix::random);
372 if (partialReassembler)
373 partialReassembler->resetJacobian(*
this);
379 void checkAssemblerState_()
const
381 if (!isStationaryProblem_ && !prevSol_)
382 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
390 template<
typename AssembleElementFunc>
391 void assemble_(AssembleElementFunc&& assembleElement)
const
394 bool succeeded =
false;
400 for (
const auto& element : elements(
gridView()))
401 assembleElement(element);
407 catch (NumericalProblem &e)
409 std::cout <<
"rank " <<
gridView().comm().rank()
410 <<
" caught an exception while assembling:" << e.what()
417 succeeded =
gridView().comm().min(succeeded);
421 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
424 template<
class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
429 if (m.first < m.second)
432 res[m.first] += res[m.second];
433 const auto end = jac[m.second].end();
434 for (
auto it = jac[m.second].begin(); it != end; ++it)
435 jac[m.first][it.index()] += (*it);
438 for (
auto it = jac[m.second].begin(); it != end; ++it)
439 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
444 template<
class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
448 std::shared_ptr<const Problem> problem_;
451 std::shared_ptr<const GridGeometry> gridGeometry_;
454 std::shared_ptr<GridVariables> gridVariables_;
457 std::shared_ptr<const TimeLoop> timeLoop_;
460 const SolutionVector* prevSol_ =
nullptr;
463 bool isStationaryProblem_;
466 std::shared_ptr<JacobianMatrix> jacobian_;
467 std::shared_ptr<SolutionVector> residual_;
An assembler for Jacobian and residual contribution per element (box method)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper function to generate Jacobian pattern for different discretization methods.
Manages the handling of time dependent problems.
The available discretization methods in Dumux.
Provides data handles for parallel communication which operate on vertices.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:123
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:52
GridGeometry FVGridGeometry
Definition: assembly/fvassembler.hh:70
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const TimeLoop > timeLoop)
The constructor for instationary problems.
Definition: assembly/fvassembler.hh:98
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:67
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:150
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:283
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:270
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:235
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:74
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:323
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:68
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:329
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: assembly/fvassembler.hh:114
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: assembly/fvassembler.hh:69
Scalar residualNorm(const SolutionVector &curSol) const
compute the residual and return it's vector norm
Definition: assembly/fvassembler.hh:185
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:266
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:287
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:295
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:170
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:71
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:310
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:262
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:163
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:72
void setJacobianPattern()
Resizes the jacobian and sets the jacobian' sparsity pattern.
Definition: assembly/fvassembler.hh:248
void setPreviousSolution(const SolutionVector &u)
Sets the solution from which to start the time integration. Has to be called prior to assembly for ti...
Definition: assembly/fvassembler.hh:317
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:279
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:343
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:291
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< SolutionVector > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: assembly/fvassembler.hh:215
void assembleJacobianAndResidual(const SolutionVector &curSol, const PartialReassembler *partialReassembler=nullptr)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: assembly/fvassembler.hh:132
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables)
The constructor for stationary problems.
Definition: assembly/fvassembler.hh:81
SolutionVector & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:299
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:303
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:335
const GridGeometry & fvGridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:275
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:424
Manages the handling of time dependent problems.
Definition: timeloop.hh:65
Data handle for parallel communication which sums up all values are attached to vertices.
Definition: vertexhandles.hh:42
Declares all properties used in Dumux.