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;
167 std::shared_ptr<const TimeLoop> timeLoop,
173 , timeLoop_(timeLoop)
175 , isStationaryProblem_(false)
177 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
186 checkAssemblerState_();
190 using namespace Dune::Hybrid;
191 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
193 auto& jacRow = (*jacobian_)[domainId];
194 auto& subRes = (*residual_)[domainId];
195 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
209 checkAssemblerState_();
214 using namespace Dune::Hybrid;
215 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
217 auto& subRes = r[domainId];
218 this->assembleResidual_(domainId, subRes, curSol);
233 return sqrt(result2);
242 std::shared_ptr<SolutionVector> r)
258 jacobian_ = std::make_shared<JacobianMatrix>();
259 residual_ = std::make_shared<SolutionVector>();
271 using namespace Dune::Hybrid;
272 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
274 forEach(jac[i], [&](
auto& jacBlock)
276 using BlockType = std::decay_t<
decltype(jacBlock)>;
277 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
278 jacBlock.setBuildMode(BlockType::BuildMode::random);
279 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
280 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
290 using namespace Dune::Hybrid;
291 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
293 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
295 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
296 pattern.exportIdx(jac[domainI][domainJ]);
306 using namespace Dune::Hybrid;
307 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
308 { res[domainId].resize(this->
numDofs(domainId)); });
316 using namespace Dune::Hybrid;
317 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
318 { this->
gridVariables(domainId).update(curSol[domainId]); });
326 using namespace Dune::Hybrid;
327 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
328 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
332 template<std::
size_t i>
333 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
334 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
337 template<std::
size_t i>
338 const auto&
problem(Dune::index_constant<i> domainId)
const
339 {
return *std::get<domainId>(problemTuple_); }
342 template<std::
size_t i>
344 {
return *std::get<domainId>(gridGeometryTuple_); }
347 template<std::
size_t i>
348 const auto&
gridView(Dune::index_constant<i> domainId)
const
352 template<std::
size_t i>
354 {
return *std::get<domainId>(gridVariablesTuple_); }
357 template<std::
size_t i>
359 {
return *std::get<domainId>(gridVariablesTuple_); }
367 {
return *jacobian_; }
371 {
return *residual_; }
375 {
return *prevSol_; }
382 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
395 {
return isStationaryProblem_; }
400 template<std::
size_t i>
402 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
410 void resetResidual_()
414 residual_ = std::make_shared<SolutionVector>();
422 void resetJacobian_()
426 jacobian_ = std::make_shared<JacobianMatrix>();
435 void checkAssemblerState_()
const
437 if (!isStationaryProblem_ && !prevSol_)
438 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
440 if (isStationaryProblem_ && prevSol_)
441 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
442 <<
" Did you forget to set the timeLoop to make this problem instationary?");
445 template<std::
size_t i,
class JacRow,
class SubRes>
446 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
449 assemble_(domainId, [&](
const auto& element)
451 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
452 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
456 template<std::
size_t i,
class SubRes>
457 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
460 assemble_(domainId, [&](
const auto& element)
462 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
463 subDomainAssembler.assembleResidual(subRes);
473 template<std::
size_t i,
class AssembleElementFunc>
474 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
477 for (
const auto& element : elements(
gridView(domainId)))
478 assembleElement(element);
482 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
483 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
484 Dune::index_constant<j> domainJ)
const
487 auto pattern = getJacobianPattern<isImplicit()>(gg);
493 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
494 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
495 Dune::index_constant<j> domainJ)
const
503 ProblemTuple problemTuple_;
506 GridGeometryTuple gridGeometryTuple_;
509 GridVariablesTuple gridVariablesTuple_;
512 std::shared_ptr<const TimeLoop> timeLoop_;
518 bool isStationaryProblem_;
521 std::shared_ptr<JacobianMatrix> jacobian_;
522 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.
The available discretization methods in Dumux.
Helper function to generate Jacobian pattern for multi domain models.
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
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...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:77
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
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:333
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:314
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:366
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:207
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:200
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:324
Scalar residualNorm(const SolutionVector &curSol)
Definition: multidomain/fvassembler.hh:224
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:358
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:388
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:241
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:343
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:304
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:256
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:348
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:163
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:79
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:362
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:381
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:374
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:353
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:406
MDTraits Traits
Definition: multidomain/fvassembler.hh:67
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:269
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:401
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:338
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:184
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:394
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:69
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:288
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:370
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:268
Declares all properties used in Dumux.
Manages the handling of time dependent problems.