12#ifndef DUMUX_FV_ASSEMBLER_HH
13#define DUMUX_FV_ASSEMBLER_HH
20#include <dune/istl/matrixindexset.hh>
45template<
class DiscretizationMethod>
51 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
58 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
65 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
72 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
76template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
79>::template type<TypeTag, Impl, diffMethod, isImplicit>;
87{
return Dune::Std::is_detected<ProblemConstraintsDetector, P>::value; }
101template<
class TypeTag, DiffMethod diffMethod,
bool isImplicit = true,
class LocalRes
idual = GetPropType<TypeTag, Properties::LocalRes
idual>>
105 using GridView =
typename GridGeo::GridView;
106 using Element =
typename GridView::template Codim<0>::Entity;
107 using ElementSeed =
typename GridView::Grid::template Codim<0>::EntitySeed;
138 , isStationaryProblem_(true)
140 static_assert(isImplicit,
"Explicit assembler for stationary problem doesn't make sense!");
144 && getParam<bool>(
"Assembly.Multithreading",
true);
146 maybeComputeColors_();
157 std::shared_ptr<const TimeLoop> timeLoop,
162 , timeLoop_(timeLoop)
164 , isStationaryProblem_(!timeLoop)
169 && getParam<bool>(
"Assembly.Multithreading",
true);
171 maybeComputeColors_();
178 template<
class PartialReassembler = DefaultPartialReassembler>
181 checkAssemblerState_();
182 resetJacobian_(partialReassembler);
185 assemble_([&](
const Element& element)
187 LocalAssembler localAssembler(*
this, element, curSol);
188 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
191 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
193 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
198 (*residual_)[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
200 auto& row = (*jacobian_)[dofIdx];
201 for (
auto col = row.begin(); col != row.end(); ++col)
202 row[col.index()][eqIdx] = 0.0;
204 (*jacobian_)[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
206 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
214 checkAssemblerState_();
217 assemble_([&](
const Element& element)
219 LocalAssembler localAssembler(*
this, element, curSol);
220 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
223 enforcePeriodicConstraints_(*jacobian_, curSol, *gridGeometry_);
225 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
230 auto& row = (*jacobian_)[dofIdx];
231 for (
auto col = row.begin(); col != row.end(); ++col)
232 row[col.index()][eqIdx] = 0.0;
234 (*jacobian_)[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
236 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
249 checkAssemblerState_();
251 assemble_([&](
const Element& element)
253 LocalAssembler localAssembler(*
this, element, curSol);
254 localAssembler.assembleResidual(r);
257 enforcePeriodicConstraints_(*residual_, curSol, *gridGeometry_);
259 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
264 r[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
266 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
275 std::shared_ptr<ResidualType> r)
281 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
282 jacobian_->setBuildMode(JacobianMatrix::random);
283 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
284 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
287 setJacobianPattern_();
296 jacobian_ = std::make_shared<JacobianMatrix>();
297 jacobian_->setBuildMode(JacobianMatrix::random);
298 residual_ = std::make_shared<ResidualType>();
301 setJacobianPattern_();
310 setJacobianPattern_();
311 maybeComputeColors_();
316 {
return gridGeometry_->numDofs(); }
320 {
return *problem_; }
324 {
return *gridGeometry_; }
332 {
return *gridVariables_; }
336 {
return *gridVariables_; }
340 {
return *jacobian_; }
344 {
return *residual_; }
348 {
return *prevSol_; }
355 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
368 {
return isStationaryProblem_; }
374 {
return LocalResidual(problem_.get(), timeLoop_.get()); }
396 void setJacobianPattern_()
403 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
406 occupationPattern.exportIdx(*jacobian_);
410 void setResidualSize_()
411 { residual_->resize(
numDofs()); }
414 void maybeComputeColors_()
416 if (enableMultithreading_)
421 void resetResidual_()
425 residual_ = std::make_shared<ResidualType>();
433 template <
class PartialReassembler = DefaultPartialReassembler>
434 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
438 jacobian_ = std::make_shared<JacobianMatrix>();
439 jacobian_->setBuildMode(JacobianMatrix::random);
440 setJacobianPattern_();
443 if (partialReassembler)
444 partialReassembler->resetJacobian(*
this);
450 void checkAssemblerState_()
const
452 if (!isStationaryProblem_ && !prevSol_)
453 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
461 template<
typename AssembleElementFunc>
462 void assemble_(AssembleElementFunc&& assembleElement)
const
465 bool succeeded =
false;
470 if (enableMultithreading_)
472 assert(elementSets_.size() > 0);
478 for (
const auto& elements : elementSets_)
482 const auto element = gridView().grid().entity(elements[i]);
483 assembleElement(element);
488 for (
const auto& element : elements(
gridView()))
489 assembleElement(element);
495 catch (NumericalProblem &e)
497 std::cout <<
"rank " <<
gridView().comm().rank()
498 <<
" caught an exception while assembling:" << e.what()
505 succeeded =
gridView().comm().min(succeeded);
509 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
515 if constexpr (Detail::hasPeriodicDofMap<GG>())
519 if (m.first < m.second)
522 res[m.first] += res[m.second];
523 const auto end = jac[m.second].end();
524 for (
auto it = jac[m.second].begin(); it != end; ++it)
525 jac[m.first][it.index()] += (*it);
528 res[m.second] = curSol[m.second] - curSol[m.first];
531 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
533 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
534 matrixBlock[eIdx][eIdx] = diagValue;
537 for (
auto it = jac[m.second].begin(); it != end; ++it)
539 auto& matrixBlock = *it;
542 assert(matrixBlock.N() == matrixBlock.M());
543 if(it.index() == m.second)
544 setMatrixBlock(matrixBlock, 1.0);
546 if(it.index() == m.first)
547 setMatrixBlock(matrixBlock, -1.0);
558 if constexpr (Detail::hasPeriodicDofMap<GG>())
562 if (m.first < m.second)
565 const auto end = jac[m.second].end();
566 for (
auto it = jac[m.second].begin(); it != end; ++it)
567 jac[m.first][it.index()] += (*it);
570 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
572 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
573 matrixBlock[eIdx][eIdx] = diagValue;
576 for (
auto it = jac[m.second].begin(); it != end; ++it)
578 auto& matrixBlock = *it;
581 assert(matrixBlock.N() == matrixBlock.M());
582 if(it.index() == m.second)
583 setMatrixBlock(matrixBlock, 1.0);
585 if(it.index() == m.first)
586 setMatrixBlock(matrixBlock, -1.0);
597 if constexpr (Detail::hasPeriodicDofMap<GG>())
601 if (m.first < m.second)
604 res[m.first] += res[m.second];
607 res[m.second] = curSol[m.second] - curSol[m.first];
613 template<
class Problem,
class GG,
typename ApplyFunction>
614 void enforceProblemConstraints_(
const Problem&
problem,
const GG&
gridGeometry,
const ApplyFunction& applyDirichletConstraint)
const
616 if constexpr (Detail::hasGlobalConstraints<Problem>())
618 for (
const auto& constraintData :
problem.constraints())
620 const auto& constraintInfo = constraintData.constraintInfo();
621 const auto& values = constraintData.values();
622 const auto dofIdx = constraintData.dofIndex();
624 for (
int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
626 if (constraintInfo.isConstraintEquation(eqIdx))
628 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
629 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
630 applyDirichletConstraint(dofIdx, values, eqIdx, pvIdx);
638 std::shared_ptr<const Problem> problem_;
641 std::shared_ptr<const GridGeometry> gridGeometry_;
644 std::shared_ptr<GridVariables> gridVariables_;
647 std::shared_ptr<const TimeLoop> timeLoop_;
653 bool isStationaryProblem_;
656 std::shared_ptr<JacobianMatrix> jacobian_;
657 std::shared_ptr<ResidualType> residual_;
660 bool enableMultithreading_ =
false;
661 std::deque<std::vector<ElementSeed>> elementSets_;
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: assembly/cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:317
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:103
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:361
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:123
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:121
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:154
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:331
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:315
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:335
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:373
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: assembly/fvassembler.hh:118
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:179
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: assembly/fvassembler.hh:307
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:319
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:367
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:323
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:247
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition: assembly/fvassembler.hh:119
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:354
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:387
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:212
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:379
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:339
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< ResidualType > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: assembly/fvassembler.hh:274
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:117
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:124
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:294
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:347
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:131
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:327
ResidualType & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:343
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:240
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:116
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:288
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:84
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper to extract native Dune vector types from particular Dumux types.
An assembler for Jacobian and residual contribution per element (cell-centered methods)
An assembler for Jacobian and residual contribution per element (CVFE methods)
An assembler for Jacobian and residual contribution per element (face-centered staggered methods)
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:45
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
Definition: cvfelocalresidual.hh:30
constexpr bool hasGlobalConstraints()
Definition: assembly/fvassembler.hh:86
decltype(std::declval< P >().constraints()) ProblemConstraintsDetector
helper struct detecting if problem has a constraints() function
Definition: assembly/fvassembler.hh:83
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod, isImplicit > LocalAssemblerChooser_t
Definition: assembly/fvassembler.hh:79
constexpr Box box
Definition: method.hh:156
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:34
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:241
Parallel for loop (multithreading)
Provides a helper class for nonoverlapping decomposition.
Type traits to detect periodicity support.
Definition: assembly/fvassembler.hh:46
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition: dunevectors.hh:59
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:294
Type traits to be used with vector types.