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!");
96 [[deprecated(
"Please use constructor taking the previous solution instead. Will be removed after release 3.2!")]]
100 std::shared_ptr<const TimeLoop> timeLoop)
104 , timeLoop_(timeLoop)
105 , isStationaryProblem_(!timeLoop)
116 std::shared_ptr<const TimeLoop> timeLoop,
121 , timeLoop_(timeLoop)
123 , isStationaryProblem_(!timeLoop)
130 template<
class PartialReassembler = DefaultPartialReassembler>
133 checkAssemblerState_();
134 resetJacobian_(partialReassembler);
137 assemble_([&](
const Element& element)
139 LocalAssembler localAssembler(*
this, element, curSol);
140 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
143 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
151 checkAssemblerState_();
154 assemble_([&](
const Element& element)
156 LocalAssembler localAssembler(*
this, element, curSol);
157 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
171 checkAssemblerState_();
174 gridVariables_->update(curSol);
176 assemble_([&](
const Element& element)
178 LocalAssembler localAssembler(*
this, element, curSol);
179 localAssembler.assembleResidual(r);
190 if (isBox &&
gridView().comm().size() > 1)
192 using VertexMapper =
typename GridGeometry::VertexMapper;
194 sumResidualHandle(gridGeometry_->vertexMapper(),
residual);
195 gridView().communicate(sumResidualHandle,
196 Dune::InteriorBorder_InteriorBorder_Interface,
197 Dune::ForwardCommunication);
203 result2 =
gridView().comm().sum(result2);
206 return sqrt(result2);
215 std::shared_ptr<SolutionVector> r)
221 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
222 jacobian_->setBuildMode(JacobianMatrix::random);
223 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
224 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
236 jacobian_ = std::make_shared<JacobianMatrix>();
237 jacobian_->setBuildMode(JacobianMatrix::random);
238 residual_ = std::make_shared<SolutionVector>();
257 occupationPattern.exportIdx(*jacobian_);
262 { residual_->resize(
numDofs()); }
266 {
return gridGeometry_->numDofs(); }
270 {
return *problem_; }
274 {
return *gridGeometry_; }
282 {
return *gridVariables_; }
286 {
return *gridVariables_; }
290 {
return *jacobian_; }
294 {
return *residual_; }
298 {
return *prevSol_; }
305 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
318 {
return isStationaryProblem_; }
324 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
344 void resetResidual_()
348 residual_ = std::make_shared<SolutionVector>();
356 template <
class PartialReassembler = DefaultPartialReassembler>
357 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
361 jacobian_ = std::make_shared<JacobianMatrix>();
362 jacobian_->setBuildMode(JacobianMatrix::random);
366 if (partialReassembler)
367 partialReassembler->resetJacobian(*
this);
373 void checkAssemblerState_()
const
375 if (!isStationaryProblem_ && !prevSol_)
376 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
384 template<
typename AssembleElementFunc>
385 void assemble_(AssembleElementFunc&& assembleElement)
const
388 bool succeeded =
false;
394 for (
const auto& element : elements(
gridView()))
395 assembleElement(element);
401 catch (NumericalProblem &e)
403 std::cout <<
"rank " <<
gridView().comm().rank()
404 <<
" caught an exception while assembling:" << e.what()
411 succeeded =
gridView().comm().min(succeeded);
415 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
418 template<
class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
423 if (m.first < m.second)
426 res[m.first] += res[m.second];
427 const auto end = jac[m.second].end();
428 for (
auto it = jac[m.second].begin(); it != end; ++it)
429 jac[m.first][it.index()] += (*it);
432 res[m.second] = curSol[m.second] - curSol[m.first];
433 for (
auto it = jac[m.second].begin(); it != end; ++it)
434 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
439 template<
class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
440 enforcePeriodicConstraints_(
JacobianMatrix& jac, SolutionVector& res,
const SolutionVector& curSol,
const GG&
gridGeometry) {}
443 std::shared_ptr<const Problem> problem_;
446 std::shared_ptr<const GridGeometry> gridGeometry_;
449 std::shared_ptr<GridVariables> gridVariables_;
452 std::shared_ptr<const TimeLoop> timeLoop_;
455 const SolutionVector* prevSol_ =
nullptr;
458 bool isStationaryProblem_;
461 std::shared_ptr<JacobianMatrix> jacobian_;
462 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).
Manages the handling of time dependent problems.
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
@ box
Definition method.hh:38
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry &gridGeometry)
Helper function to generate Jacobian pattern for the box method.
Definition jacobianpattern.hh:39
VectorCommDataHandle< Mapper, Vector, codim, Detail::Sum > VectorCommDataHandleSum
Definition vectorcommdatahandle.hh:136
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
const Scalar PengRobinsonMixture< Scalar, StaticParameters >::u
Definition pengrobinsonmixture.hh:167
An assembler for Jacobian and residual contribution per element (cell-centered methods).
Definition cclocalassembler.hh:123
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:97
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:149
const GridView & gridView() const
The gridview.
Definition assembly/fvassembler.hh:277
const Problem & problem() const
Definition assembly/fvassembler.hh:269
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition assembly/fvassembler.hh:234
SolutionVector ResidualType
Definition assembly/fvassembler.hh:73
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition assembly/fvassembler.hh:317
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:323
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:113
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:184
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition assembly/fvassembler.hh:265
GridVariables & gridVariables()
Definition assembly/fvassembler.hh:281
JacobianMatrix & jacobian()
The jacobian matrix.
Definition assembly/fvassembler.hh:289
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition assembly/fvassembler.hh:169
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:304
void setResidualSize()
Resizes the residual.
Definition assembly/fvassembler.hh:261
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition assembly/fvassembler.hh:162
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:247
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:311
const GridGeometry & gridGeometry() const
Definition assembly/fvassembler.hh:273
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition assembly/fvassembler.hh:337
const GridVariables & gridVariables() const
The global grid variables.
Definition assembly/fvassembler.hh:285
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:214
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:131
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:293
const SolutionVector & prevSol() const
Definition assembly/fvassembler.hh:297
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition assembly/fvassembler.hh:329
detects which entries in the Jacobian have to be recomputed
Definition partialreassembler.hh:424
Manages the handling of time dependent problems.
Definition timeloop.hh:72
Declares all properties used in Dumux.