13#ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FFPM_FFMASSPM_COUPLINGMANAGER_HH
14#define DUMUX_MULTIDOMAIN_BOUNDARY_FFPM_FFMASSPM_COUPLINGMANAGER_HH
19#include <dune/common/float_cmp.hh>
20#include <dune/common/exceptions.hh>
32template<
class MDTraits>
36 using Scalar =
typename MDTraits::Scalar;
40 static constexpr auto freeFlowMassIndex =
typename MDTraits::template SubDomain<0>::Index();
41 static constexpr auto porousMediumIndex =
typename MDTraits::template SubDomain<1>::Index();
47 using FreeFlowMassTypeTag =
typename MDTraits::template SubDomain<freeFlowMassIndex>::TypeTag;
48 using PorousMediumTypeTag =
typename MDTraits::template SubDomain<porousMediumIndex>::TypeTag;
50 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
51 using CouplingStencil = CouplingStencils::mapped_type;
54 template<std::
size_t id>
55 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
63 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
65 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
67 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
68 template<std::
size_t id>
using SubControlVolume =
typename FVElementGeometry<id>::SubControlVolume;
70 using VelocityVector =
typename Element<freeFlowMassIndex>::Geometry::GlobalCoordinate;
72 struct FreeFlowMassCouplingContext
74 Element<porousMediumIndex> element;
75 VolumeVariables<porousMediumIndex> volVars;
76 FVElementGeometry<porousMediumIndex> fvGeometry;
78 std::size_t freeFlowMassScvfIdx;
79 std::size_t porousMediumScvfIdx;
80 mutable VelocityVector velocity;
83 struct PorousMediumCouplingContext
85 Element<freeFlowMassIndex> element;
86 VolumeVariables<freeFlowMassIndex> volVars;
87 FVElementGeometry<freeFlowMassIndex> fvGeometry;
89 std::size_t porousMediumScvfIdx;
90 std::size_t freeFlowMassScvfIdx;
91 mutable VelocityVector velocity;
94 using CouplingMapper = DarcyDarcyBoundaryCouplingMapper<MDTraits>;
104 void init(std::shared_ptr<Problem<freeFlowMassIndex>> freeFlowMassProblem,
105 std::shared_ptr<Problem<porousMediumIndex>> darcyProblem,
108 this->
setSubProblems(std::make_tuple(freeFlowMassProblem, darcyProblem));
111 couplingMapper_.
update(*
this);
124 template<std::
size_t i,
class Assembler>
125 void bindCouplingContext(Dune::index_constant<i> domainI,
const Element<i>& element,
const Assembler& assembler)
const
127 bindCouplingContext_(domainI, element);
133 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
135 const LocalAssemblerI& localAssemblerI,
136 Dune::index_constant<j> domainJ,
137 std::size_t dofIdxGlobalJ,
138 const PrimaryVariables<j>& priVarsJ,
141 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
148 auto& context = std::get<1-j>(couplingContext_);
149 for (
auto& c : context)
151 if (c.dofIdx == dofIdxGlobalJ)
153 const auto elemSol =
elementSolution(c.element, this->curSol(domainJ), this->problem(domainJ).gridGeometry());
154 const auto& scv = *scvs(c.fvGeometry).begin();
155 c.volVars.update(elemSol, this->
problem(domainJ), c.element, scv);
165 template<std::
size_t i>
167 const FVElementGeometry<i>& fvGeometry,
168 const SubControlVolumeFace<i> scvf)
const
170 auto& contexts = std::get<i>(couplingContext_);
172 if (contexts.empty() || couplingContextBoundForElement_[i] != scvf.insideScvIdx())
173 bindCouplingContext_(domainI, fvGeometry);
175 for (
const auto& context : contexts)
177 const auto expectedScvfIdx = domainI ==
freeFlowMassIndex ? context.freeFlowMassScvfIdx : context.porousMediumScvfIdx;
178 if (scvf.index() == expectedScvfIdx)
182 DUNE_THROW(Dune::InvalidStateException,
"No coupling context found at scvf " << scvf.center());
193 template<std::
size_t i, std::
size_t j>
195 const Element<i>& element,
196 Dune::index_constant<j> domainJ)
const
198 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
207 template<std::
size_t i>
208 bool isCoupled(Dune::index_constant<i> domainI,
const SubControlVolumeFace<i>& scvf)
const
209 {
return couplingMapper_.
isCoupled(domainI, scvf); }
215 template<std::
size_t i>
216 bool isCoupledElement_(Dune::index_constant<i> domainI, std::size_t eIdx)
const
220 template<std::
size_t i>
221 VolumeVariables<i> volVars_(Dune::index_constant<i> domainI,
222 const Element<i>& element,
223 const SubControlVolume<i>& scv)
const
225 VolumeVariables<i> volVars;
227 element, this->
curSol(domainI), this->
problem(domainI).gridGeometry()
229 volVars.update(elemSol, this->
problem(domainI), element, scv);
236 template<std::
size_t i>
237 void bindCouplingContext_(Dune::index_constant<i> domainI,
const Element<i>& element)
const
239 const auto fvGeometry =
localView(this->
problem(domainI).gridGeometry()).bindElement(element);
240 bindCouplingContext_(domainI, fvGeometry);
246 template<std::
size_t i>
247 void bindCouplingContext_(Dune::index_constant<i> domainI,
const FVElementGeometry<i>& fvGeometry)
const
249 auto& context = std::get<domainI>(couplingContext_);
252 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(fvGeometry.element());
255 if (!isCoupledElement_(domainI, eIdx))
258 couplingContextBoundForElement_[domainI] = eIdx;
260 for (
const auto& scvf : scvfs(fvGeometry))
265 constexpr auto domainJ = Dune::index_constant<1-domainI>();
266 const auto& otherGridGeometry = this->
problem(domainJ).gridGeometry();
267 const auto& otherElement = otherGridGeometry.element(otherElementIdx);
268 auto otherFvGeometry =
localView(otherGridGeometry).bindElement(otherElement);
273 volVars_(domainJ, otherElement, *std::begin(scvs(otherFvGeometry))),
274 std::move(otherFvGeometry),
284 mutable std::tuple<std::vector<FreeFlowMassCouplingContext>, std::vector<PorousMediumCouplingContext>> couplingContext_;
285 mutable std::array<std::size_t, 2> couplingContextBoundForElement_;
287 CouplingMapper couplingMapper_;
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:322
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:287
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:338
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:71
bool isCoupledElement(Dune::index_constant< i >, std::size_t eIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/darcydarcy/couplingmapper.hh:178
void update(const CouplingManager &couplingManager)
Main update routine.
Definition: boundary/darcydarcy/couplingmapper.hh:69
const std::vector< std::size_t > & couplingStencil(Dune::index_constant< i > domainI, const std::size_t eIdxI, Dune::index_constant< j > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: boundary/darcydarcy/couplingmapper.hh:163
std::size_t outsideElementIndex(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Return the outside element index (the element index of the other domain)
Definition: boundary/darcydarcy/couplingmapper.hh:211
bool isCoupled(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
If the boundary entity is on a coupling boundary.
Definition: boundary/darcydarcy/couplingmapper.hh:187
std::size_t flipScvfIndex(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Return the scvf index of the flipped scvf in the other domain.
Definition: boundary/darcydarcy/couplingmapper.hh:199
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:35
void bindCouplingContext(Dune::index_constant< i > domainI, const Element< i > &element, const Assembler &assembler) const
Methods to be accessed by the assembly.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:125
void init(std::shared_ptr< Problem< freeFlowMassIndex > > freeFlowMassProblem, std::shared_ptr< Problem< porousMediumIndex > > darcyProblem, SolutionVectorStorage &curSol)
Methods to be accessed by main.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:104
static constexpr auto porousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:41
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &element, Dune::index_constant< j > domainJ) const
The coupling stencils.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:194
bool isCoupled(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Returns whether a given scvf is coupled to the other domain.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:208
static constexpr auto freeFlowMassIndex
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:40
const auto & couplingContext(Dune::index_constant< i > domainI, const FVElementGeometry< i > &fvGeometry, const SubControlVolumeFace< i > scvf) const
Access the coupling context needed for the Stokes domain.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:166
typename ParentType::SolutionVectorStorage SolutionVectorStorage
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:43
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
Update the coupling context.
Definition: multidomain/boundary/freeflowporousmedium/ffmasspm/couplingmanager.hh:134
Defines all properties used in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:267
The interface of the coupling manager for multi domain problems.
The local element solution class for staggered methods.