14#ifndef DUMUX_MULTIDOMAIN_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_ASSEMBLER_HH
22#include <dune/common/hybridutilities.hh>
23#include <dune/istl/matrixindexset.hh>
39#include "subdomaincvfelocalassembler_.hh"
44#if HAVE_DUMUX_OLD_STAGGERED
45#include <dumux/multidomain/subdomainstaggeredlocalassembler.hh>
61template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
64 template<std::
size_t id>
65 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
70 using Scalar =
typename MDTraits::Scalar;
72 template<std::
size_t id>
73 using LocalResidual =
typename MDTraits::template SubDomain<id>::LocalResidual;
75 template<std::
size_t id>
76 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
78 template<std::
size_t id>
79 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
81 template<std::
size_t id>
82 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
94 {
return useImplicitAssembly; }
98 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
99 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
100 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
105 template<std::
size_t id>
108 template<
class DiscretizationMethod, std::
size_t id>
109 struct SubDomainAssemblerType;
111 template<std::
size_t id>
112 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
117 template<std::
size_t id>
118 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
123 template<std::
size_t id,
class DM>
124 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
129 template<std::
size_t id>
130 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
132 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
135 template<std::
size_t id>
136 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
138 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
141 template<std::
size_t id>
142 using SubDomainAssembler =
typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod,
id>::type;
156 , problemTuple_(std::move(
problem))
160 , isStationaryProblem_(true)
161 , warningIssued_(false)
163 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
164 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
169 && getParam<bool>(
"Assembly.Multithreading",
true);
171 maybeComputeColors_();
183 std::shared_ptr<const TimeLoop> timeLoop,
186 , problemTuple_(std::move(
problem))
189 , timeLoop_(timeLoop)
191 , isStationaryProblem_(false)
192 , warningIssued_(false)
194 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
199 && getParam<bool>(
"Assembly.Multithreading",
true);
201 maybeComputeColors_();
210 checkAssemblerState_();
214 using namespace Dune::Hybrid;
215 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
217 auto& jacRow = (*jacobian_)[domainId];
218 auto& subRes = (*residual_)[domainId];
219 const auto& subSol = curSol[domainId];
220 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
222 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
223 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, subSol);
224 enforceProblemConstraints_(domainId, jacRow, subRes, *
gridGeometry, subSol);
240 checkAssemblerState_();
245 using namespace Dune::Hybrid;
246 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
248 auto& subRes = r[domainId];
249 const auto& subSol = curSol[domainId];
250 this->assembleResidual_(domainId, subRes, curSol);
252 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
253 enforcePeriodicConstraints_(domainId, subRes, *
gridGeometry, subSol);
254 enforceProblemConstraints_(domainId, subRes, *
gridGeometry, subSol);
264 std::shared_ptr<ResidualType> r)
270 setJacobianPattern_(*jacobian_);
271 setResidualSize_(*residual_);
280 jacobian_ = std::make_shared<JacobianMatrix>();
281 residual_ = std::make_shared<ResidualType>();
284 setJacobianPattern_(*jacobian_);
285 setResidualSize_(*residual_);
293 using namespace Dune::Hybrid;
294 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
296 forEach(jac[i], [&](
auto& jacBlock)
298 using BlockType = std::decay_t<
decltype(jacBlock)>;
299 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
300 jacBlock.setBuildMode(BlockType::BuildMode::random);
301 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
302 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
312 setJacobianPattern_(*jacobian_);
313 setResidualSize_(*residual_);
314 maybeComputeColors_();
322 using namespace Dune::Hybrid;
323 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
324 { this->
gridVariables(domainId).update(curSol[domainId]); });
332 using namespace Dune::Hybrid;
333 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
334 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
338 template<std::
size_t i>
339 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
340 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
343 template<std::
size_t i>
344 const auto&
problem(Dune::index_constant<i> domainId)
const
345 {
return *std::get<domainId>(problemTuple_); }
348 template<std::
size_t i>
350 {
return *std::get<domainId>(gridGeometryTuple_); }
353 template<std::
size_t i>
354 const auto&
gridView(Dune::index_constant<i> domainId)
const
358 template<std::
size_t i>
360 {
return *std::get<domainId>(gridVariablesTuple_); }
363 template<std::
size_t i>
365 {
return *std::get<domainId>(gridVariablesTuple_); }
373 {
return *jacobian_; }
377 {
return *residual_; }
381 {
return *prevSol_; }
388 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
401 {
return isStationaryProblem_; }
406 template<std::
size_t i>
408 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
420 using namespace Dune::Hybrid;
421 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
423 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
425 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
426 pattern.exportIdx(jac[domainI][domainJ]);
436 using namespace Dune::Hybrid;
437 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
438 { res[domainId].resize(this->
numDofs(domainId)); });
442 void resetResidual_()
446 residual_ = std::make_shared<ResidualType>();
447 setResidualSize_(*residual_);
454 void resetJacobian_()
458 jacobian_ = std::make_shared<JacobianMatrix>();
460 setJacobianPattern_(*jacobian_);
467 void maybeComputeColors_()
469 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
470 if (enableMultithreading_)
475 void checkAssemblerState_()
const
477 if (!isStationaryProblem_ && !prevSol_)
478 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
480 if (isStationaryProblem_ && prevSol_)
481 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
482 <<
" Did you forget to set the timeLoop to make this problem instationary?");
485 template<std::
size_t i,
class JacRow,
class SubRes>
486 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
489 assemble_(domainId, [&](
const auto& element)
491 MultiDomainAssemblerSubDomainView view{*
this, domainId};
492 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
493 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
497 template<std::
size_t i,
class SubRes>
498 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
501 assemble_(domainId, [&](
const auto& element)
503 MultiDomainAssemblerSubDomainView view{*
this, domainId};
504 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
505 subDomainAssembler.assembleResidual(subRes);
514 template<std::
size_t i,
class AssembleElementFunc>
515 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
518 bool succeeded =
false;
523 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
525 if (enableMultithreading_)
528 domainId, std::forward<AssembleElementFunc>(assembleElement)
537 for (
const auto& element : elements(
gridView(domainId)))
538 assembleElement(element);
544 catch (NumericalProblem &e)
546 std::cout <<
"rank " <<
gridView(domainId).comm().rank()
547 <<
" caught an exception while assembling:" << e.what()
553 if (
gridView(domainId).comm().size() > 1)
554 succeeded =
gridView(domainId).comm().min(succeeded);
558 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
562 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
563 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
564 Dune::index_constant<j> domainJ)
const
567 auto pattern = getJacobianPattern<isImplicit()>(gg);
573 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
574 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
575 Dune::index_constant<j> domainJ)
const
583 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
584 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
586 if constexpr (Dumux::Detail::hasPeriodicDofMap<GG>())
590 if (m.first < m.second)
592 auto& jac = jacRow[domainI];
595 res[m.first] += res[m.second];
597 const auto end = jac[m.second].end();
598 for (
auto it = jac[m.second].begin(); it != end; ++it)
599 jac[m.first][it.index()] += (*it);
603 res[m.second] = curSol[m.second] - curSol[m.first];
606 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
608 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
609 matrixBlock[eIdx][eIdx] = diagValue;
612 for (
auto it = jac[m.second].begin(); it != end; ++it)
614 auto& matrixBlock = *it;
617 assert(matrixBlock.N() == matrixBlock.M());
618 if(it.index() == m.second)
619 setMatrixBlock(matrixBlock, 1.0);
621 if(it.index() == m.first)
622 setMatrixBlock(matrixBlock, -1.0);
626 using namespace Dune::Hybrid;
629 auto& jacCoupling = jacRow[couplingDomainId];
631 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
632 jacCoupling[m.first][it.index()] += (*it);
634 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
643 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
644 void enforceProblemConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
646 if constexpr (Dumux::Detail::hasSubProblemGlobalConstraints<Problem<domainI>>())
648 auto& jac = jacRow[domainI];
649 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
654 (res)[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
656 auto& row = jac[dofIdx];
657 for (
auto col = row.begin(); col != row.end(); ++col)
658 row[col.index()][eqIdx] = 0.0;
660 jac[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
662 using namespace Dune::Hybrid;
665 auto& jacCoupling = jacRow[couplingDomainId];
667 auto& rowCoupling = jacCoupling[dofIdx];
668 for (
auto colCoupling = rowCoupling.begin(); colCoupling != rowCoupling.end(); ++colCoupling)
669 rowCoupling[colCoupling.index()][eqIdx] = 0.0;
674 for (
const auto& constraintData : this->
problem(domainI).constraints())
676 const auto& constraintInfo = constraintData.constraintInfo();
677 const auto& values = constraintData.values();
678 const auto dofIdx = constraintData.dofIndex();
680 for (
int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
682 if (constraintInfo.isConstraintEquation(eqIdx))
684 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
685 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
686 applyDirichletConstraint(dofIdx, values, eqIdx, pvIdx);
694 template<std::
size_t i,
class Res,
class GG,
class Sol>
695 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
697 if constexpr (Dumux::Detail::hasPeriodicDofMap<GG>())
701 if (m.first < m.second)
704 res[m.first] += res[m.second];
707 res[m.second] = curSol[m.second] - curSol[m.first];
714 template<std::
size_t i,
class Res,
class GG,
class Sol>
715 void enforceProblemConstraints_(Dune::index_constant<i> domainI, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
717 if constexpr (Dumux::Detail::hasSubProblemGlobalConstraints<Problem<domainI>>())
719 for (
const auto& constraintData : this->
problem(domainI).constraints())
721 const auto& constraintInfo = constraintData.constraintInfo();
722 const auto& values = constraintData.values();
723 const auto dofIdx = constraintData.dofIndex();
724 for (
int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
726 if (constraintInfo.isConstraintEquation(eqIdx))
728 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
729 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
730 res[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
738 ProblemTuple problemTuple_;
741 GridGeometryTuple gridGeometryTuple_;
744 GridVariablesTuple gridVariablesTuple_;
747 std::shared_ptr<const TimeLoop> timeLoop_;
753 bool isStationaryProblem_;
756 std::shared_ptr<JacobianMatrix> jacobian_;
757 std::shared_ptr<ResidualType> residual_;
760 mutable bool warningIssued_;
763 bool enableMultithreading_ =
false;
Subdomain-specific views on multidomain assemblers.
A linear system assembler (residual and Jacobian) for general schemes (FV, CVFE, FE,...
Definition: multidomain/assembler.hh:63
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/assembler.hh:86
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/assembler.hh:400
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/assembler.hh:76
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/assembler.hh:291
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/assembler.hh:407
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/assembler.hh:320
ResidualType & residual()
the full residual vector
Definition: multidomain/assembler.hh:376
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/assembler.hh:310
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/assembler.hh:79
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/assembler.hh:263
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/assembler.hh:364
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/assembler.hh:339
MultiDomainAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/assembler.hh:151
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/assembler.hh:387
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/assembler.hh:354
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/assembler.hh:236
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/assembler.hh:85
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/assembler.hh:412
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/assembler.hh:394
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/assembler.hh:372
typename MDTraits::Scalar Scalar
Definition: multidomain/assembler.hh:70
MDTraits Traits
Definition: multidomain/assembler.hh:68
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/assembler.hh:84
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/assembler.hh:278
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/assembler.hh:359
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/assembler.hh:93
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/assembler.hh:380
typename MDTraits::template SubDomain< id >::LocalResidual LocalResidual
Definition: multidomain/assembler.hh:73
CMType CouplingManager
Definition: multidomain/assembler.hh:88
MultiDomainAssembler(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/assembler.hh:179
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/assembler.hh:330
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/assembler.hh:344
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/assembler.hh:368
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/assembler.hh:229
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/assembler.hh:82
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/assembler.hh:349
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/assembler.hh:208
The cell-centered scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincclocalassembler.hh:195
The CVFE scheme multidomain local assembler.
Definition: experimental/assembly/subdomaincvfelocalassembler.hh:207
Subdomain-specific view on a multidomain assembler. Allows retrieval of sub-domain specific objects w...
Definition: assemblerview.hh:31
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
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
A linear system assembler (residual and Jacobian) for general grid-based schemes with multiple domain...
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
Definition: assembly/assembler.hh:44
bool allGridsSupportsMultithreading(const std::tuple< GG... > &gridGeometries)
Definition: multistagemultidomainfvassembler.hh:60
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 traits to detect periodicity support.
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...
Utilities for template meta programming.