26#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
32#include <dune/common/hybridutilities.hh>
33#include <dune/istl/matrixindexset.hh>
58namespace Grid::Capabilities {
62template<
class T, std::size_t... I>
96template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
99 template<std::
size_t id>
100 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
105 using Scalar =
typename MDTraits::Scalar;
108 template<std::
size_t id>
111 template<std::
size_t id>
112 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
114 template<std::
size_t id>
115 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
117 template<std::
size_t id>
118 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
130 {
return useImplicitAssembly; }
134 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
135 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
136 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
141 template<
class DiscretizationMethod, std::
size_t id>
142 struct SubDomainAssemblerType;
144 template<std::
size_t id>
145 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
150 template<std::
size_t id>
151 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
153 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
156 template<std::
size_t id>
157 struct SubDomainAssemblerType<DiscretizationMethods::
Box, id>
159 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
162 template<std::
size_t id>
163 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
165 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
168 template<std::
size_t id>
169 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
171 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
174 template<std::
size_t id>
175 struct SubDomainAssemblerType<DiscretizationMethods::
FCDiamond, id>
177 using type = SubDomainFaceCenteredDiamondLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
180 template<std::
size_t id>
181 struct SubDomainAssemblerType<DiscretizationMethods::
PQ1Bubble, id>
183 using type = SubDomainPQ1BubbleLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
186 template<std::
size_t id>
187 using SubDomainAssembler =
typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod,
id>::type;
202 , problemTuple_(std::move(
problem))
206 , isStationaryProblem_(true)
207 , warningIssued_(false)
209 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
210 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
215 && getParam<bool>(
"Assembly.Multithreading",
true);
217 maybeComputeColors_();
229 std::shared_ptr<const TimeLoop> timeLoop,
232 , problemTuple_(std::move(
problem))
235 , timeLoop_(timeLoop)
237 , isStationaryProblem_(false)
238 , warningIssued_(false)
240 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
245 && getParam<bool>(
"Assembly.Multithreading",
true);
247 maybeComputeColors_();
256 checkAssemblerState_();
260 using namespace Dune::Hybrid;
261 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
263 auto& jacRow = (*jacobian_)[domainId];
264 auto& subRes = (*residual_)[domainId];
265 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
267 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
268 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, curSol[domainId]);
284 checkAssemblerState_();
289 using namespace Dune::Hybrid;
290 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
292 auto& subRes = r[domainId];
293 this->assembleResidual_(domainId, subRes, curSol);
305 Scalar resultSquared = 0.0;
308 using namespace Dune::Hybrid;
309 forEach(integralRange(Dune::Hybrid::size(
residual)), [&](
const auto domainId)
311 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
324 vectorHelper.makeNonOverlappingConsistent(
residual[domainId]);
327 else if (!warningIssued_)
330 std::cout <<
"\nWarning: norm calculation adds entries corresponding to\n"
331 <<
"overlapping entities multiple times. Please use the norm\n"
332 <<
"function provided by a linear solver instead." << std::endl;
334 warningIssued_ =
true;
341 localNormSquared =
gridView.comm().sum(localNormSquared);
344 resultSquared += localNormSquared;
348 return sqrt(resultSquared);
357 std::shared_ptr<SolutionVector> r)
363 setJacobianPattern_(*jacobian_);
364 setResidualSize_(*residual_);
373 jacobian_ = std::make_shared<JacobianMatrix>();
374 residual_ = std::make_shared<SolutionVector>();
377 setJacobianPattern_(*jacobian_);
378 setResidualSize_(*residual_);
386 using namespace Dune::Hybrid;
387 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
389 forEach(jac[i], [&](
auto& jacBlock)
391 using BlockType = std::decay_t<
decltype(jacBlock)>;
392 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
393 jacBlock.setBuildMode(BlockType::BuildMode::random);
394 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
395 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
406 setJacobianPattern_();
407 maybeComputeColors_();
415 using namespace Dune::Hybrid;
416 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
417 { this->
gridVariables(domainId).update(curSol[domainId]); });
425 using namespace Dune::Hybrid;
426 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
427 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
431 template<std::
size_t i>
432 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
433 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
436 template<std::
size_t i>
437 const auto&
problem(Dune::index_constant<i> domainId)
const
438 {
return *std::get<domainId>(problemTuple_); }
441 template<std::
size_t i>
443 {
return *std::get<domainId>(gridGeometryTuple_); }
446 template<std::
size_t i>
447 const auto&
gridView(Dune::index_constant<i> domainId)
const
451 template<std::
size_t i>
453 {
return *std::get<domainId>(gridVariablesTuple_); }
456 template<std::
size_t i>
458 {
return *std::get<domainId>(gridVariablesTuple_); }
466 {
return *jacobian_; }
470 {
return *residual_; }
474 {
return *prevSol_; }
481 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
494 {
return isStationaryProblem_; }
499 template<std::
size_t i>
501 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
513 using namespace Dune::Hybrid;
514 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
516 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
518 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
519 pattern.exportIdx(jac[domainI][domainJ]);
529 using namespace Dune::Hybrid;
530 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
531 { res[domainId].resize(this->
numDofs(domainId)); });
535 void resetResidual_()
539 residual_ = std::make_shared<SolutionVector>();
540 setResidualSize_(*residual_);
547 void resetJacobian_()
551 jacobian_ = std::make_shared<JacobianMatrix>();
553 setJacobianPattern_(*jacobian_);
560 void maybeComputeColors_()
562 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
563 if (enableMultithreading_)
568 void checkAssemblerState_()
const
570 if (!isStationaryProblem_ && !prevSol_)
571 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
573 if (isStationaryProblem_ && prevSol_)
574 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
575 <<
" Did you forget to set the timeLoop to make this problem instationary?");
578 template<std::
size_t i,
class JacRow,
class SubRes>
579 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
582 assemble_(domainId, [&](
const auto& element)
584 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
585 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
589 template<std::
size_t i,
class SubRes>
590 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
593 assemble_(domainId, [&](
const auto& element)
595 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
596 subDomainAssembler.assembleResidual(subRes);
605 template<std::
size_t i,
class AssembleElementFunc>
606 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
609 bool succeeded =
false;
614 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
616 if (enableMultithreading_)
619 domainId, std::forward<AssembleElementFunc>(assembleElement)
628 for (
const auto& element : elements(
gridView(domainId)))
629 assembleElement(element);
635 catch (NumericalProblem &e)
637 std::cout <<
"rank " <<
gridView(domainId).comm().rank()
638 <<
" caught an exception while assembling:" << e.what()
644 if (
gridView(domainId).comm().size() > 1)
645 succeeded =
gridView(domainId).comm().min(succeeded);
649 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
653 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
654 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
655 Dune::index_constant<j> domainJ)
const
658 auto pattern = getJacobianPattern<isImplicit()>(gg);
664 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
665 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
666 Dune::index_constant<j> domainJ)
const
674 template<std::
size_t i,
class JacRow,
class Sol,
class GG>
675 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Sol& res,
const GG&
gridGeometry,
const Sol& curSol)
681 if (m.first < m.second)
683 auto& jac = jacRow[domainI];
686 res[m.first] += res[m.second];
689 res[m.second] = curSol[m.second] - curSol[m.first];
691 const auto end = jac[m.second].end();
692 for (
auto it = jac[m.second].begin(); it != end; ++it)
693 jac[m.first][it.index()] += (*it);
696 for (
auto it = jac[m.second].begin(); it != end; ++it)
697 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
699 using namespace Dune::Hybrid;
702 auto& jacCoupling = jacRow[couplingDomainId];
704 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
705 jacCoupling[m.first][it.index()] += (*it);
707 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
716 ProblemTuple problemTuple_;
719 GridGeometryTuple gridGeometryTuple_;
722 GridVariablesTuple gridVariablesTuple_;
725 std::shared_ptr<const TimeLoop> timeLoop_;
731 bool isStationaryProblem_;
734 std::shared_ptr<JacobianMatrix> jacobian_;
735 std::shared_ptr<SolutionVector> residual_;
741 bool enableMultithreading_ =
false;
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.
Utilities for template meta programming.
Some exceptions thrown in DuMux
dune-grid capabilities compatibility layer
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
Helper function to generate Jacobian pattern for multi domain models.
An assembler for Jacobian and residual contribution per element (box methods) for multidomain problem...
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
An assembler for Jacobian and residual contribution per element (face-centered diamond methods) for m...
An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for...
An assembler for Jacobian and residual contribution per element for multidomain problems.
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:57
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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:71
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
bool allGridsSupportsMultithreading(const std::tuple< GG... > &gridGeometries)
Definition: multidomain/fvassembler.hh:71
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:86
CVFE< CVFEMethods::CR_RT > FCDiamond
Definition: method.hh:90
constexpr Box box
Definition: method.hh:136
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
CVFE< CVFEMethods::PQ1Bubble > PQ1Bubble
Definition: method.hh:97
constexpr FCStaggered fcstaggered
Definition: method.hh:140
bool allGridsSupportsMultithreadingImpl(const T &gridGeometries, std::index_sequence< I... >)
Definition: multidomain/fvassembler.hh:63
Definition: common/properties.hh:72
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:68
The default time loop for instationary simulations.
Definition: common/timeloop.hh:113
Definition: parallelhelpers.hh:473
trait that is specialized for coupling manager supporting multithreaded assembly
Definition: multidomain/fvassembler.hh:85
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:98
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:432
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:413
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:403
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:121
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:122
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:118
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:465
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:129
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:280
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:120
CMType CouplingManager
Definition: multidomain/fvassembler.hh:124
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:273
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:423
Scalar residualNorm(const SolutionVector &curSol)
compute the residual and return it's vector norm
Definition: multidomain/fvassembler.hh:298
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:457
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:487
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:112
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: multidomain/fvassembler.hh:356
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:442
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:371
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:447
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:115
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:461
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:480
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:473
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:452
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:505
MDTraits Traits
Definition: multidomain/fvassembler.hh:103
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:225
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:384
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:109
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:500
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:437
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:197
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:254
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:493
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:105
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:469
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Declares all properties used in Dumux.
Manages the handling of time dependent problems.