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 static bool warningIssued =
false;
177 using DM =
typename GridGeometry::VertexMapper;
180 PVHelper vectorHelper(
gridView(), gridGeometry_->vertexMapper());
182 vectorHelper.makeNonOverlappingConsistent(
residual);
185 else if (!warningIssued)
188 std::cout <<
"\nWarning: norm calculation adds entries corresponding to\n"
189 <<
"overlapping entities multiple times. Please use the norm\n"
190 <<
"function provided by a linear solver instead." << std::endl;
192 warningIssued =
true;
198 result2 =
gridView().comm().sum(result2);
201 return sqrt(result2);
210 std::shared_ptr<SolutionVector> r)
216 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
217 jacobian_->setBuildMode(JacobianMatrix::random);
218 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
219 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
231 jacobian_ = std::make_shared<JacobianMatrix>();
232 jacobian_->setBuildMode(JacobianMatrix::random);
233 residual_ = std::make_shared<SolutionVector>();
249 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
252 occupationPattern.exportIdx(*jacobian_);
257 { residual_->resize(
numDofs()); }
261 {
return gridGeometry_->numDofs(); }
265 {
return *problem_; }
269 {
return *gridGeometry_; }
277 {
return *gridVariables_; }
281 {
return *gridVariables_; }
285 {
return *jacobian_; }
289 {
return *residual_; }
293 {
return *prevSol_; }
300 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
313 {
return isStationaryProblem_; }
319 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
339 void resetResidual_()
343 residual_ = std::make_shared<SolutionVector>();
351 template <
class PartialReassembler = DefaultPartialReassembler>
352 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
356 jacobian_ = std::make_shared<JacobianMatrix>();
357 jacobian_->setBuildMode(JacobianMatrix::random);
361 if (partialReassembler)
362 partialReassembler->resetJacobian(*
this);
368 void checkAssemblerState_()
const
370 if (!isStationaryProblem_ && !prevSol_)
371 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
379 template<
typename AssembleElementFunc>
380 void assemble_(AssembleElementFunc&& assembleElement)
const
383 bool succeeded =
false;
389 for (
const auto& element : elements(
gridView()))
390 assembleElement(element);
396 catch (NumericalProblem &e)
398 std::cout <<
"rank " <<
gridView().comm().rank()
399 <<
" caught an exception while assembling:" << e.what()
406 succeeded =
gridView().comm().min(succeeded);
410 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
413 template<
class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
418 if (m.first < m.second)
421 res[m.first] += res[m.second];
422 const auto end = jac[m.second].end();
423 for (
auto it = jac[m.second].begin(); it != end; ++it)
424 jac[m.first][it.index()] += (*it);
427 res[m.second] = curSol[m.second] - curSol[m.first];
428 for (
auto it = jac[m.second].begin(); it != end; ++it)
429 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
434 template<
class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
435 enforcePeriodicConstraints_(
JacobianMatrix& jac, SolutionVector& res,
const SolutionVector& curSol,
const GG&
gridGeometry) {}
438 std::shared_ptr<const Problem> problem_;
441 std::shared_ptr<const GridGeometry> gridGeometry_;
444 std::shared_ptr<GridVariables> gridVariables_;
447 std::shared_ptr<const TimeLoop> timeLoop_;
450 const SolutionVector* prevSol_ =
nullptr;
453 bool isStationaryProblem_;
456 std::shared_ptr<JacobianMatrix> jacobian_;
457 std::shared_ptr<SolutionVector> residual_;
An enum class to define various differentiation methods available in order to compute the derivatives...
An assembler for Jacobian and residual contribution per element (box method)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
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
Definition: propertysystem.hh:150
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:136
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:272
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:264
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:229
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:73
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:312
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:318
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:260
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:276
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:284
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:299
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:256
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:242
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:306
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:268
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:332
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:280
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:209
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:288
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:292
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:324
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
Definition: parallelhelpers.hh:430
Declares all properties used in Dumux.
Manages the handling of time dependent problems.