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>
46namespace Grid::Capabilities {
50template<
class T, std::size_t... I>
72struct CouplingManagerSupportsMultithreadedAssembly :
public std::false_type
84template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
87 template<std::
size_t id>
88 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
93 using Scalar =
typename MDTraits::Scalar;
96 template<std::
size_t id>
99 template<std::
size_t id>
100 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
102 template<std::
size_t id>
103 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
105 template<std::
size_t id>
106 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
118 {
return useImplicitAssembly; }
122 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
123 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
124 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
129 template<std::
size_t id>
132 template<
class DiscretizationMethod, std::
size_t id>
133 struct SubDomainAssemblerType;
135 template<std::
size_t id>
136 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
141 template<std::
size_t id>
142 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
144 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
147 template<std::
size_t id,
class DM>
148 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
150 using type = SubDomainCVFELocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
153 template<std::
size_t id>
154 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
156 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
159 template<std::
size_t id>
160 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
162 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
165 template<std::
size_t id>
166 using SubDomainAssembler =
typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod,
id>::type;
181 , problemTuple_(std::move(
problem))
185 , isStationaryProblem_(true)
186 , warningIssued_(false)
188 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
189 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
194 && getParam<bool>(
"Assembly.Multithreading",
true);
196 maybeComputeColors_();
208 std::shared_ptr<const TimeLoop> timeLoop,
211 , problemTuple_(std::move(
problem))
214 , timeLoop_(timeLoop)
216 , isStationaryProblem_(false)
217 , warningIssued_(false)
219 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
224 && getParam<bool>(
"Assembly.Multithreading",
true);
226 maybeComputeColors_();
235 checkAssemblerState_();
239 using namespace Dune::Hybrid;
240 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
242 auto& jacRow = (*jacobian_)[domainId];
243 auto& subRes = (*residual_)[domainId];
244 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
246 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
247 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, curSol[domainId]);
263 checkAssemblerState_();
268 using namespace Dune::Hybrid;
269 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
271 auto& subRes = r[domainId];
272 this->assembleResidual_(domainId, subRes, curSol);
282 std::shared_ptr<ResidualType> r)
288 setJacobianPattern_(*jacobian_);
289 setResidualSize_(*residual_);
298 jacobian_ = std::make_shared<JacobianMatrix>();
299 residual_ = std::make_shared<ResidualType>();
302 setJacobianPattern_(*jacobian_);
303 setResidualSize_(*residual_);
311 using namespace Dune::Hybrid;
312 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
314 forEach(jac[i], [&](
auto& jacBlock)
316 using BlockType = std::decay_t<
decltype(jacBlock)>;
317 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
318 jacBlock.setBuildMode(BlockType::BuildMode::random);
319 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
320 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
330 setJacobianPattern_(*jacobian_);
331 setResidualSize_(*residual_);
332 maybeComputeColors_();
340 using namespace Dune::Hybrid;
341 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
342 { this->
gridVariables(domainId).update(curSol[domainId]); });
350 using namespace Dune::Hybrid;
351 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
352 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
356 template<std::
size_t i>
357 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
358 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
361 template<std::
size_t i>
362 const auto&
problem(Dune::index_constant<i> domainId)
const
363 {
return *std::get<domainId>(problemTuple_); }
366 template<std::
size_t i>
368 {
return *std::get<domainId>(gridGeometryTuple_); }
371 template<std::
size_t i>
372 const auto&
gridView(Dune::index_constant<i> domainId)
const
376 template<std::
size_t i>
378 {
return *std::get<domainId>(gridVariablesTuple_); }
381 template<std::
size_t i>
383 {
return *std::get<domainId>(gridVariablesTuple_); }
391 {
return *jacobian_; }
395 {
return *residual_; }
399 {
return *prevSol_; }
406 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
419 {
return isStationaryProblem_; }
424 template<std::
size_t i>
426 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
438 using namespace Dune::Hybrid;
439 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
441 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
443 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
444 pattern.exportIdx(jac[domainI][domainJ]);
454 using namespace Dune::Hybrid;
455 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
456 { res[domainId].resize(this->
numDofs(domainId)); });
460 void resetResidual_()
464 residual_ = std::make_shared<ResidualType>();
465 setResidualSize_(*residual_);
472 void resetJacobian_()
476 jacobian_ = std::make_shared<JacobianMatrix>();
478 setJacobianPattern_(*jacobian_);
485 void maybeComputeColors_()
487 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
488 if (enableMultithreading_)
493 void checkAssemblerState_()
const
495 if (!isStationaryProblem_ && !prevSol_)
496 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
498 if (isStationaryProblem_ && prevSol_)
499 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
500 <<
" Did you forget to set the timeLoop to make this problem instationary?");
503 template<std::
size_t i,
class JacRow,
class SubRes>
504 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
507 assemble_(domainId, [&](
const auto& element)
509 MultiDomainAssemblerSubDomainView view{*
this, domainId};
510 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
511 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
515 template<std::
size_t i,
class SubRes>
516 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
519 assemble_(domainId, [&](
const auto& element)
521 MultiDomainAssemblerSubDomainView view{*
this, domainId};
522 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
523 subDomainAssembler.assembleResidual(subRes);
532 template<std::
size_t i,
class AssembleElementFunc>
533 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
536 bool succeeded =
false;
541 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
543 if (enableMultithreading_)
546 domainId, std::forward<AssembleElementFunc>(assembleElement)
555 for (
const auto& element : elements(
gridView(domainId)))
556 assembleElement(element);
562 catch (NumericalProblem &e)
564 std::cout <<
"rank " <<
gridView(domainId).comm().rank()
565 <<
" caught an exception while assembling:" << e.what()
571 if (
gridView(domainId).comm().size() > 1)
572 succeeded =
gridView(domainId).comm().min(succeeded);
576 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
580 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
581 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
582 Dune::index_constant<j> domainJ)
const
585 auto pattern = getJacobianPattern<isImplicit()>(gg);
591 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
592 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
593 Dune::index_constant<j> domainJ)
const
601 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
602 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
604 if constexpr (Detail::hasPeriodicDofMap<GG>())
608 if (m.first < m.second)
610 auto& jac = jacRow[domainI];
613 res[m.first] += res[m.second];
615 const auto end = jac[m.second].end();
616 for (
auto it = jac[m.second].begin(); it != end; ++it)
617 jac[m.first][it.index()] += (*it);
621 res[m.second] = curSol[m.second] - curSol[m.first];
624 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
626 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
627 matrixBlock[eIdx][eIdx] = diagValue;
630 for (
auto it = jac[m.second].begin(); it != end; ++it)
632 auto& matrixBlock = *it;
635 assert(matrixBlock.N() == matrixBlock.M());
636 if(it.index() == m.second)
637 setMatrixBlock(matrixBlock, 1.0);
639 if(it.index() == m.first)
640 setMatrixBlock(matrixBlock, -1.0);
644 using namespace Dune::Hybrid;
647 auto& jacCoupling = jacRow[couplingDomainId];
649 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
650 jacCoupling[m.first][it.index()] += (*it);
652 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
661 ProblemTuple problemTuple_;
664 GridGeometryTuple gridGeometryTuple_;
667 GridVariablesTuple gridVariablesTuple_;
670 std::shared_ptr<const TimeLoop> timeLoop_;
676 bool isStationaryProblem_;
679 std::shared_ptr<JacobianMatrix> jacobian_;
680 std::shared_ptr<ResidualType> residual_;
683 mutable bool warningIssued_;
686 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:86
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:357
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:338
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:328
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:109
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:106
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:390
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:117
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:259
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:108
CMType CouplingManager
Definition: multidomain/fvassembler.hh:112
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:252
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:348
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:382
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:412
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:100
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:367
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:296
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:372
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:103
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:386
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:405
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:398
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:377
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:430
MDTraits Traits
Definition: multidomain/fvassembler.hh:91
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:204
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:309
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:97
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:425
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:362
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/fvassembler.hh:110
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:176
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:233
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:418
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:93
ResidualType & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:394
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:281
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...
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 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...
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
Utilities for template meta programming.