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>,
87 , isStationaryProblem_(true)
89 static_assert(isImplicit,
"Explicit assembler for stationary problem doesn't make sense!");
100 std::shared_ptr<const TimeLoop> timeLoop,
105 , timeLoop_(timeLoop)
107 , isStationaryProblem_(!timeLoop)
114 template<
class PartialReassembler = DefaultPartialReassembler>
117 checkAssemblerState_();
118 resetJacobian_(partialReassembler);
121 assemble_([&](
const Element& element)
123 LocalAssembler localAssembler(*
this, element, curSol);
124 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
127 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
135 checkAssemblerState_();
138 assemble_([&](
const Element& element)
140 LocalAssembler localAssembler(*
this, element, curSol);
141 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
155 checkAssemblerState_();
157 assemble_([&](
const Element& element)
159 LocalAssembler localAssembler(*
this, element, curSol);
160 localAssembler.assembleResidual(r);
171 if (isBox &&
gridView().comm().size() > 1)
173 using VertexMapper =
typename GridGeometry::VertexMapper;
175 sumResidualHandle(gridGeometry_->vertexMapper(),
residual);
176 gridView().communicate(sumResidualHandle,
177 Dune::InteriorBorder_InteriorBorder_Interface,
178 Dune::ForwardCommunication);
184 result2 =
gridView().comm().sum(result2);
187 return sqrt(result2);
196 std::shared_ptr<SolutionVector> r)
202 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
203 jacobian_->setBuildMode(JacobianMatrix::random);
204 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
205 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
217 jacobian_ = std::make_shared<JacobianMatrix>();
218 jacobian_->setBuildMode(JacobianMatrix::random);
219 residual_ = std::make_shared<SolutionVector>();
235 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
238 occupationPattern.exportIdx(*jacobian_);
243 { residual_->resize(
numDofs()); }
247 {
return gridGeometry_->numDofs(); }
251 {
return *problem_; }
255 {
return *gridGeometry_; }
263 {
return *gridVariables_; }
267 {
return *gridVariables_; }
271 {
return *jacobian_; }
275 {
return *residual_; }
279 {
return *prevSol_; }
286 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
299 {
return isStationaryProblem_; }
305 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
325 void resetResidual_()
329 residual_ = std::make_shared<SolutionVector>();
337 template <
class PartialReassembler = DefaultPartialReassembler>
338 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
342 jacobian_ = std::make_shared<JacobianMatrix>();
343 jacobian_->setBuildMode(JacobianMatrix::random);
347 if (partialReassembler)
348 partialReassembler->resetJacobian(*
this);
354 void checkAssemblerState_()
const
356 if (!isStationaryProblem_ && !prevSol_)
357 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
365 template<
typename AssembleElementFunc>
366 void assemble_(AssembleElementFunc&& assembleElement)
const
369 bool succeeded =
false;
375 for (
const auto& element : elements(
gridView()))
376 assembleElement(element);
382 catch (NumericalProblem &e)
384 std::cout <<
"rank " <<
gridView().comm().rank()
385 <<
" caught an exception while assembling:" << e.what()
392 succeeded =
gridView().comm().min(succeeded);
396 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
399 template<
class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
404 if (m.first < m.second)
407 res[m.first] += res[m.second];
408 const auto end = jac[m.second].end();
409 for (
auto it = jac[m.second].begin(); it != end; ++it)
410 jac[m.first][it.index()] += (*it);
413 res[m.second] = curSol[m.second] - curSol[m.first];
414 for (
auto it = jac[m.second].begin(); it != end; ++it)
415 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
420 template<
class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
421 enforcePeriodicConstraints_(
JacobianMatrix& jac, SolutionVector& res,
const SolutionVector& curSol,
const GG&
gridGeometry) {}
424 std::shared_ptr<const Problem> problem_;
427 std::shared_ptr<const GridGeometry> gridGeometry_;
430 std::shared_ptr<GridVariables> gridVariables_;
433 std::shared_ptr<const TimeLoop> timeLoop_;
436 const SolutionVector* prevSol_ =
nullptr;
439 bool isStationaryProblem_;
442 std::shared_ptr<JacobianMatrix> jacobian_;
443 std::shared_ptr<SolutionVector> residual_;
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.
An assembler for Jacobian and residual contribution per element (box method)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
The available discretization methods in Dumux.
Contains a class to exchange entries of a vector.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:121
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:52
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:133
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:258
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:250
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:215
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:73
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:298
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:304
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:97
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:165
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:246
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:262
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:270
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:153
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:70
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:285
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:242
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:146
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:71
void setJacobianPattern()
Resizes the jacobian and sets the jacobian' sparsity pattern.
Definition: assembly/fvassembler.hh:228
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:292
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:254
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:318
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:266
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:195
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:115
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:80
SolutionVector & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:274
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:278
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:310
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:425
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78
Declares all properties used in Dumux.
Manages the handling of time dependent problems.