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>
60template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
63 template<std::
size_t id>
64 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
69 using Scalar =
typename MDTraits::Scalar;
72 template<std::
size_t id>
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<DiscretizationMethod discMethod, std::
size_t id>
106 struct SubDomainAssemblerType;
108 template<std::
size_t id>
114 template<std::
size_t id>
117 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
120 template<std::
size_t id>
123 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
126 template<std::
size_t id>
129 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
132 template<std::
size_t id>
133 using SubDomainAssembler =
typename SubDomainAssemblerType<GridGeometry<id>::discMethod,
id>::type;
152 , isStationaryProblem_(true)
154 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
155 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
162 [[deprecated(
"Please use constructor taking the previous solution instead. Will be removed after release 3.2!")]]
167 std::shared_ptr<const TimeLoop> timeLoop)
172 , timeLoop_(timeLoop)
173 , isStationaryProblem_(false)
175 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
187 std::shared_ptr<const TimeLoop> timeLoop,
193 , timeLoop_(timeLoop)
195 , isStationaryProblem_(false)
197 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
200 template<
class DeprecatedGridGeometryTuple,
class DeprecatedGridVariablesTuple,
202 bool isDeprecated = !std::is_convertible<DeprecatedGridVariablesTuple, GridVariablesTuple>::value,
203 std::enable_if_t<isStaggered && isDeprecated, int> = 0>
204 [[deprecated(
"Please change the order within your tuples for gridGeometry and gridVariables. Will be removed after release 3.2!")]]
213 template<
class DeprecatedGridGeometryTuple,
class DeprecatedGridVariablesTuple,
215 bool isDeprecated = !std::is_convertible<DeprecatedGridVariablesTuple, GridVariablesTuple>::value,
216 std::enable_if_t<isStaggered && isDeprecated, int> = 0>
217 [[deprecated(
"Please change the order within your tuples for gridGeometry and gridVariables. Will be removed after release 3.2!")]]
222 std::shared_ptr<const TimeLoop> timeLoop,
228 template<
class DeprecatedGridGeometryTuple,
class DeprecatedGridVariablesTuple,
230 bool isDeprecated = !std::is_convertible<DeprecatedGridVariablesTuple, GridVariablesTuple>::value,
231 std::enable_if_t<isStaggered && isDeprecated, int> = 0>
232 [[deprecated(
"Please change the order within your tuples for gridGeometry and gridVariables. Will be removed after release 3.2!")]]
237 std::shared_ptr<const TimeLoop> timeLoop)
248 checkAssemblerState_();
252 using namespace Dune::Hybrid;
253 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
255 auto& jacRow = (*jacobian_)[domainId];
256 auto& subRes = (*residual_)[domainId];
257 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
271 checkAssemblerState_();
276 using namespace Dune::Hybrid;
277 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
279 auto& subRes = r[domainId];
280 this->assembleResidual_(domainId, subRes, curSol);
295 return sqrt(result2);
304 std::shared_ptr<SolutionVector> r)
320 jacobian_ = std::make_shared<JacobianMatrix>();
321 residual_ = std::make_shared<SolutionVector>();
333 using namespace Dune::Hybrid;
334 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
336 forEach(jac[i], [&](
auto& jacBlock)
338 using BlockType = std::decay_t<
decltype(jacBlock)>;
339 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
340 jacBlock.setBuildMode(BlockType::BuildMode::random);
341 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
342 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
352 using namespace Dune::Hybrid;
353 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
355 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
357 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
358 pattern.exportIdx(jac[domainI][domainJ]);
368 using namespace Dune::Hybrid;
369 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
370 { res[domainId].resize(this->
numDofs(domainId)); });
378 using namespace Dune::Hybrid;
379 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
380 { this->
gridVariables(domainId).update(curSol[domainId]); });
388 using namespace Dune::Hybrid;
389 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
390 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
394 template<std::
size_t i>
395 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
396 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
399 template<std::
size_t i>
400 const auto&
problem(Dune::index_constant<i> domainId)
const
401 {
return *std::get<domainId>(problemTuple_); }
404 template<std::
size_t i>
406 {
return *std::get<domainId>(gridGeometryTuple_); }
409 template<std::
size_t i>
410 const auto&
gridView(Dune::index_constant<i> domainId)
const
414 template<std::
size_t i>
416 {
return *std::get<domainId>(gridVariablesTuple_); }
419 template<std::
size_t i>
421 {
return *std::get<domainId>(gridVariablesTuple_); }
429 {
return *jacobian_; }
433 {
return *residual_; }
437 {
return *prevSol_; }
444 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
457 {
return isStationaryProblem_; }
462 template<std::
size_t i>
464 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
472 void resetResidual_()
476 residual_ = std::make_shared<SolutionVector>();
484 void resetJacobian_()
488 jacobian_ = std::make_shared<JacobianMatrix>();
497 void checkAssemblerState_()
const
499 if (!isStationaryProblem_ && !prevSol_)
500 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
502 if (isStationaryProblem_ && prevSol_)
503 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
504 <<
" Did you forget to set the timeLoop to make this problem instationary?");
507 template<std::
size_t i,
class JacRow,
class SubRes>
508 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
511 assemble_(domainId, [&](
const auto& element)
513 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
514 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
518 template<std::
size_t i,
class SubRes>
519 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
522 assemble_(domainId, [&](
const auto& element)
524 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
525 subDomainAssembler.assembleResidual(subRes);
535 template<std::
size_t i,
class AssembleElementFunc>
536 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
539 for (
const auto& element : elements(
gridView(domainId)))
540 assembleElement(element);
544 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
545 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
546 Dune::index_constant<j> domainJ)
const
549 auto pattern = getJacobianPattern<isImplicit()>(gg);
555 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
556 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
557 Dune::index_constant<j> domainJ)
const
564 template<
class DeprecatedGridGeometryTuple,
class DeprecatedGridVariablesTuple, std::size_t... I>
565 [[deprecated(
"\n\n Please change the order within your tuples for gridGeometry and gridVariables. Usage:\n\n"
566 "auto assembler = std::make_shared<Assembler>(std::make_tuple(ffProblem, ffProblem, otherProblem, ...),\n"
567 " std::make_tuple(ffGridGeometry->faceFVGridGeometryPtr(),\n"
568 " ffFvGridGeometry->cellCenterFVGridGeometryPtr(),\n"
569 " otherFvGridGeometry, ...),\n"
570 " std::make_tuple(ffGridVariables->faceGridVariablesPtr(),\n"
571 " ffGridVariables->cellCenterGridVariablesPtr(),\n"
572 " otherGridVariables, ...),\n"
573 " couplingManager);\n\n"
574 "Will be removed after release 3.2!")]]
579 std::index_sequence<I...>)
582 , gridGeometryTuple_(std::make_tuple(std::move(std::get<1>(
gridGeometry)),
585 , gridVariablesTuple_(std::make_tuple(std::move(std::get<1>(
gridVariables)),
589 , isStationaryProblem_(true)
591 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
592 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
595 template<
class DeprecatedGridGeometryTuple,
class DeprecatedGridVariablesTuple, std::size_t... I>
596 [[deprecated(
"\n\n Please change the order within your tuples for gridGeometry and gridVariables. Usage:\n\n"
597 "auto assembler = std::make_shared<Assembler>(std::make_tuple(ffProblem, ffProblem, otherProblem, ...),\n"
598 " std::make_tuple(ffGridGeometry->faceFVGridGeometryPtr(),\n"
599 " ffFvGridGeometry->cellCenterFVGridGeometryPtr(),\n"
600 " otherFvGridGeometry, ...),\n"
601 " std::make_tuple(ffGridVariables->faceGridVariablesPtr(),\n"
602 " ffGridVariables->cellCenterGridVariablesPtr(),\n"
603 " otherGridVariables, ...),\n"
604 " couplingManager,\n"
605 " timeLoop, solOld);\n\n"
606 "Will be removed after release 3.2!")]]
611 std::shared_ptr<const TimeLoop> timeLoop,
613 std::index_sequence<I...>)
616 , gridGeometryTuple_(std::make_tuple(std::move(std::get<1>(
gridGeometry)),
619 , gridVariablesTuple_(std::make_tuple(std::move(std::get<1>(
gridVariables)),
622 , timeLoop_(timeLoop)
624 , isStationaryProblem_(false)
626 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
630 ProblemTuple problemTuple_;
633 GridGeometryTuple gridGeometryTuple_;
636 GridVariablesTuple gridVariablesTuple_;
639 std::shared_ptr<const TimeLoop> timeLoop_;
645 bool isStationaryProblem_;
648 std::shared_ptr<JacobianMatrix> jacobian_;
649 std::shared_ptr<SolutionVector> residual_;
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.
Manages the handling of time dependent problems.
The available discretization methods in Dumux.
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...
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
constexpr auto makeIndexSequenceWithOffset()
Definition: utility.hh:96
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: common/properties.hh:89
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
The default time loop for instationary simulations.
Definition: timeloop.hh:117
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:62
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:395
MultiDomainFVAssembler(ProblemTuple &&problem, DeprecatedGridGeometryTuple &&gridGeometry, DeprecatedGridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
Definition: multidomain/fvassembler.hh:218
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:376
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:85
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:86
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:82
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:428
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:93
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:269
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:84
CMType CouplingManager
Definition: multidomain/fvassembler.hh:88
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:262
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:386
Scalar residualNorm(const SolutionVector &curSol)
Definition: multidomain/fvassembler.hh:286
MultiDomainFVAssembler(ProblemTuple &&problem, DeprecatedGridGeometryTuple &&gridGeometry, DeprecatedGridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop)
Definition: multidomain/fvassembler.hh:233
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:420
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:450
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:76
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:303
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:405
void setResidualSize(SolutionVector &res) const
Resizes the residual.
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:318
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:410
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:183
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:79
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop)
The constructor for instationary problems.
Definition: multidomain/fvassembler.hh:163
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:424
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:443
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:436
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:415
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:468
MDTraits Traits
Definition: multidomain/fvassembler.hh:67
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:331
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:73
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:143
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:463
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:400
MultiDomainFVAssembler(ProblemTuple &&problem, DeprecatedGridGeometryTuple &&gridGeometry, DeprecatedGridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
Definition: multidomain/fvassembler.hh:205
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:246
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:456
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:69
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:350
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:432
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:268
Declares all properties used in Dumux.