14#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
20#include <dune/common/hybridutilities.hh>
21#include <dune/istl/matrixindexset.hh>
45namespace Grid::Capabilities {
49template<
class T, std::size_t... I>
71struct CouplingManagerSupportsMultithreadedAssembly :
public std::false_type
83template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
86 template<std::
size_t id>
87 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
92 using Scalar =
typename MDTraits::Scalar;
95 template<std::
size_t id>
98 template<std::
size_t id>
99 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
101 template<std::
size_t id>
102 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
104 template<std::
size_t id>
105 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
117 {
return useImplicitAssembly; }
121 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
122 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
123 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
128 template<std::
size_t id>
131 template<
class DiscretizationMethod, std::
size_t id>
132 struct SubDomainAssemblerType;
134 template<std::
size_t id>
135 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
140 template<std::
size_t id>
141 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
143 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
146 template<std::
size_t id,
class DM>
147 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
149 using type = SubDomainCVFELocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
152 template<std::
size_t id>
153 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
155 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
158 template<std::
size_t id>
159 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
161 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
164 template<std::
size_t id>
165 using SubDomainAssembler =
typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod,
id>::type;
180 , problemTuple_(std::move(
problem))
184 , isStationaryProblem_(true)
185 , warningIssued_(false)
187 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
188 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
193 && getParam<bool>(
"Assembly.Multithreading",
true);
195 maybeComputeColors_();
207 std::shared_ptr<const TimeLoop> timeLoop,
210 , problemTuple_(std::move(
problem))
213 , timeLoop_(timeLoop)
215 , isStationaryProblem_(false)
216 , warningIssued_(false)
218 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
223 && getParam<bool>(
"Assembly.Multithreading",
true);
225 maybeComputeColors_();
234 checkAssemblerState_();
238 using namespace Dune::Hybrid;
239 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
241 auto& jacRow = (*jacobian_)[domainId];
242 auto& subRes = (*residual_)[domainId];
243 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
245 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
246 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, curSol[domainId]);
262 checkAssemblerState_();
267 using namespace Dune::Hybrid;
268 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
270 auto& subRes = r[domainId];
271 this->assembleResidual_(domainId, subRes, curSol);
281 std::shared_ptr<ResidualType> r)
287 setJacobianPattern_(*jacobian_);
288 setResidualSize_(*residual_);
297 jacobian_ = std::make_shared<JacobianMatrix>();
298 residual_ = std::make_shared<ResidualType>();
301 setJacobianPattern_(*jacobian_);
302 setResidualSize_(*residual_);
310 using namespace Dune::Hybrid;
311 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
313 forEach(jac[i], [&](
auto& jacBlock)
315 using BlockType = std::decay_t<
decltype(jacBlock)>;
316 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
317 jacBlock.setBuildMode(BlockType::BuildMode::random);
318 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
319 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
329 setJacobianPattern_(*jacobian_);
330 setResidualSize_(*residual_);
331 maybeComputeColors_();
339 using namespace Dune::Hybrid;
340 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
341 { this->
gridVariables(domainId).update(curSol[domainId]); });
349 using namespace Dune::Hybrid;
350 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
351 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
355 template<std::
size_t i>
356 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
357 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
360 template<std::
size_t i>
361 const auto&
problem(Dune::index_constant<i> domainId)
const
362 {
return *std::get<domainId>(problemTuple_); }
365 template<std::
size_t i>
367 {
return *std::get<domainId>(gridGeometryTuple_); }
370 template<std::
size_t i>
371 const auto&
gridView(Dune::index_constant<i> domainId)
const
375 template<std::
size_t i>
377 {
return *std::get<domainId>(gridVariablesTuple_); }
380 template<std::
size_t i>
382 {
return *std::get<domainId>(gridVariablesTuple_); }
390 {
return *jacobian_; }
394 {
return *residual_; }
398 {
return *prevSol_; }
405 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
418 {
return isStationaryProblem_; }
423 template<std::
size_t i>
425 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
437 using namespace Dune::Hybrid;
438 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
440 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
442 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
443 pattern.exportIdx(jac[domainI][domainJ]);
453 using namespace Dune::Hybrid;
454 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
455 { res[domainId].resize(this->
numDofs(domainId)); });
459 void resetResidual_()
463 residual_ = std::make_shared<ResidualType>();
464 setResidualSize_(*residual_);
471 void resetJacobian_()
475 jacobian_ = std::make_shared<JacobianMatrix>();
477 setJacobianPattern_(*jacobian_);
484 void maybeComputeColors_()
486 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
487 if (enableMultithreading_)
492 void checkAssemblerState_()
const
494 if (!isStationaryProblem_ && !prevSol_)
495 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
497 if (isStationaryProblem_ && prevSol_)
498 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
499 <<
" Did you forget to set the timeLoop to make this problem instationary?");
502 template<std::
size_t i,
class JacRow,
class SubRes>
503 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
506 assemble_(domainId, [&](
const auto& element)
508 MultiDomainAssemblerSubDomainView view{*
this, domainId};
509 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
510 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
514 template<std::
size_t i,
class SubRes>
515 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
518 assemble_(domainId, [&](
const auto& element)
520 MultiDomainAssemblerSubDomainView view{*
this, domainId};
521 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
522 subDomainAssembler.assembleResidual(subRes);
531 template<std::
size_t i,
class AssembleElementFunc>
532 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
535 bool succeeded =
false;
540 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
542 if (enableMultithreading_)
545 domainId, std::forward<AssembleElementFunc>(assembleElement)
554 for (
const auto& element : elements(
gridView(domainId)))
555 assembleElement(element);
561 catch (NumericalProblem &e)
563 std::cout <<
"rank " <<
gridView(domainId).comm().rank()
564 <<
" caught an exception while assembling:" << e.what()
570 if (
gridView(domainId).comm().size() > 1)
571 succeeded =
gridView(domainId).comm().min(succeeded);
575 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
579 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
580 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
581 Dune::index_constant<j> domainJ)
const
584 auto pattern = getJacobianPattern<isImplicit()>(gg);
590 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
591 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
592 Dune::index_constant<j> domainJ)
const
600 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
601 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
607 if (m.first < m.second)
609 auto& jac = jacRow[domainI];
612 res[m.first] += res[m.second];
615 res[m.second] = curSol[m.second] - curSol[m.first];
617 const auto end = jac[m.second].end();
618 for (
auto it = jac[m.second].begin(); it != end; ++it)
619 jac[m.first][it.index()] += (*it);
622 for (
auto it = jac[m.second].begin(); it != end; ++it)
623 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
625 using namespace Dune::Hybrid;
628 auto& jacCoupling = jacRow[couplingDomainId];
630 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
631 jacCoupling[m.first][it.index()] += (*it);
633 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
642 ProblemTuple problemTuple_;
645 GridGeometryTuple gridGeometryTuple_;
648 GridVariablesTuple gridVariablesTuple_;
651 std::shared_ptr<const TimeLoop> timeLoop_;
657 bool isStationaryProblem_;
660 std::shared_ptr<JacobianMatrix> jacobian_;
661 std::shared_ptr<ResidualType> residual_;
664 mutable bool warningIssued_;
667 bool enableMultithreading_ =
false;
Subdomain-specific views on multidomain assemblers.
Subdomain-specific view on a multidomain assembler. Allows retrieval of sub-domain specific objects w...
Definition: assemblerview.hh:31
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:85
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:356
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:337
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:327
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:108
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:105
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:389
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:116
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:258
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:107
CMType CouplingManager
Definition: multidomain/fvassembler.hh:111
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:251
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:347
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:381
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: multidomain/fvassembler.hh:411
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:99
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:366
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:295
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:371
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:102
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:385
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:404
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:397
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:376
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:429
MDTraits Traits
Definition: multidomain/fvassembler.hh:90
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: multidomain/fvassembler.hh:203
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:308
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:96
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:424
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:361
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/fvassembler.hh:109
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:175
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:232
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:417
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:92
ResidualType & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:393
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: multidomain/fvassembler.hh:280
The cell-centered scheme multidomain local assembler.
Definition: multidomain/subdomaincclocalassembler.hh:270
The default time loop for instationary simulations.
Definition: common/timeloop.hh:139
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
Helper function to generate Jacobian pattern for multi domain models.
An enum class to define various differentiation methods available in order to compute the derivatives...
Some exceptions thrown in DuMux
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:45
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.
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
An assembler for Jacobian and residual contribution per element (CVFE methods) for multidomain proble...
constexpr Box box
Definition: method.hh:147
constexpr FCStaggered fcstaggered
Definition: method.hh:151
bool allGridsSupportsMultithreadingImpl(const T &gridGeometries, std::index_sequence< I... >)
Definition: multistagemultidomainfvassembler.hh:52
bool allGridsSupportsMultithreading(const std::tuple< GG... > &gridGeometries)
Definition: multistagemultidomainfvassembler.hh:60
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:74
typename Detail::ConcatSeq< decltype(std::make_index_sequence< e >{}), e+1, decltype(std::make_index_sequence<(n > e) ?(n - e - 1) :0 >{})>::type makeIncompleteIntegerSequence
Definition: utility.hh:59
Provides a helper class for nonoverlapping decomposition.
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for...
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
Utilities for template meta programming.