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>
47template<
class MDTraits>
51 using Scalar =
typename MDTraits::Scalar;
55 static constexpr auto stokesFaceIdx =
typename MDTraits::template SubDomain<0>::Index();
58 static constexpr auto darcyIdx =
typename MDTraits::template SubDomain<2>::Index();
62 using SolutionVector =
typename MDTraits::SolutionVector;
65 using StokesTypeTag =
typename MDTraits::template SubDomain<0>::TypeTag;
66 using DarcyTypeTag =
typename MDTraits::template SubDomain<2>::TypeTag;
68 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
69 using CouplingStencil = CouplingStencils::mapped_type;
72 template<std::
size_t id>
73 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
79 template<std::
size_t id>
using PrimaryVariables =
typename MDTraits::template SubDomain<id>::PrimaryVariables;
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 SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
94 using VelocityVector =
typename Element<stokesIdx>::Geometry::GlobalCoordinate;
98 struct StationaryStokesCouplingContext
100 Element<darcyIdx> element;
101 FVElementGeometry<darcyIdx> fvGeometry;
102 std::size_t darcyScvfIdx;
103 std::size_t stokesScvfIdx;
104 VolumeVariables<darcyIdx> volVars;
107 struct StationaryDarcyCouplingContext
109 Element<stokesIdx> element;
110 FVElementGeometry<stokesIdx> fvGeometry;
111 std::size_t stokesScvfIdx;
112 std::size_t darcyScvfIdx;
113 VelocityVector velocity;
114 VolumeVariables<stokesIdx> volVars;
124 std::shared_ptr<
const GridGeometry<darcyIdx>> darcyFvGridGeometry)
133 void init(std::shared_ptr<
const Problem<stokesIdx>> stokesProblem,
134 std::shared_ptr<
const Problem<darcyIdx>> darcyProblem,
135 const SolutionVector&
curSol)
137 if (Dune::FloatCmp::ne(stokesProblem->gravity(), darcyProblem->spatialParams().gravity({})))
138 DUNE_THROW(Dune::InvalidStateException,
"Both models must use the same gravity vector");
140 this->
setSubProblems(std::make_tuple(stokesProblem, stokesProblem, darcyProblem));
142 couplingData_ = std::make_shared<CouplingData>(*
this);
152 darcyToStokesCellCenterCouplingStencils_,
153 darcyToStokesFaceCouplingStencils_,
154 stokesCellCenterCouplingStencils_,
155 stokesFaceCouplingStencils_);
157 for(
auto&& stencil : darcyToStokesCellCenterCouplingStencils_)
159 for(
auto&& stencil : darcyToStokesFaceCouplingStencils_)
161 for(
auto&& stencil : stokesCellCenterCouplingStencils_)
163 for(
auto&& stencil : stokesFaceCouplingStencils_)
177 template<std::
size_t i,
class Assembler, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
178 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<stokesCellCenterIdx>& element,
const Assembler& assembler)
const
184 template<std::
size_t i, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
185 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<stokesCellCenterIdx>& element)
const
187 stokesCouplingContext_.clear();
189 const auto stokesElementIdx = this->
problem(
stokesIdx).gridGeometry().elementMapper().index(element);
190 boundStokesElemIdx_ = stokesElementIdx;
199 for(
auto&& indices : darcyIndices)
201 const auto& darcyElement = this->
problem(
darcyIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
202 darcyFvGeometry.bindElement(darcyElement);
203 const auto& scv = (*scvs(darcyFvGeometry).begin());
206 VolumeVariables<darcyIdx> darcyVolVars;
207 darcyVolVars.update(darcyElemSol, this->
problem(
darcyIdx), darcyElement, scv);
210 stokesCouplingContext_.push_back({darcyElement, darcyFvGeometry, indices.scvfIdx, indices.flipScvfIdx, darcyVolVars});
217 template<
class Assembler>
218 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI,
const Element<darcyIdx>& element,
const Assembler& assembler)
const
226 darcyCouplingContext_.clear();
228 const auto darcyElementIdx = this->
problem(
darcyIdx).gridGeometry().elementMapper().index(element);
229 boundDarcyElemIdx_ = darcyElementIdx;
239 for(
auto&& indices : stokesElementIndices)
241 const auto& stokesElement = this->
problem(
stokesIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
242 stokesFvGeometry.bindElement(stokesElement);
244 VelocityVector faceVelocity(0.0);
246 for(
const auto& scvf : scvfs(stokesFvGeometry))
248 if(scvf.index() == indices.scvfIdx)
252 using PriVarsType =
typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
254 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(cellCenterPriVars);
256 VolumeVariables<stokesIdx> stokesVolVars;
257 for(
const auto& scv : scvs(stokesFvGeometry))
258 stokesVolVars.update(elemSol, this->
problem(
stokesIdx), stokesElement, scv);
261 darcyCouplingContext_.push_back({stokesElement, stokesFvGeometry, indices.scvfIdx, indices.flipScvfIdx, faceVelocity, stokesVolVars});
268 template<
class LocalAssemblerI>
270 const LocalAssemblerI& localAssemblerI,
271 Dune::index_constant<darcyIdx> domainJ,
272 std::size_t dofIdxGlobalJ,
273 const PrimaryVariables<darcyIdx>& priVarsJ,
276 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
282 template<
class LocalAssemblerI>
284 const LocalAssemblerI& localAssemblerI,
285 Dune::index_constant<stokesCellCenterIdx> domainJ,
286 const std::size_t dofIdxGlobalJ,
287 const PrimaryVariables<stokesCellCenterIdx>& priVars,
290 this->
curSol(domainJ)[dofIdxGlobalJ] = priVars;
292 for (
auto& data : darcyCouplingContext_)
294 const auto stokesElemIdx = this->
problem(
stokesIdx).gridGeometry().elementMapper().index(data.element);
296 if(stokesElemIdx != dofIdxGlobalJ)
299 using PriVarsType =
typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
300 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(priVars);
302 for(
const auto& scv : scvs(data.fvGeometry))
310 template<
class LocalAssemblerI>
312 const LocalAssemblerI& localAssemblerI,
313 Dune::index_constant<stokesFaceIdx> domainJ,
314 const std::size_t dofIdxGlobalJ,
315 const PrimaryVariables<stokesFaceIdx>& priVars,
318 this->
curSol(domainJ)[dofIdxGlobalJ] = priVars;
320 for (
auto& data : darcyCouplingContext_)
322 for(
const auto& scvf : scvfs(data.fvGeometry))
324 if(scvf.dofIndex() == dofIdxGlobalJ)
325 data.velocity[scvf.directionIndex()] = priVars;
333 template<std::
size_t i,
class LocalAssemblerI, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx),
int> = 0>
335 const LocalAssemblerI& localAssemblerI,
336 Dune::index_constant<darcyIdx> domainJ,
337 const std::size_t dofIdxGlobalJ,
338 const PrimaryVariables<darcyIdx>& priVars,
341 this->
curSol(domainJ)[dofIdxGlobalJ] = priVars;
343 for (
auto& data : stokesCouplingContext_)
345 const auto darcyElemIdx = this->
problem(
darcyIdx).gridGeometry().elementMapper().index(data.element);
347 if(darcyElemIdx != dofIdxGlobalJ)
352 for(
const auto& scv : scvs(data.fvGeometry))
353 data.volVars.update(darcyElemSol, this->
problem(
darcyIdx), data.element, scv);
364 return *couplingData_;
370 const auto&
stokesCouplingContext(
const Element<stokesIdx>& element,
const SubControlVolumeFace<stokesIdx>& scvf)
const
372 if (stokesCouplingContext_.empty() || boundStokesElemIdx_ != scvf.insideScvIdx())
375 for(
const auto& context : stokesCouplingContext_)
377 if(scvf.index() == context.stokesScvfIdx)
381 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
387 const auto&
darcyCouplingContext(
const Element<darcyIdx>& element,
const SubControlVolumeFace<darcyIdx>& scvf)
const
389 if (darcyCouplingContext_.empty() || boundDarcyElemIdx_ != scvf.insideScvIdx())
392 for(
const auto& context : darcyCouplingContext_)
394 if(scvf.index() == context.darcyScvfIdx)
398 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
409 const CouplingStencil&
couplingStencil(Dune::index_constant<stokesCellCenterIdx> domainI,
410 const Element<stokesIdx>& element,
411 Dune::index_constant<darcyIdx> domainJ)
const
413 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
414 if(stokesCellCenterCouplingStencils_.count(eIdx))
415 return stokesCellCenterCouplingStencils_.at(eIdx);
417 return emptyStencil_;
424 template<std::
size_t i, std::
size_t j>
426 const Element<i>& element,
427 Dune::index_constant<j> domainJ)
const
428 {
return emptyStencil_; }
435 const Element<darcyIdx>& element,
436 Dune::index_constant<stokesCellCenterIdx> domainJ)
const
438 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
439 if(darcyToStokesCellCenterCouplingStencils_.count(eIdx))
440 return darcyToStokesCellCenterCouplingStencils_.at(eIdx);
442 return emptyStencil_;
450 const Element<darcyIdx>& element,
451 Dune::index_constant<stokesFaceIdx> domainJ)
const
453 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
454 if (darcyToStokesFaceCouplingStencils_.count(eIdx))
455 return darcyToStokesFaceCouplingStencils_.at(eIdx);
457 return emptyStencil_;
464 template<std::
size_t i, std::
size_t j>
466 const SubControlVolumeFace<stokesIdx>& scvf,
467 Dune::index_constant<j> domainJ)
const
468 {
return emptyStencil_; }
474 const SubControlVolumeFace<stokesIdx>& scvf,
475 Dune::index_constant<darcyIdx> domainJ)
const
477 const auto faceDofIdx = scvf.dofIndex();
478 if(stokesFaceCouplingStencils_.count(faceDofIdx))
479 return stokesFaceCouplingStencils_.at(faceDofIdx);
481 return emptyStencil_;
489 template<
class IdType>
491 {
return emptyStencil_; }
496 template<
class IdType>
498 {
return emptyStencil_; }
503 bool isCoupledEntity(Dune::index_constant<stokesIdx>,
const SubControlVolumeFace<stokesFaceIdx>& scvf)
const
505 return stokesFaceCouplingStencils_.count(scvf.dofIndex());
511 bool isCoupledEntity(Dune::index_constant<darcyIdx>,
const SubControlVolumeFace<darcyIdx>& scvf)
const
520 {
return emptyStencil_; }
524 std::sort(stencil.begin(), stencil.end());
525 stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
530 std::vector<bool> isCoupledDarcyDof_;
531 std::shared_ptr<CouplingData> couplingData_;
533 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesCellCenterCouplingStencils_;
534 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesFaceCouplingStencils_;
535 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesCellCenterCouplingStencils_;
536 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesFaceCouplingStencils_;
537 std::vector<std::size_t> emptyStencil_;
542 mutable std::vector<StationaryStokesCouplingContext> stokesCouplingContext_;
543 mutable std::vector<StationaryDarcyCouplingContext> darcyCouplingContext_;
545 mutable std::size_t boundStokesElemIdx_;
546 mutable std::size_t boundDarcyElemIdx_;
548 CouplingMapper couplingMapper_;
A helper to deduce a vector with the same size as numbers of equations.
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==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:99
Definition: common/properties.hh:102
The type for a global container for the volume variables.
Definition: common/properties.hh:109
The global vector of flux variable containers.
Definition: common/properties.hh:119
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:123
Definition: couplingdata.hh:206
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:50
void computeStencils()
Prepare the coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:149
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:185
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:425
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:269
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:473
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:519
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:409
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:334
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:490
static constexpr auto darcyIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:58
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:370
StokesDarcyCouplingManager(std::shared_ptr< const GridGeometry< stokesIdx > > stokesFvGridGeometry, std::shared_ptr< const GridGeometry< darcyIdx > > darcyFvGridGeometry)
Constructor.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:123
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:224
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:387
void removeDuplicates_(std::vector< std::size_t > &stencil)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:522
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:133
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:434
const auto & couplingData() const
Access the coupling data.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:362
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:283
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:449
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:503
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:218
static constexpr auto stokesCellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:56
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:497
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:178
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:311
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:465
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:511
static constexpr auto stokesFaceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:55
static constexpr auto stokesIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:57
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:131
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:139
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:147
decltype(auto) curSol()
the solution vector of the coupled problem
Definition: multidomain/couplingmanager.hh:370
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:299
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:321
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:231
Base coupling manager for the staggered discretization.
Definition: staggeredcouplingmanager.hh:43
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:150
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:94
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:69
Declares all properties used in Dumux.
The local element solution class for staggered methods.