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>
44namespace Grid::Capabilities {
48template<
class T, std::size_t... I>
82template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
85 template<std::
size_t id>
86 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
91 using Scalar =
typename MDTraits::Scalar;
94 template<std::
size_t id>
97 template<std::
size_t id>
98 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
100 template<std::
size_t id>
101 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
103 template<std::
size_t id>
104 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
116 {
return useImplicitAssembly; }
120 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
121 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
122 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
127 template<
class DiscretizationMethod, std::
size_t id>
128 struct SubDomainAssemblerType;
130 template<std::
size_t id>
131 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
136 template<std::
size_t id>
137 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
139 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
142 template<std::
size_t id,
class DM>
143 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
145 using type = SubDomainCVFELocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
148 template<std::
size_t id>
149 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
151 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
154 template<std::
size_t id>
155 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
157 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
160 template<std::
size_t id>
161 using SubDomainAssembler =
typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod,
id>::type;
176 , problemTuple_(std::move(
problem))
180 , isStationaryProblem_(true)
181 , warningIssued_(false)
183 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
184 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
189 && getParam<bool>(
"Assembly.Multithreading",
true);
191 maybeComputeColors_();
203 std::shared_ptr<const TimeLoop> timeLoop,
206 , problemTuple_(std::move(
problem))
209 , timeLoop_(timeLoop)
211 , isStationaryProblem_(false)
212 , warningIssued_(false)
214 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
219 && getParam<bool>(
"Assembly.Multithreading",
true);
221 maybeComputeColors_();
230 checkAssemblerState_();
234 using namespace Dune::Hybrid;
235 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
237 auto& jacRow = (*jacobian_)[domainId];
238 auto& subRes = (*residual_)[domainId];
239 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
241 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
242 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, curSol[domainId]);
258 checkAssemblerState_();
263 using namespace Dune::Hybrid;
264 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
266 auto& subRes = r[domainId];
267 this->assembleResidual_(domainId, subRes, curSol);
272 [[deprecated(
"Use the linear solver's norm. Will be deleted after 3.7")]]
276 Scalar resultSquared = 0.0;
279 using namespace Dune::Hybrid;
280 forEach(integralRange(Dune::Hybrid::size(
residual)), [&](
const auto domainId)
282 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
295 vectorHelper.makeNonOverlappingConsistent(
residual[domainId]);
298 else if (!warningIssued_)
301 std::cout <<
"\nWarning: norm calculation adds entries corresponding to\n"
302 <<
"overlapping entities multiple times. Please use the norm\n"
303 <<
"function provided by a linear solver instead." << std::endl;
305 warningIssued_ =
true;
312 localNormSquared =
gridView.comm().sum(localNormSquared);
315 resultSquared += localNormSquared;
319 return sqrt(resultSquared);
323 [[deprecated(
"Use norm(curSol) provided by the linear solver class instead. Will be deleted after 3.7")]]
338 std::shared_ptr<ResidualType> r)
344 setJacobianPattern_(*jacobian_);
345 setResidualSize_(*residual_);
354 jacobian_ = std::make_shared<JacobianMatrix>();
355 residual_ = std::make_shared<ResidualType>();
358 setJacobianPattern_(*jacobian_);
359 setResidualSize_(*residual_);
367 using namespace Dune::Hybrid;
368 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
370 forEach(jac[i], [&](
auto& jacBlock)
372 using BlockType = std::decay_t<
decltype(jacBlock)>;
373 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
374 jacBlock.setBuildMode(BlockType::BuildMode::random);
375 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
376 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
387 setJacobianPattern_();
388 maybeComputeColors_();
396 using namespace Dune::Hybrid;
397 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
398 { this->
gridVariables(domainId).update(curSol[domainId]); });
406 using namespace Dune::Hybrid;
407 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
408 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
412 template<std::
size_t i>
413 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
414 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
417 template<std::
size_t i>
418 const auto&
problem(Dune::index_constant<i> domainId)
const
419 {
return *std::get<domainId>(problemTuple_); }
422 template<std::
size_t i>
424 {
return *std::get<domainId>(gridGeometryTuple_); }
427 template<std::
size_t i>
428 const auto&
gridView(Dune::index_constant<i> domainId)
const
432 template<std::
size_t i>
434 {
return *std::get<domainId>(gridVariablesTuple_); }
437 template<std::
size_t i>
439 {
return *std::get<domainId>(gridVariablesTuple_); }
447 {
return *jacobian_; }
451 {
return *residual_; }
455 {
return *prevSol_; }
462 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
475 {
return isStationaryProblem_; }
480 template<std::
size_t i>
482 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
494 using namespace Dune::Hybrid;
495 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
497 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
499 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
500 pattern.exportIdx(jac[domainI][domainJ]);
510 using namespace Dune::Hybrid;
511 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
512 { res[domainId].resize(this->
numDofs(domainId)); });
516 void resetResidual_()
520 residual_ = std::make_shared<ResidualType>();
521 setResidualSize_(*residual_);
528 void resetJacobian_()
532 jacobian_ = std::make_shared<JacobianMatrix>();
534 setJacobianPattern_(*jacobian_);
541 void maybeComputeColors_()
543 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
544 if (enableMultithreading_)
549 void checkAssemblerState_()
const
551 if (!isStationaryProblem_ && !prevSol_)
552 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
554 if (isStationaryProblem_ && prevSol_)
555 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
556 <<
" Did you forget to set the timeLoop to make this problem instationary?");
559 template<std::
size_t i,
class JacRow,
class SubRes>
560 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
563 assemble_(domainId, [&](
const auto& element)
565 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
566 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
570 template<std::
size_t i,
class SubRes>
571 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
574 assemble_(domainId, [&](
const auto& element)
576 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
577 subDomainAssembler.assembleResidual(subRes);
586 template<std::
size_t i,
class AssembleElementFunc>
587 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
590 bool succeeded =
false;
595 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
597 if (enableMultithreading_)
600 domainId, std::forward<AssembleElementFunc>(assembleElement)
609 for (
const auto& element : elements(
gridView(domainId)))
610 assembleElement(element);
616 catch (NumericalProblem &e)
618 std::cout <<
"rank " <<
gridView(domainId).comm().rank()
619 <<
" caught an exception while assembling:" << e.what()
625 if (
gridView(domainId).comm().size() > 1)
626 succeeded =
gridView(domainId).comm().min(succeeded);
630 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
634 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
635 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
636 Dune::index_constant<j> domainJ)
const
639 auto pattern = getJacobianPattern<isImplicit()>(gg);
645 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
646 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
647 Dune::index_constant<j> domainJ)
const
655 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
656 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
662 if (m.first < m.second)
664 auto& jac = jacRow[domainI];
667 res[m.first] += res[m.second];
670 res[m.second] = curSol[m.second] - curSol[m.first];
672 const auto end = jac[m.second].end();
673 for (
auto it = jac[m.second].begin(); it != end; ++it)
674 jac[m.first][it.index()] += (*it);
677 for (
auto it = jac[m.second].begin(); it != end; ++it)
678 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
680 using namespace Dune::Hybrid;
683 auto& jacCoupling = jacRow[couplingDomainId];
685 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
686 jacCoupling[m.first][it.index()] += (*it);
688 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
697 ProblemTuple problemTuple_;
700 GridGeometryTuple gridGeometryTuple_;
703 GridVariablesTuple gridVariablesTuple_;
706 std::shared_ptr<const TimeLoop> timeLoop_;
712 bool isStationaryProblem_;
715 std::shared_ptr<JacobianMatrix> jacobian_;
716 std::shared_ptr<ResidualType> residual_;
719 mutable bool warningIssued_;
722 bool enableMultithreading_ =
false;
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:84
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:413
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:394
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:384
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:107
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:104
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:446
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:115
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:254
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:106
CMType CouplingManager
Definition: multidomain/fvassembler.hh:110
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:247
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:404
Scalar residualNorm(const SolutionVector &curSol)
compute the residual and return it's vector norm
Definition: multidomain/fvassembler.hh:324
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:438
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:468
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:98
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:423
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:352
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:428
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:101
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:442
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:461
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:454
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:433
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:486
MDTraits Traits
Definition: multidomain/fvassembler.hh:89
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:199
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:365
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:95
Scalar normOfResidual(ResidualType &residual) const
compute a residual's vector norm (this is a temporary interface introduced during the deprecation per...
Definition: multidomain/fvassembler.hh:273
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:481
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:418
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/fvassembler.hh:108
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:171
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:228
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:474
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:91
ResidualType & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:450
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:337
Definition: parallelhelpers.hh:476
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:270
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:56
The default time loop for instationary simulations.
Definition: common/timeloop.hh:101
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:267
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
constexpr FCStaggered fcstaggered
Definition: method.hh:151
bool allGridsSupportsMultithreadingImpl(const T &gridGeometries, std::index_sequence< I... >)
Definition: multidomain/fvassembler.hh:49
bool allGridsSupportsMultithreading(const std::tuple< GG... > &gridGeometries)
Definition: multidomain/fvassembler.hh:57
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.
trait that is specialized for coupling manager supporting multithreaded assembly
Definition: multidomain/fvassembler.hh:71
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...
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.