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>
57template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
60 template<std::
size_t id>
61 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
66 using Scalar =
typename MDTraits::Scalar;
69 template<std::
size_t id>
72 template<std::
size_t id>
73 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
75 template<std::
size_t id>
76 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
78 template<std::
size_t id>
79 using FVGridGeometry [[deprecated(
"Use GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] =
GridGeometry<id>;
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;
206 checkAssemblerState_();
210 using namespace Dune::Hybrid;
211 forEach(integralRange(Dune::Hybrid::size(*jacobian_)), [&](
const auto domainId)
213 auto& jacRow = (*jacobian_)[domainId];
214 auto& subRes = (*residual_)[domainId];
215 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
229 checkAssemblerState_();
234 using namespace Dune::Hybrid;
235 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
237 auto& subRes = r[domainId];
238 this->assembleResidual_(domainId, subRes, curSol);
253 return sqrt(result2);
262 std::shared_ptr<SolutionVector> r)
278 jacobian_ = std::make_shared<JacobianMatrix>();
279 residual_ = std::make_shared<SolutionVector>();
291 using namespace Dune::Hybrid;
292 forEach(jac, [](
auto& jacRow)
294 forEach(jacRow, [](
auto& jacBlock)
296 using BlockType = std::decay_t<
decltype(jacBlock)>;
297 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
298 jacBlock.setBuildMode(BlockType::BuildMode::random);
299 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
300 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
310 using namespace Dune::Hybrid;
311 forEach(integralRange(Dune::Hybrid::size(jac)), [&](
const auto domainI)
313 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
315 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
316 pattern.exportIdx(jac[domainI][domainJ]);
326 using namespace Dune::Hybrid;
327 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
328 { res[domainId].resize(this->
numDofs(domainId)); });
336 using namespace Dune::Hybrid;
337 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
338 { this->
gridVariables(domainId).update(curSol[domainId]); });
346 using namespace Dune::Hybrid;
347 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
348 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
352 template<std::
size_t i>
353 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
354 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
357 template<std::
size_t i>
358 const auto&
problem(Dune::index_constant<i> domainId)
const
359 {
return *std::get<domainId>(problemTuple_); }
362 template<std::
size_t i>
363 [[deprecated(
"Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
368 template<std::
size_t i>
370 {
return *std::get<domainId>(gridGeometryTuple_); }
373 template<std::
size_t i>
374 const auto&
gridView(Dune::index_constant<i> domainId)
const
378 template<std::
size_t i>
380 {
return *std::get<domainId>(gridVariablesTuple_); }
383 template<std::
size_t i>
385 {
return *std::get<domainId>(gridVariablesTuple_); }
393 {
return *jacobian_; }
397 {
return *residual_; }
401 {
return *prevSol_; }
408 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
421 {
return isStationaryProblem_; }
426 template<std::
size_t i>
428 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
436 void resetResidual_()
440 residual_ = std::make_shared<SolutionVector>();
448 void resetJacobian_()
452 jacobian_ = std::make_shared<JacobianMatrix>();
461 void checkAssemblerState_()
const
463 if (!isStationaryProblem_ && !prevSol_)
464 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
466 if (isStationaryProblem_ && prevSol_)
467 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
468 <<
" Did you forget to set the timeLoop to make this problem instationary?");
471 template<std::
size_t i,
class JacRow,
class SubRes>
472 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
475 assemble_(domainId, [&](
const auto& element)
477 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
478 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
482 template<std::
size_t i,
class SubRes>
483 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
486 assemble_(domainId, [&](
const auto& element)
488 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
489 subDomainAssembler.assembleResidual(subRes);
499 template<std::
size_t i,
class AssembleElementFunc>
500 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
503 for (
const auto& element : elements(
gridView(domainId)))
504 assembleElement(element);
508 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
509 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
510 Dune::index_constant<j> domainJ)
const
519 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
520 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
521 Dune::index_constant<j> domainJ)
const
529 ProblemTuple problemTuple_;
532 GridGeometryTuple gridGeometryTuple_;
535 GridVariablesTuple gridVariablesTuple_;
538 std::shared_ptr<const TimeLoop> timeLoop_;
544 bool isStationaryProblem_;
547 std::shared_ptr<JacobianMatrix> jacobian_;
548 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.
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
@ box
Definition method.hh:38
@ ccmpfa
Definition method.hh:38
@ staggered
Definition method.hh:38
@ cctpfa
Definition method.hh:38
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry &gridGeometry)
Helper function to generate Jacobian pattern for the box method.
Definition jacobianpattern.hh:39
Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager &couplingManager, Dune::index_constant< i > domainI, const GridGeometryI &gridGeometryI, Dune::index_constant< j > domainJ, const GridGeometryJ &gridGeometryJ)
Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for cell-centered schemes...
Definition couplingjacobianpattern.hh:42
make the local view function available whenever we use the grid geometry
Definition adapt.hh:29
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
const Scalar PengRobinsonMixture< Scalar, StaticParameters >::u
Definition pengrobinsonmixture.hh:167
Definition common/properties.hh:91
Manages the handling of time dependent problems.
Definition timeloop.hh:65
The default time loop for instationary simulations.
Definition timeloop.hh:110
Definition multidomain/couplingmanager.hh:46
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition multidomain/fvassembler.hh:59
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition multidomain/fvassembler.hh:353
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition multidomain/fvassembler.hh:334
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::SolutionVector SolutionVector
Definition multidomain/fvassembler.hh:85
SolutionVector ResidualType
Definition multidomain/fvassembler.hh:86
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::template SubDomain< id >::Problem Problem
Definition multidomain/fvassembler.hh:82
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition multidomain/fvassembler.hh:392
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:227
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::JacobianMatrix JacobianMatrix
Definition multidomain/fvassembler.hh:84
StaggeredCouplingManager< StaggeredMultiDomainTraits< TypeTag, TypeTag > > CouplingManager
Definition multidomain/fvassembler.hh:88
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition multidomain/fvassembler.hh:220
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition multidomain/fvassembler.hh:344
Scalar residualNorm(const SolutionVector &curSol)
Definition multidomain/fvassembler.hh:244
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition multidomain/fvassembler.hh:384
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:414
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::template SubDomain< id >::GridVariables GridVariables
Definition multidomain/fvassembler.hh:73
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:261
const auto & fvGridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition multidomain/fvassembler.hh:364
const auto & gridGeometry(Dune::index_constant< i > domainId) const
Definition multidomain/fvassembler.hh:369
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition multidomain/fvassembler.hh:324
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition multidomain/fvassembler.hh:276
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition multidomain/fvassembler.hh:374
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 StaggeredMultiDomainTraits< TypeTag, TypeTag >::template SubDomain< id >::GridGeometry GridGeometry
Definition multidomain/fvassembler.hh:76
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
Definition multidomain/fvassembler.hh:388
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition multidomain/fvassembler.hh:407
const SolutionVector & prevSol() const
Definition multidomain/fvassembler.hh:400
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
Definition multidomain/fvassembler.hh:379
std::shared_ptr< CouplingManager > couplingManager_
Definition multidomain/fvassembler.hh:432
StaggeredMultiDomainTraits< TypeTag, TypeTag > Traits
Definition multidomain/fvassembler.hh:64
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition multidomain/fvassembler.hh:289
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
Definition multidomain/fvassembler.hh:70
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:427
const auto & problem(Dune::index_constant< i > domainId) const
Definition multidomain/fvassembler.hh:358
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition multidomain/fvassembler.hh:204
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition multidomain/fvassembler.hh:420
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::Scalar Scalar
Definition multidomain/fvassembler.hh:66
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition multidomain/fvassembler.hh:308
SolutionVector & residual()
Definition multidomain/fvassembler.hh:396
The cell-centered scheme multidomain local assembler.
Definition subdomaincclocalassembler.hh:268
Declares all properties used in Dumux.