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;
55 static constexpr auto stokesFaceIdx =
typename MDTraits::template SubDomain<1>::Index();
56 static constexpr auto cellCenterIdx =
typename MDTraits::template SubDomain<0>::Index();
57 static constexpr auto faceIdx =
typename MDTraits::template SubDomain<1>::Index();
59 static constexpr auto darcyIdx =
typename MDTraits::template SubDomain<2>::Index();
63 using SolutionVector =
typename MDTraits::SolutionVector;
66 using StokesTypeTag =
typename MDTraits::template SubDomain<0>::TypeTag;
67 using DarcyTypeTag =
typename MDTraits::template SubDomain<2>::TypeTag;
69 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
70 using CouplingStencil = CouplingStencils::mapped_type;
73 template<std::
size_t id>
74 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
85 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
89 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
90 template<std::
size_t id>
using PrimaryVariables =
typename MDTraits::template SubDomain<id>::PrimaryVariables;
91 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
95 using VelocityVector =
typename Element<stokesIdx>::Geometry::GlobalCoordinate;
99 struct StationaryStokesCouplingContext
101 Element<darcyIdx> element;
102 FVElementGeometry<darcyIdx> fvGeometry;
103 std::size_t darcyScvfIdx;
104 std::size_t stokesScvfIdx;
105 VolumeVariables<darcyIdx> volVars;
108 struct StationaryDarcyCouplingContext
110 Element<stokesIdx> element;
111 FVElementGeometry<stokesIdx> fvGeometry;
112 std::size_t stokesScvfIdx;
113 std::size_t darcyScvfIdx;
114 VelocityVector velocity;
115 VolumeVariables<stokesIdx> volVars;
125 std::shared_ptr<
const GridGeometry<darcyIdx>> darcyFvGridGeometry) : couplingMapper_(*this)
134 void init(std::shared_ptr<
const Problem<stokesIdx>> stokesProblem,
135 std::shared_ptr<
const Problem<darcyIdx>> darcyProblem,
136 const SolutionVector&
curSol)
138 if (Dune::FloatCmp::ne(stokesProblem->gravity(), darcyProblem->spatialParams().gravity({})))
139 DUNE_THROW(Dune::InvalidStateException,
"Both models must use the same gravity vector");
141 this->
setSubProblems(std::make_tuple(stokesProblem, stokesProblem, darcyProblem));
143 couplingData_ = std::make_shared<CouplingData>(*
this);
161 darcyToStokesFaceCouplingStencils_,
162 stokesCellCenterCouplingStencils_,
163 stokesFaceCouplingStencils_);
165 for(
auto&& stencil : darcyToStokesCellCenterCouplingStencils_)
167 for(
auto&& stencil : darcyToStokesFaceCouplingStencils_)
169 for(
auto&& stencil : stokesCellCenterCouplingStencils_)
171 for(
auto&& stencil : stokesFaceCouplingStencils_)
185 template<std::
size_t i,
class Assembler, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
186 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<stokesCellCenterIdx>& element,
const Assembler& assembler)
const
192 template<std::
size_t i, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
193 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<stokesCellCenterIdx>& element)
const
195 stokesCouplingContext_.clear();
197 const auto stokesElementIdx = this->
problem(
stokesIdx).gridGeometry().elementMapper().index(element);
198 boundStokesElemIdx_ = stokesElementIdx;
208 for(
auto&& indices : darcyIndices)
210 const auto& darcyElement = this->
problem(
darcyIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
211 darcyFvGeometry.bindElement(darcyElement);
212 const auto& scv = (*scvs(darcyFvGeometry).begin());
215 VolumeVariables<darcyIdx> darcyVolVars;
216 darcyVolVars.update(darcyElemSol, this->
problem(
darcyIdx), darcyElement, scv);
219 stokesCouplingContext_.push_back({darcyElement, darcyFvGeometry, indices.scvfIdx, indices.flipScvfIdx, darcyVolVars});
226 template<
class Assembler>
227 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI,
const Element<darcyIdx>& element,
const Assembler& assembler)
const
235 darcyCouplingContext_.clear();
237 const auto darcyElementIdx = this->
problem(
darcyIdx).gridGeometry().elementMapper().index(element);
238 boundDarcyElemIdx_ = darcyElementIdx;
248 for(
auto&& indices : stokesElementIndices)
250 const auto& stokesElement = this->
problem(
stokesIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
251 stokesFvGeometry.bindElement(stokesElement);
253 VelocityVector faceVelocity(0.0);
255 for(
const auto& scvf : scvfs(stokesFvGeometry))
257 if(scvf.index() == indices.scvfIdx)
261 using PriVarsType =
typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
263 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(cellCenterPriVars);
265 VolumeVariables<stokesIdx> stokesVolVars;
266 for(
const auto& scv : scvs(stokesFvGeometry))
267 stokesVolVars.update(elemSol, this->
problem(
stokesIdx), stokesElement, scv);
270 darcyCouplingContext_.push_back({stokesElement, stokesFvGeometry, indices.scvfIdx, indices.flipScvfIdx, faceVelocity, stokesVolVars});
277 template<
class LocalAssemblerI>
279 const LocalAssemblerI& localAssemblerI,
280 Dune::index_constant<darcyIdx> domainJ,
281 std::size_t dofIdxGlobalJ,
282 const PrimaryVariables<darcyIdx>& priVarsJ,
285 this->
curSol()[domainJ][dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
291 template<
class LocalAssemblerI>
293 const LocalAssemblerI& localAssemblerI,
294 Dune::index_constant<stokesCellCenterIdx> domainJ,
295 const std::size_t dofIdxGlobalJ,
296 const PrimaryVariables<stokesCellCenterIdx>& priVars,
299 this->
curSol()[domainJ][dofIdxGlobalJ] = priVars;
301 for (
auto& data : darcyCouplingContext_)
303 const auto stokesElemIdx = this->
problem(
stokesIdx).gridGeometry().elementMapper().index(data.element);
305 if(stokesElemIdx != dofIdxGlobalJ)
308 using PriVarsType =
typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
309 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(priVars);
311 for(
const auto& scv : scvs(data.fvGeometry))
319 template<
class LocalAssemblerI>
321 const LocalAssemblerI& localAssemblerI,
322 Dune::index_constant<stokesFaceIdx> domainJ,
323 const std::size_t dofIdxGlobalJ,
324 const PrimaryVariables<1>& priVars,
327 this->
curSol()[domainJ][dofIdxGlobalJ] = priVars;
329 for (
auto& data : darcyCouplingContext_)
331 for(
const auto& scvf : scvfs(data.fvGeometry))
333 if(scvf.dofIndex() == dofIdxGlobalJ)
334 data.velocity[scvf.directionIndex()] = priVars;
342 template<std::
size_t i,
class LocalAssemblerI, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
344 const LocalAssemblerI& localAssemblerI,
345 Dune::index_constant<darcyIdx> domainJ,
346 const std::size_t dofIdxGlobalJ,
347 const PrimaryVariables<darcyIdx>& priVars,
350 this->
curSol()[domainJ][dofIdxGlobalJ] = priVars;
352 for (
auto& data : stokesCouplingContext_)
354 const auto darcyElemIdx = this->
problem(
darcyIdx).gridGeometry().elementMapper().index(data.element);
356 if(darcyElemIdx != dofIdxGlobalJ)
361 for(
const auto& scv : scvs(data.fvGeometry))
362 data.volVars.update(darcyElemSol, this->
problem(
darcyIdx), data.element, scv);
373 return *couplingData_;
379 const auto&
stokesCouplingContext(
const Element<stokesIdx>& element,
const SubControlVolumeFace<stokesIdx>& scvf)
const
381 if (stokesCouplingContext_.empty() || boundStokesElemIdx_ != scvf.insideScvIdx())
384 for(
const auto& context : stokesCouplingContext_)
386 if(scvf.index() == context.stokesScvfIdx)
390 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
396 const auto&
darcyCouplingContext(
const Element<darcyIdx>& element,
const SubControlVolumeFace<darcyIdx>& scvf)
const
398 if (darcyCouplingContext_.empty() || boundDarcyElemIdx_ != scvf.insideScvIdx())
401 for(
const auto& context : darcyCouplingContext_)
403 if(scvf.index() == context.darcyScvfIdx)
407 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
418 const CouplingStencil&
couplingStencil(Dune::index_constant<stokesCellCenterIdx> domainI,
419 const Element<stokesIdx>& element,
420 Dune::index_constant<darcyIdx> domainJ)
const
422 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
423 if(stokesCellCenterCouplingStencils_.count(eIdx))
424 return stokesCellCenterCouplingStencils_.at(eIdx);
426 return emptyStencil_;
433 template<std::
size_t i, std::
size_t j>
435 const Element<i>& element,
436 Dune::index_constant<j> domainJ)
const
437 {
return emptyStencil_; }
444 const Element<darcyIdx>& element,
445 Dune::index_constant<stokesCellCenterIdx> domainJ)
const
447 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
448 if(darcyToStokesCellCenterCouplingStencils_.count(eIdx))
449 return darcyToStokesCellCenterCouplingStencils_.at(eIdx);
451 return emptyStencil_;
459 const Element<darcyIdx>& element,
460 Dune::index_constant<stokesFaceIdx> domainJ)
const
462 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
463 if (darcyToStokesFaceCouplingStencils_.count(eIdx))
464 return darcyToStokesFaceCouplingStencils_.at(eIdx);
466 return emptyStencil_;
473 template<std::
size_t i, std::
size_t j>
475 const SubControlVolumeFace<stokesIdx>& scvf,
476 Dune::index_constant<j> domainJ)
const
477 {
return emptyStencil_; }
483 const SubControlVolumeFace<stokesIdx>& scvf,
484 Dune::index_constant<darcyIdx> domainJ)
const
486 const auto faceDofIdx = scvf.dofIndex();
487 if(stokesFaceCouplingStencils_.count(faceDofIdx))
488 return stokesFaceCouplingStencils_.at(faceDofIdx);
490 return emptyStencil_;
498 template<
class IdType>
500 {
return emptyStencil_; }
505 template<
class IdType>
507 {
return emptyStencil_; }
512 bool isCoupledEntity(Dune::index_constant<stokesIdx>,
const SubControlVolumeFace<stokesFaceIdx>& scvf)
const
514 return stokesFaceCouplingStencils_.count(scvf.dofIndex());
520 bool isCoupledEntity(Dune::index_constant<darcyIdx>,
const SubControlVolumeFace<darcyIdx>& scvf)
const
529 {
return emptyStencil_; }
533 std::sort(stencil.begin(), stencil.end());
534 stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
539 std::vector<bool> isCoupledDarcyDof_;
540 std::shared_ptr<CouplingData> couplingData_;
542 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesCellCenterCouplingStencils_;
543 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesFaceCouplingStencils_;
544 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesCellCenterCouplingStencils_;
545 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesFaceCouplingStencils_;
546 std::vector<std::size_t> emptyStencil_;
551 mutable std::vector<StationaryStokesCouplingContext> stokesCouplingContext_;
552 mutable std::vector<StationaryDarcyCouplingContext> darcyCouplingContext_;
554 mutable std::size_t boundStokesElemIdx_;
555 mutable std::size_t boundDarcyElemIdx_;
557 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
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
A vector of size number equations that can be used for Neumann fluxes, sources, residuals,...
Definition: common/properties.hh:61
The type of the grid view according to the grid type.
Definition: common/properties.hh:63
Traits class encapsulating model specifications.
Definition: common/properties.hh:65
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Stores the boundary types on an element.
Definition: common/properties.hh:112
Definition: common/properties.hh:150
The type for a global container for the volume variables.
Definition: common/properties.hh:176
The global vector of flux variable containers.
Definition: common/properties.hh:186
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:190
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:158
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:193
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:434
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:278
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:482
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:528
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:418
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:343
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:499
static constexpr auto darcyIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:59
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:379
StokesDarcyCouplingManager(std::shared_ptr< const GridGeometry< stokesIdx > > stokesFvGridGeometry, std::shared_ptr< const GridGeometry< darcyIdx > > darcyFvGridGeometry)
Constructor.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:124
void update()
Update after the grid has changed.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:148
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:233
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:396
void removeDuplicates_(std::vector< std::size_t > &stencil)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:531
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:134
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:443
const auto & couplingData() const
Access the coupling data.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:371
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesFaceIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< 1 > &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:320
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:292
static constexpr auto faceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:57
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:458
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:512
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:227
static constexpr auto stokesCellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:54
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:506
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:186
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:474
static constexpr auto cellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:56
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:520
void updateSolution(const SolutionVector &curSol)
Update the solution vector before assembly.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:154
static constexpr auto stokesFaceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:55
static constexpr auto stokesIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:58
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:49
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:167
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:175
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:183
void computeCouplingMapsAndStencils(Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:94
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.