25#ifndef DUMUX_STOKES_DARCY_COUPLINGMANAGER_HH
26#define DUMUX_STOKES_DARCY_COUPLINGMANAGER_HH
31#include <dune/common/float_cmp.hh>
32#include <dune/common/exceptions.hh>
46template<
class MDTraits>
50 using Scalar =
typename MDTraits::Scalar;
54 static constexpr auto stokesFaceIdx =
typename MDTraits::template SubDomain<0>::Index();
57 static constexpr auto darcyIdx =
typename MDTraits::template SubDomain<2>::Index();
61 using SolutionVector =
typename MDTraits::SolutionVector;
64 using StokesTypeTag =
typename MDTraits::template SubDomain<0>::TypeTag;
65 using DarcyTypeTag =
typename MDTraits::template SubDomain<2>::TypeTag;
67 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
68 using CouplingStencil = CouplingStencils::mapped_type;
71 template<std::
size_t id>
72 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
83 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
87 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
88 template<std::
size_t id>
using PrimaryVariables =
typename MDTraits::template SubDomain<id>::PrimaryVariables;
89 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
93 using VelocityVector =
typename Element<stokesIdx>::Geometry::GlobalCoordinate;
97 struct StationaryStokesCouplingContext
99 Element<darcyIdx> element;
100 FVElementGeometry<darcyIdx> fvGeometry;
101 std::size_t darcyScvfIdx;
102 std::size_t stokesScvfIdx;
103 VolumeVariables<darcyIdx> volVars;
106 struct StationaryDarcyCouplingContext
108 Element<stokesIdx> element;
109 FVElementGeometry<stokesIdx> fvGeometry;
110 std::size_t stokesScvfIdx;
111 std::size_t darcyScvfIdx;
112 VelocityVector velocity;
113 VolumeVariables<stokesIdx> volVars;
123 std::shared_ptr<
const GridGeometry<darcyIdx>> darcyFvGridGeometry)
132 void init(std::shared_ptr<
const Problem<stokesIdx>> stokesProblem,
133 std::shared_ptr<
const Problem<darcyIdx>> darcyProblem,
134 const SolutionVector&
curSol)
136 if (Dune::FloatCmp::ne(stokesProblem->gravity(), darcyProblem->spatialParams().gravity({})))
137 DUNE_THROW(Dune::InvalidStateException,
"Both models must use the same gravity vector");
139 this->
setSubProblems(std::make_tuple(stokesProblem, stokesProblem, darcyProblem));
141 couplingData_ = std::make_shared<CouplingData>(*
this);
159 darcyToStokesCellCenterCouplingStencils_,
160 darcyToStokesFaceCouplingStencils_,
161 stokesCellCenterCouplingStencils_,
162 stokesFaceCouplingStencils_);
164 for(
auto&& stencil : darcyToStokesCellCenterCouplingStencils_)
166 for(
auto&& stencil : darcyToStokesFaceCouplingStencils_)
168 for(
auto&& stencil : stokesCellCenterCouplingStencils_)
170 for(
auto&& stencil : stokesFaceCouplingStencils_)
184 template<std::
size_t i,
class Assembler, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
185 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<stokesCellCenterIdx>& element,
const Assembler& assembler)
const
191 template<std::
size_t i, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
192 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<stokesCellCenterIdx>& element)
const
194 stokesCouplingContext_.clear();
196 const auto stokesElementIdx = this->
problem(
stokesIdx).gridGeometry().elementMapper().index(element);
197 boundStokesElemIdx_ = stokesElementIdx;
207 for(
auto&& indices : darcyIndices)
209 const auto& darcyElement = this->
problem(
darcyIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
210 darcyFvGeometry.bindElement(darcyElement);
211 const auto& scv = (*scvs(darcyFvGeometry).begin());
214 VolumeVariables<darcyIdx> darcyVolVars;
215 darcyVolVars.update(darcyElemSol, this->
problem(
darcyIdx), darcyElement, scv);
218 stokesCouplingContext_.push_back({darcyElement, darcyFvGeometry, indices.scvfIdx, indices.flipScvfIdx, darcyVolVars});
225 template<
class Assembler>
226 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI,
const Element<darcyIdx>& element,
const Assembler& assembler)
const
234 darcyCouplingContext_.clear();
236 const auto darcyElementIdx = this->
problem(
darcyIdx).gridGeometry().elementMapper().index(element);
237 boundDarcyElemIdx_ = darcyElementIdx;
247 for(
auto&& indices : stokesElementIndices)
249 const auto& stokesElement = this->
problem(
stokesIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
250 stokesFvGeometry.bindElement(stokesElement);
252 VelocityVector faceVelocity(0.0);
254 for(
const auto& scvf : scvfs(stokesFvGeometry))
256 if(scvf.index() == indices.scvfIdx)
260 using PriVarsType =
typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
262 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(cellCenterPriVars);
264 VolumeVariables<stokesIdx> stokesVolVars;
265 for(
const auto& scv : scvs(stokesFvGeometry))
266 stokesVolVars.update(elemSol, this->
problem(
stokesIdx), stokesElement, scv);
269 darcyCouplingContext_.push_back({stokesElement, stokesFvGeometry, indices.scvfIdx, indices.flipScvfIdx, faceVelocity, stokesVolVars});
276 template<
class LocalAssemblerI>
278 const LocalAssemblerI& localAssemblerI,
279 Dune::index_constant<darcyIdx> domainJ,
280 std::size_t dofIdxGlobalJ,
281 const PrimaryVariables<darcyIdx>& priVarsJ,
284 this->
curSol()[domainJ][dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
290 template<
class LocalAssemblerI>
292 const LocalAssemblerI& localAssemblerI,
293 Dune::index_constant<stokesCellCenterIdx> domainJ,
294 const std::size_t dofIdxGlobalJ,
295 const PrimaryVariables<stokesCellCenterIdx>& priVars,
298 this->
curSol()[domainJ][dofIdxGlobalJ] = priVars;
300 for (
auto& data : darcyCouplingContext_)
302 const auto stokesElemIdx = this->
problem(
stokesIdx).gridGeometry().elementMapper().index(data.element);
304 if(stokesElemIdx != dofIdxGlobalJ)
307 using PriVarsType =
typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
308 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(priVars);
310 for(
const auto& scv : scvs(data.fvGeometry))
318 template<
class LocalAssemblerI>
320 const LocalAssemblerI& localAssemblerI,
321 Dune::index_constant<stokesFaceIdx> domainJ,
322 const std::size_t dofIdxGlobalJ,
323 const PrimaryVariables<stokesFaceIdx>& priVars,
326 this->
curSol()[domainJ][dofIdxGlobalJ] = priVars;
328 for (
auto& data : darcyCouplingContext_)
330 for(
const auto& scvf : scvfs(data.fvGeometry))
332 if(scvf.dofIndex() == dofIdxGlobalJ)
333 data.velocity[scvf.directionIndex()] = priVars;
341 template<std::
size_t i,
class LocalAssemblerI, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
343 const LocalAssemblerI& localAssemblerI,
344 Dune::index_constant<darcyIdx> domainJ,
345 const std::size_t dofIdxGlobalJ,
346 const PrimaryVariables<darcyIdx>& priVars,
349 this->
curSol()[domainJ][dofIdxGlobalJ] = priVars;
351 for (
auto& data : stokesCouplingContext_)
353 const auto darcyElemIdx = this->
problem(
darcyIdx).gridGeometry().elementMapper().index(data.element);
355 if(darcyElemIdx != dofIdxGlobalJ)
360 for(
const auto& scv : scvs(data.fvGeometry))
361 data.volVars.update(darcyElemSol, this->
problem(
darcyIdx), data.element, scv);
372 return *couplingData_;
378 const auto&
stokesCouplingContext(
const Element<stokesIdx>& element,
const SubControlVolumeFace<stokesIdx>& scvf)
const
380 if (stokesCouplingContext_.empty() || boundStokesElemIdx_ != scvf.insideScvIdx())
383 for(
const auto& context : stokesCouplingContext_)
385 if(scvf.index() == context.stokesScvfIdx)
389 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
395 const auto&
darcyCouplingContext(
const Element<darcyIdx>& element,
const SubControlVolumeFace<darcyIdx>& scvf)
const
397 if (darcyCouplingContext_.empty() || boundDarcyElemIdx_ != scvf.insideScvIdx())
400 for(
const auto& context : darcyCouplingContext_)
402 if(scvf.index() == context.darcyScvfIdx)
406 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
417 const CouplingStencil&
couplingStencil(Dune::index_constant<stokesCellCenterIdx> domainI,
418 const Element<stokesIdx>& element,
419 Dune::index_constant<darcyIdx> domainJ)
const
421 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
422 if(stokesCellCenterCouplingStencils_.count(eIdx))
423 return stokesCellCenterCouplingStencils_.at(eIdx);
425 return emptyStencil_;
432 template<std::
size_t i, std::
size_t j>
434 const Element<i>& element,
435 Dune::index_constant<j> domainJ)
const
436 {
return emptyStencil_; }
443 const Element<darcyIdx>& element,
444 Dune::index_constant<stokesCellCenterIdx> domainJ)
const
446 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
447 if(darcyToStokesCellCenterCouplingStencils_.count(eIdx))
448 return darcyToStokesCellCenterCouplingStencils_.at(eIdx);
450 return emptyStencil_;
458 const Element<darcyIdx>& element,
459 Dune::index_constant<stokesFaceIdx> domainJ)
const
461 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
462 if (darcyToStokesFaceCouplingStencils_.count(eIdx))
463 return darcyToStokesFaceCouplingStencils_.at(eIdx);
465 return emptyStencil_;
472 template<std::
size_t i, std::
size_t j>
474 const SubControlVolumeFace<stokesIdx>& scvf,
475 Dune::index_constant<j> domainJ)
const
476 {
return emptyStencil_; }
482 const SubControlVolumeFace<stokesIdx>& scvf,
483 Dune::index_constant<darcyIdx> domainJ)
const
485 const auto faceDofIdx = scvf.dofIndex();
486 if(stokesFaceCouplingStencils_.count(faceDofIdx))
487 return stokesFaceCouplingStencils_.at(faceDofIdx);
489 return emptyStencil_;
497 template<
class IdType>
499 {
return emptyStencil_; }
504 template<
class IdType>
506 {
return emptyStencil_; }
511 bool isCoupledEntity(Dune::index_constant<stokesIdx>,
const SubControlVolumeFace<stokesFaceIdx>& scvf)
const
513 return stokesFaceCouplingStencils_.count(scvf.dofIndex());
519 bool isCoupledEntity(Dune::index_constant<darcyIdx>,
const SubControlVolumeFace<darcyIdx>& scvf)
const
528 {
return emptyStencil_; }
532 std::sort(stencil.begin(), stencil.end());
533 stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
538 std::vector<bool> isCoupledDarcyDof_;
539 std::shared_ptr<CouplingData> couplingData_;
541 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesCellCenterCouplingStencils_;
542 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesFaceCouplingStencils_;
543 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesCellCenterCouplingStencils_;
544 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesFaceCouplingStencils_;
545 std::vector<std::size_t> emptyStencil_;
550 mutable std::vector<StationaryStokesCouplingContext> stokesCouplingContext_;
551 mutable std::vector<StationaryDarcyCouplingContext> darcyCouplingContext_;
553 mutable std::size_t boundStokesElemIdx_;
554 mutable std::size_t boundDarcyElemIdx_;
556 CouplingMapper couplingMapper_;
The interface of the coupling manager for multi domain problems.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
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
A vector of size number equations that can be used for Neumann fluxes, sources, residuals,...
Definition: common/properties.hh:51
Traits class encapsulating model specifications.
Definition: common/properties.hh:53
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:57
Stores the boundary types on an element.
Definition: common/properties.hh:98
Definition: common/properties.hh:101
The type for a global container for the volume variables.
Definition: common/properties.hh:108
The global vector of flux variable containers.
Definition: common/properties.hh:118
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:122
Definition: couplingdata.hh:206
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:49
void computeStencils()
Prepare the coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:156
void bindCouplingContext(Dune::index_constant< i > domainI, const Element< stokesCellCenterIdx > &element) const
prepares all data and variables that are necessary to evaluate the residual of an Darcy element (i....
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:192
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &element, Dune::index_constant< j > domainJ) const
The coupling stencil of domain I, i.e. which domain J DOFs the given domain I element's residual depe...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:433
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< darcyIdx > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< darcyIdx > &priVarsJ, int pvIdxJ)
Update the coupling context for the Darcy residual w.r.t. Darcy DOFs.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:277
const CouplingStencil & couplingStencil(Dune::index_constant< stokesFaceIdx > domainI, const SubControlVolumeFace< stokesIdx > &scvf, Dune::index_constant< darcyIdx > domainJ) const
The coupling stencil of a Stokes face w.r.t. Darcy DOFs.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:481
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:527
const CouplingStencil & couplingStencil(Dune::index_constant< stokesCellCenterIdx > domainI, const Element< stokesIdx > &element, Dune::index_constant< darcyIdx > domainJ) const
The coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:417
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< darcyIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< darcyIdx > &priVars, int pvIdxJ)
Update the coupling context for the Stokes cc residual w.r.t. the Darcy DOFs (FaceToDarcy)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:342
const std::vector< std::size_t > & getAdditionalDofDependencies(IdType id, std::size_t stokesElementIdx) const
There are no additional dof dependencies.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:498
static constexpr auto darcyIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:57
const auto & stokesCouplingContext(const Element< stokesIdx > &element, const SubControlVolumeFace< stokesIdx > &scvf) const
Access the coupling context needed for the Stokes domain.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:378
StokesDarcyCouplingManager(std::shared_ptr< const GridGeometry< stokesIdx > > stokesFvGridGeometry, std::shared_ptr< const GridGeometry< darcyIdx > > darcyFvGridGeometry)
Constructor.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:122
void update()
Update after the grid has changed.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:146
void bindCouplingContext(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &element) const
prepares all data and variables that are necessary to evaluate the residual of an Darcy element (i....
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:232
const auto & darcyCouplingContext(const Element< darcyIdx > &element, const SubControlVolumeFace< darcyIdx > &scvf) const
Access the coupling context needed for the Darcy domain.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:395
void removeDuplicates_(std::vector< std::size_t > &stencil)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:530
void init(std::shared_ptr< const Problem< stokesIdx > > stokesProblem, std::shared_ptr< const Problem< darcyIdx > > darcyProblem, const SolutionVector &curSol)
Methods to be accessed by main.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:132
const CouplingStencil & couplingStencil(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &element, Dune::index_constant< stokesCellCenterIdx > domainJ) const
The coupling stencil of domain I, i.e. which domain J dofs the given domain I element's residual depe...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:442
const auto & couplingData() const
Access the coupling data.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:370
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesCellCenterIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< stokesCellCenterIdx > &priVars, int pvIdxJ)
Update the coupling context for the Darcy residual w.r.t. the Stokes cell-center DOFs (DarcyToCC)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:291
const CouplingStencil & couplingStencil(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &element, Dune::index_constant< stokesFaceIdx > domainJ) const
The coupling stencil of domain I, i.e. which domain J dofs the given domain I element's residual depe...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:457
bool isCoupledEntity(Dune::index_constant< stokesIdx >, const SubControlVolumeFace< stokesFaceIdx > &scvf) const
Returns whether a given free flow scvf is coupled to the other domain.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:511
void bindCouplingContext(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &element, const Assembler &assembler) const
prepares all data and variables that are necessary to evaluate the residual of an Darcy element (i....
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:226
static constexpr auto stokesCellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:55
const std::vector< std::size_t > & getAdditionalDofDependenciesInverse(IdType id, std::size_t darcyElementIdx) const
There are no additional dof dependencies.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:505
void bindCouplingContext(Dune::index_constant< i > domainI, const Element< stokesCellCenterIdx > &element, const Assembler &assembler) const
prepares all data and variables that are necessary to evaluate the residual of an Darcy element (i....
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:185
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesFaceIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< stokesFaceIdx > &priVars, int pvIdxJ)
Update the coupling context for the Darcy residual w.r.t. the Stokes face DOFs (DarcyToFace)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:319
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const SubControlVolumeFace< stokesIdx > &scvf, Dune::index_constant< j > domainJ) const
The coupling stencil of domain I, i.e. which domain J DOFs the given domain I element's residual depe...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:473
bool isCoupledEntity(Dune::index_constant< darcyIdx >, const SubControlVolumeFace< darcyIdx > &scvf) const
Returns whether a given free flow scvf is coupled to the other domain.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:519
void updateSolution(const SolutionVector &curSol)
Update the solution vector before assembly.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:152
static constexpr auto stokesFaceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:54
static constexpr auto stokesIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:56
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:42
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:133
void computeCouplingMapsAndStencils(const CouplingManager &couplingManager, Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:56
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:141
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:149
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:247
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264
SolutionVector & curSol()
the solution vector of the coupled problem
Definition: multidomain/couplingmanager.hh:278
Base coupling manager for the staggered discretization.
Definition: staggeredcouplingmanager.hh:42
decltype(auto) evalCouplingResidual(Dune::index_constant< cellCenterIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluates the element residual of a coupled element of domain i which depends on the variables at the...
Definition: staggeredcouplingmanager.hh:149
const CouplingStencil & couplingStencil(Dune::index_constant< cellCenterIdx > domainI, const Element &elementI, Dune::index_constant< faceIdx > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: staggeredcouplingmanager.hh:93
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, const std::size_t dofIdxGlobalJ, const PriVarsJ &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: staggeredcouplingmanager.hh:68
Declares all properties used in Dumux.
The local element solution class for staggered methods.