12#ifndef DUMUX_ASSEMBLER_HH
13#define DUMUX_ASSEMBLER_HH
22#include <dune/common/std/type_traits.hh>
23#include <dune/grid/common/rangegenerators.hh>
42#include "cvfelocalassembler_.hh"
46template<
class DiscretizationMethod>
52 template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
56template<
class TypeTag,
class Impl, DiffMethod diffMethod,
bool isImplicit>
59>::template type<TypeTag, Impl, diffMethod, isImplicit>;
67{
return Dune::Std::is_detected<ProblemConstraintsDetector, P>::value; }
81template<
class TypeTag, DiffMethod diffMethod,
bool isImplicit = true,
class LocalRes
idual = GetPropType<TypeTag, Properties::LocalRes
idual>>
85 using GridView =
typename GridGeo::GridView;
86 using Element =
typename GridView::template Codim<0>::Entity;
87 using ElementSeed =
typename GridView::Grid::template Codim<0>::EntitySeed;
116 , isStationaryProblem_(true)
118 static_assert(isImplicit,
"Explicit assembler for stationary problem doesn't make sense!");
122 && getParam<bool>(
"Assembly.Multithreading",
true);
124 maybeComputeColors_();
135 std::shared_ptr<const TimeLoop> timeLoop,
140 , timeLoop_(timeLoop)
142 , isStationaryProblem_(!timeLoop)
147 && getParam<bool>(
"Assembly.Multithreading",
true);
149 maybeComputeColors_();
156 template<
class PartialReassembler = DefaultPartialReassembler>
159 checkAssemblerState_();
160 resetJacobian_(partialReassembler);
163 assemble_([&](
const Element& element)
165 LocalAssembler localAssembler(*
this, element, curSol);
166 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
169 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
171 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
176 (*residual_)[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
178 auto& row = (*jacobian_)[dofIdx];
179 for (
auto col = row.begin(); col != row.end(); ++col)
180 row[col.index()][eqIdx] = 0.0;
182 (*jacobian_)[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
184 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
192 checkAssemblerState_();
195 assemble_([&](
const Element& element)
197 LocalAssembler localAssembler(*
this, element, curSol);
198 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
201 enforcePeriodicConstraints_(*jacobian_, curSol, *gridGeometry_);
203 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
208 auto& row = (*jacobian_)[dofIdx];
209 for (
auto col = row.begin(); col != row.end(); ++col)
210 row[col.index()][eqIdx] = 0.0;
212 (*jacobian_)[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
214 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
227 checkAssemblerState_();
229 assemble_([&](
const Element& element)
231 LocalAssembler localAssembler(*
this, element, curSol);
232 localAssembler.assembleResidual(r);
235 enforcePeriodicConstraints_(r, curSol, *gridGeometry_);
237 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
242 r[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
244 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
253 std::shared_ptr<ResidualType> r)
259 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
260 jacobian_->setBuildMode(JacobianMatrix::random);
261 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
262 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
265 setJacobianPattern_();
274 jacobian_ = std::make_shared<JacobianMatrix>();
275 jacobian_->setBuildMode(JacobianMatrix::random);
276 residual_ = std::make_shared<ResidualType>();
279 setJacobianPattern_();
288 setJacobianPattern_();
289 maybeComputeColors_();
294 {
return gridGeometry_->numDofs(); }
298 {
return *problem_; }
302 {
return *gridGeometry_; }
310 {
return *gridVariables_; }
314 {
return *gridVariables_; }
318 {
return *jacobian_; }
322 {
return *residual_; }
326 {
return *prevSol_; }
333 { timeLoop_ = timeLoop; isStationaryProblem_ = !
static_cast<bool>(timeLoop); }
346 {
return isStationaryProblem_; }
374 void setJacobianPattern_()
381 const auto occupationPattern = getJacobianPattern<isImplicit>(
gridGeometry());
384 occupationPattern.exportIdx(*jacobian_);
388 void setResidualSize_()
389 { residual_->resize(
numDofs()); }
392 void maybeComputeColors_()
394 if (enableMultithreading_)
399 void resetResidual_()
403 residual_ = std::make_shared<ResidualType>();
411 template <
class PartialReassembler = DefaultPartialReassembler>
412 void resetJacobian_(
const PartialReassembler *partialReassembler =
nullptr)
416 jacobian_ = std::make_shared<JacobianMatrix>();
417 jacobian_->setBuildMode(JacobianMatrix::random);
418 setJacobianPattern_();
421 if (partialReassembler)
422 partialReassembler->resetJacobian(*
this);
428 void checkAssemblerState_()
const
430 if (!isStationaryProblem_ && !prevSol_)
431 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
439 template<
typename AssembleElementFunc>
440 void assemble_(AssembleElementFunc&& assembleElement)
const
443 bool succeeded =
false;
448 if (enableMultithreading_)
450 assert(elementSets_.size() > 0);
456 for (
const auto& elements : elementSets_)
460 const auto element = gridView().grid().entity(elements[i]);
461 assembleElement(element);
466 for (
const auto& element : elements(
gridView()))
467 assembleElement(element);
473 catch (
const NumericalProblem& e)
475 std::cout <<
"rank " <<
gridView().comm().rank()
476 <<
" caught an exception while assembling:" << e.what()
483 succeeded =
gridView().comm().min(succeeded);
487 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
493 if constexpr (Dumux::Detail::hasPeriodicDofMap<GG>())
497 if (m.first < m.second)
500 res[m.first] += res[m.second];
501 const auto end = jac[m.second].end();
502 for (
auto it = jac[m.second].begin(); it != end; ++it)
503 jac[m.first][it.index()] += (*it);
506 res[m.second] = curSol[m.second] - curSol[m.first];
509 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
511 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
512 matrixBlock[eIdx][eIdx] = diagValue;
515 for (
auto it = jac[m.second].begin(); it != end; ++it)
517 auto& matrixBlock = *it;
520 assert(matrixBlock.N() == matrixBlock.M());
521 if(it.index() == m.second)
522 setMatrixBlock(matrixBlock, 1.0);
524 if(it.index() == m.first)
525 setMatrixBlock(matrixBlock, -1.0);
536 if constexpr (Dumux::Detail::hasPeriodicDofMap<GG>())
540 if (m.first < m.second)
543 const auto end = jac[m.second].end();
544 for (
auto it = jac[m.second].begin(); it != end; ++it)
545 jac[m.first][it.index()] += (*it);
548 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
550 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
551 matrixBlock[eIdx][eIdx] = diagValue;
554 for (
auto it = jac[m.second].begin(); it != end; ++it)
556 auto& matrixBlock = *it;
559 assert(matrixBlock.N() == matrixBlock.M());
560 if(it.index() == m.second)
561 setMatrixBlock(matrixBlock, 1.0);
563 if(it.index() == m.first)
564 setMatrixBlock(matrixBlock, -1.0);
575 if constexpr (Dumux::Detail::hasPeriodicDofMap<GG>())
579 if (m.first < m.second)
582 res[m.first] += res[m.second];
585 res[m.second] = curSol[m.second] - curSol[m.first];
591 template<
class Problem,
class GG,
typename ApplyFunction>
592 void enforceProblemConstraints_(
const Problem&
problem,
const GG&,
const ApplyFunction& applyDirichletConstraint)
const
594 if constexpr (Detail::hasGlobalConstraints<Problem>())
596 for (
const auto& constraintData :
problem.constraints())
598 const auto& constraintInfo = constraintData.constraintInfo();
599 const auto& values = constraintData.values();
600 const auto dofIdx = constraintData.dofIndex();
602 for (
int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
604 if (constraintInfo.isConstraintEquation(eqIdx))
606 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
607 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
608 applyDirichletConstraint(dofIdx, values, eqIdx, pvIdx);
616 std::shared_ptr<const Problem> problem_;
619 std::shared_ptr<const GridGeometry> gridGeometry_;
622 std::shared_ptr<GridVariables> gridVariables_;
625 std::shared_ptr<const TimeLoop> timeLoop_;
631 bool isStationaryProblem_;
634 std::shared_ptr<JacobianMatrix> jacobian_;
635 std::shared_ptr<ResidualType> residual_;
638 bool enableMultithreading_ =
false;
639 std::deque<std::vector<ElementSeed>> elementSets_;
A linear system assembler (residual and Jacobian) for general discretization schemes.
Definition: assembly/assembler.hh:83
Assembler(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/assembler.hh:132
typename Dumux::Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition: assembly/assembler.hh:97
void updateGridVariables(const SolutionVector &curSol)
Update the grid variables.
Definition: assembly/assembler.hh:357
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/assembler.hh:99
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/assembler.hh:218
Assembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables)
The constructor for stationary problems.
Definition: assembly/assembler.hh:109
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/assembler.hh:351
ResidualType & residual()
The residual vector (rhs)
Definition: assembly/assembler.hh:321
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/assembler.hh:157
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/assembler.hh:301
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/assembler.hh:94
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: assembly/assembler.hh:285
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/assembler.hh:325
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/assembler.hh:293
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/assembler.hh:225
GridGeo GridGeometry
Definition: assembly/assembler.hh:101
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: assembly/assembler.hh:96
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/assembler.hh:332
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/assembler.hh:252
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/assembler.hh:339
const GridView & gridView() const
The gridview.
Definition: assembly/assembler.hh:305
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/assembler.hh:102
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/assembler.hh:190
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/assembler.hh:317
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/assembler.hh:345
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/assembler.hh:313
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/assembler.hh:272
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/assembler.hh:95
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/assembler.hh:309
void resetTimeStep(const SolutionVector &curSol)
Reset the gridVariables.
Definition: assembly/assembler.hh:365
const Problem & problem() const
The problem.
Definition: assembly/assembler.hh:297
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:328
The element-wise residual for grid-based discretization schemes.
Definition: assembly/localresidual.hh:37
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.
Some exceptions thrown in DuMux
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: assembly/assembler.hh:44
constexpr bool hasGlobalConstraints()
Definition: assembly/assembler.hh:66
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod, isImplicit > LocalAssemblerChooser_t
Definition: assembly/assembler.hh:59
decltype(std::declval< P >().constraints()) ProblemConstraintsDetector
helper struct detecting if problem has a constraints() function
Definition: assembly/assembler.hh:63
Definition: assembly/assembler.hh:44
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)
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Type traits to detect periodicity support.
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition: dunevectors.hh:59
Definition: assembly/assembler.hh:47
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:294