version 3.7
multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FFPM_FFMMOMENTUMPM_COUPLINGMANAGER_HH
14#define DUMUX_MULTIDOMAIN_BOUNDARY_FFPM_FFMMOMENTUMPM_COUPLINGMANAGER_HH
15
16#include <utility>
17#include <memory>
18
19#include <dune/common/float_cmp.hh>
20#include <dune/common/exceptions.hh>
24
25#include "couplingmapper.hh"
26
27namespace Dumux {
28
33template<class MDTraits>
35: public CouplingManager<MDTraits>
36{
37 using Scalar = typename MDTraits::Scalar;
39
40public:
41 static constexpr auto freeFlowMomentumIndex = typename MDTraits::template SubDomain<0>::Index();
42 static constexpr auto porousMediumIndex = typename MDTraits::template SubDomain<1>::Index();
43
45private:
46 // obtain the type tags of the sub problems
47 using FreeFlowMomentumTypeTag = typename MDTraits::template SubDomain<freeFlowMomentumIndex>::TypeTag;
48 using PorousMediumTypeTag = typename MDTraits::template SubDomain<porousMediumIndex>::TypeTag;
49
50 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
51 using CouplingStencil = CouplingStencils::mapped_type;
52
53 // the sub domain type tags
54 template<std::size_t id>
55 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
56
57 template<std::size_t id> using GridView = typename GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>::GridView;
58 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
59 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
60 template<std::size_t id> using GridVolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>;
61 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
62 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
63 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
64 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
65 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
66 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
67 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
68 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
69
70 using VelocityVector = typename Element<freeFlowMomentumIndex>::Geometry::GlobalCoordinate;
71
72 struct FreeFlowMomentumCouplingContext
73 {
74 Element<porousMediumIndex> element;
75 VolumeVariables<porousMediumIndex> volVars;
76 FVElementGeometry<porousMediumIndex> fvGeometry;
77 std::size_t freeFlowMomentumScvfIdx;
78 std::size_t porousMediumScvfIdx;
79 std::size_t porousMediumDofIdx;
80 VelocityVector gravity;
81 };
82
83 struct PorousMediumCouplingContext
84 {
85 Element<freeFlowMomentumIndex> element;
86 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
87 std::size_t porousMediumScvfIdx;
88 std::size_t freeFlowMomentumScvfIdx;
89 std::size_t freeFlowMomentumDofIdx;
90 VelocityVector faceVelocity;
91 };
92
93 using CouplingMapper = FreeFlowMomentumPorousMediumCouplingMapper<MDTraits, FreeFlowMomentumPorousMediumCouplingManager<MDTraits>>;
94
95public:
96
100 // \{
101
103 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> freeFlowMomentumProblem,
104 std::shared_ptr<Problem<porousMediumIndex>> porousMediumProblem,
106 {
107 this->setSubProblems(std::make_tuple(freeFlowMomentumProblem, porousMediumProblem));
108 this->attachSolution(curSol);
109
110 couplingMapper_.update(*this);
111 }
112
113 // \}
114
115
119 // \{
120
124 template<std::size_t i, class Assembler>
125 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<i>& element, const Assembler& assembler) const
126 {
127 bindCouplingContext_(domainI, element);
128 }
129
133 template<std::size_t i, std::size_t j, class LocalAssemblerI>
134 void updateCouplingContext(Dune::index_constant<i> domainI,
135 const LocalAssemblerI& localAssemblerI,
136 Dune::index_constant<j> domainJ,
137 std::size_t dofIdxGlobalJ,
138 const PrimaryVariables<j>& priVarsJ,
139 int pvIdxJ)
140 {
141 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
142
143 // we need to update all solution-depenent components of the coupling context
144 // the dof of domain J has been deflected
145
146 // update the faceVelocity in the PorousMediumCouplingContext
147 if constexpr (domainJ == freeFlowMomentumIndex)
148 {
149 // we only need to update if we are assembling the porous medium domain
150 // since the freeflow domain will not use the velocity from the context
151 if constexpr (domainI == porousMediumIndex)
152 {
153 auto& context = std::get<porousMediumIndex>(couplingContext_);
154 for (auto& c : context)
155 {
156 if (c.freeFlowMomentumDofIdx == dofIdxGlobalJ)
157 {
158 const auto& scvf = localAssemblerI.fvGeometry().scvf(c.porousMediumScvfIdx);
159 c.faceVelocity = faceVelocity(localAssemblerI.element(), scvf);
160 }
161 }
162 }
163 }
164
165 // update volVars in the FreeFlowMomentumCouplingContext
166 else if (domainJ == porousMediumIndex)
167 {
168 auto& context = std::get<freeFlowMomentumIndex>(couplingContext_);
169 for (auto& c : context)
170 {
171 if (c.porousMediumDofIdx == dofIdxGlobalJ)
172 {
173 const auto& ggJ = c.fvGeometry.gridGeometry();
174 const auto& scv = *scvs(c.fvGeometry).begin();
175 const auto elemSol = elementSolution(c.element, this->curSol(domainJ), ggJ);
176 c.volVars.update(elemSol, this->problem(domainJ), c.element, scv);
177 }
178 }
179 }
180 }
181
182 // \}
183
187 template<std::size_t i>
188 const auto& couplingContext(Dune::index_constant<i> domainI,
189 const FVElementGeometry<i>& fvGeometry,
190 const SubControlVolumeFace<i> scvf) const
191 {
192 auto& contexts = std::get<i>(couplingContext_);
193
194 if (contexts.empty() || couplingContextBoundForElement_[i] != scvf.insideScvIdx())
195 bindCouplingContext_(domainI, fvGeometry);
196
197 for (const auto& context : contexts)
198 {
199 const auto expectedScvfIdx = domainI == freeFlowMomentumIndex ? context.freeFlowMomentumScvfIdx : context.porousMediumScvfIdx;
200 if (scvf.index() == expectedScvfIdx)
201 return context;
202 }
203
204 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
205 }
206
210 // \{
211
215 const CouplingStencil& couplingStencil(Dune::index_constant<porousMediumIndex> domainI,
216 const Element<porousMediumIndex>& element,
217 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
218 {
219 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
220 return couplingMapper_.couplingStencil(domainI, eIdx, domainJ);
221 }
222
232 const CouplingStencil& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
233 const Element<freeFlowMomentumIndex>& elementI,
234 const SubControlVolume<freeFlowMomentumIndex>& scvI,
235 Dune::index_constant<porousMediumIndex> domainJ) const
236 {
237 return couplingMapper_.couplingStencil(domainI, elementI, scvI, domainJ);
238 }
239
241
246 template<class LocalAssemblerI, std::size_t j>
247 decltype(auto) evalCouplingResidual(Dune::index_constant<freeFlowMomentumIndex> domainI,
248 const LocalAssemblerI& localAssemblerI,
249 const SubControlVolume<freeFlowMomentumIndex>& scvI,
250 Dune::index_constant<j> domainJ,
251 std::size_t dofIdxGlobalJ) const
252 {
253 return localAssemblerI.evalLocalResidual();
254 }
255
256 // \}
257
261 bool isCoupledLateralScvf(Dune::index_constant<freeFlowMomentumIndex> domainI, const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
262 { return couplingMapper_.isCoupledLateralScvf(domainI, scvf); }
263
267 template<std::size_t i>
268 bool isCoupled(Dune::index_constant<i> domainI, const SubControlVolumeFace<i>& scvf) const
269 {
270 if constexpr (i == freeFlowMomentumIndex)
271 return couplingMapper_.isCoupled(domainI, scvf) || couplingMapper_.isCoupledLateralScvf(domainI, scvf);
272 else
273 return couplingMapper_.isCoupled(domainI, scvf);
274 }
275
281 bool isCoupled(Dune::index_constant<freeFlowMomentumIndex> domainI,
282 const SubControlVolume<freeFlowMomentumIndex>& scv) const
283 { return couplingMapper_.isCoupled(domainI, scv); }
284
288 auto faceVelocity(const Element<porousMediumIndex>& element,
289 const SubControlVolumeFace<porousMediumIndex>& scvf) const
290 {
291 // create a unit normal vector oriented in positive coordinate direction
292 auto velocity = scvf.unitOuterNormal();
293 using std::abs;
294 std::for_each(velocity.begin(), velocity.end(), [](auto& v){ v = abs(v); });
295
296 // create the actual velocity vector
297 velocity *= this->curSol(freeFlowMomentumIndex)[
298 couplingMapper_.outsideDofIndex(Dune::index_constant<porousMediumIndex>(), scvf)
299 ];
300
301 return velocity;
302 }
303
304private:
306 template<std::size_t i>
307 VolumeVariables<i> volVars_(Dune::index_constant<i> domainI,
308 const Element<i>& element,
309 const SubControlVolume<i>& scv) const
310 {
311 VolumeVariables<i> volVars;
312 const auto elemSol = elementSolution(
313 element, this->curSol(domainI), this->problem(domainI).gridGeometry()
314 );
315 volVars.update(elemSol, this->problem(domainI), element, scv);
316 return volVars;
317 }
318
322 template<std::size_t i>
323 bool isCoupledElement_(Dune::index_constant<i> domainI, std::size_t eIdx) const
324 { return couplingMapper_.isCoupledElement(domainI, eIdx); }
325
329 template<std::size_t i>
330 void bindCouplingContext_(Dune::index_constant<i> domainI, const Element<i>& element) const
331 {
332 const auto fvGeometry = localView(this->problem(domainI).gridGeometry()).bindElement(element);
333 bindCouplingContext_(domainI, fvGeometry);
334 }
335
339 template<std::size_t i>
340 void bindCouplingContext_(Dune::index_constant<i> domainI, const FVElementGeometry<i>& fvGeometry) const
341 {
342 auto& context = std::get<domainI>(couplingContext_);
343 context.clear();
344
345 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(fvGeometry.element());
346
347 // do nothing if the element is not coupled to the other domain
348 if (!isCoupledElement_(domainI, eIdx))
349 return;
350
351 couplingContextBoundForElement_[domainI] = eIdx;
352
353 for (const auto& scvf : scvfs(fvGeometry))
354 {
355 if (couplingMapper_.isCoupled(domainI, scvf))
356 {
357 if constexpr (domainI == freeFlowMomentumIndex)
358 {
359 const auto otherElementIdx = couplingMapper_.outsideElementIndex(domainI, scvf);
360 constexpr auto domainJ = Dune::index_constant<1-i>();
361 const auto& otherGridGeometry = this->problem(domainJ).gridGeometry();
362 const auto& otherElement = otherGridGeometry.element(otherElementIdx);
363 auto otherFvGeometry = localView(otherGridGeometry).bindElement(otherElement);
364
365 // there is only one scv for TPFA
366 context.push_back({
367 otherElement,
368 volVars_(domainJ, otherElement, *std::begin(scvs(otherFvGeometry))),
369 std::move(otherFvGeometry),
370 scvf.index(),
371 couplingMapper_.flipScvfIndex(domainI, scvf),
372 otherElementIdx,
373 this->problem(domainJ).spatialParams().gravity(scvf.center())
374 });
375 }
376
377 else if constexpr (domainI == porousMediumIndex)
378 {
379 const auto otherElementIdx = couplingMapper_.outsideElementIndex(domainI, scvf);
380 constexpr auto domainJ = Dune::index_constant<1-i>();
381 const auto& otherGridGeometry = this->problem(domainJ).gridGeometry();
382 const auto& otherElement = otherGridGeometry.element(otherElementIdx);
383 auto otherFvGeometry = localView(otherGridGeometry).bindElement(otherElement);
384
385 const auto otherScvfIdx = couplingMapper_.flipScvfIndex(domainI, scvf);
386 const auto& otherScvf = otherFvGeometry.scvf(otherScvfIdx);
387 const auto& otherScv = otherFvGeometry.scv(otherScvf.insideScvIdx());
388
389 context.push_back({
390 otherElement,
391 std::move(otherFvGeometry),
392 scvf.index(),
393 otherScvfIdx,
394 otherScv.dofIndex(),
395 faceVelocity(fvGeometry.element(), scvf)
396 });
397 }
398 }
399 }
400 }
401
402 mutable std::tuple<std::vector<FreeFlowMomentumCouplingContext>, std::vector<PorousMediumCouplingContext>> couplingContext_;
403 mutable std::array<std::size_t, 2> couplingContextBoundForElement_;
404
405 CouplingMapper couplingMapper_;
406};
407
408} // end namespace Dumux
409
410#endif
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
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:36
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/ffmomentumpm/couplingmanager.hh:134
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/ffmomentumpm/couplingmanager.hh:188
bool isCoupledLateralScvf(Dune::index_constant< freeFlowMomentumIndex > domainI, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns whether a given scvf is coupled to the other domain.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:261
auto faceVelocity(const Element< porousMediumIndex > &element, const SubControlVolumeFace< porousMediumIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:288
bool isCoupled(Dune::index_constant< freeFlowMomentumIndex > domainI, const SubControlVolume< freeFlowMomentumIndex > &scv) const
If the boundary entity is on a coupling boundary.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:281
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:41
const CouplingStencil & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< porousMediumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:232
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/ffmomentumpm/couplingmanager.hh:268
const CouplingStencil & couplingStencil(Dune::index_constant< porousMediumIndex > domainI, const Element< porousMediumIndex > &element, Dune::index_constant< freeFlowMomentumIndex > domainJ) const
The coupling stencils.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:215
decltype(auto) evalCouplingResidual(Dune::index_constant< freeFlowMomentumIndex > domainI, const LocalAssemblerI &localAssemblerI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluate the coupling residual special interface for fcstaggered methods
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:247
static constexpr auto porousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:42
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > freeFlowMomentumProblem, std::shared_ptr< Problem< porousMediumIndex > > porousMediumProblem, SolutionVectorStorage &curSol)
Methods to be accessed by main.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:103
typename ParentType::SolutionVectorStorage SolutionVectorStorage
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:44
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/ffmomentumpm/couplingmanager.hh:125
bool isCoupled(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
If the boundary entity is on a coupling boundary.
Definition: boundary/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:260
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/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:246
bool isCoupledLateralScvf(Dune::index_constant< CouplingManager::freeFlowMomentumIndex > domainI, const SubControlVolumeFace< CouplingManager::freeFlowMomentumIndex > &scvf) const
If the boundary entity is on a coupling boundary.
Definition: boundary/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:271
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/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:290
const std::vector< std::size_t > & couplingStencil(Dune::index_constant< CouplingManager::porousMediumIndex > domainI, const std::size_t eIdxI, Dune::index_constant< CouplingManager::freeFlowMomentumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: boundary/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:206
std::size_t outsideDofIndex(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
Return the outside element index (the element index of the other domain)
Definition: boundary/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:314
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/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:302
void update(const CouplingManager &couplingManager)
Main update routine.
Definition: boundary/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:81
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
decltype(auto) evalCouplingResidual(Dune::index_constant< i > 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: multidomain/couplingmanager.hh:249
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.
Definition: adapt.hh:17
The local element solution class for staggered methods.