3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_BOUNDARY_FFPM_FFMMOMENTUMPM_COUPLINGMANAGER_HH
26#define DUMUX_MULTIDOMAIN_BOUNDARY_FFPM_FFMMOMENTUMPM_COUPLINGMANAGER_HH
27
28#include <utility>
29#include <memory>
30
31#include <dune/common/float_cmp.hh>
32#include <dune/common/exceptions.hh>
36
37#include "couplingmapper.hh"
38
39namespace Dumux {
40
45template<class MDTraits>
47: public CouplingManager<MDTraits>
48{
49 using Scalar = typename MDTraits::Scalar;
51
52public:
53 static constexpr auto freeFlowMomentumIndex = typename MDTraits::template SubDomain<0>::Index();
54 static constexpr auto porousMediumIndex = typename MDTraits::template SubDomain<1>::Index();
55
57private:
58 // obtain the type tags of the sub problems
59 using FreeFlowMomentumTypeTag = typename MDTraits::template SubDomain<freeFlowMomentumIndex>::TypeTag;
60 using PorousMediumTypeTag = typename MDTraits::template SubDomain<porousMediumIndex>::TypeTag;
61
62 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
63 using CouplingStencil = CouplingStencils::mapped_type;
64
65 // the sub domain type tags
66 template<std::size_t id>
67 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
68
69 template<std::size_t id> using GridView = typename GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>::GridView;
70 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
71 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
72 template<std::size_t id> using GridVolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>;
73 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
74 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
75 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
76 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
77 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
78 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
79 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
80 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
81
82 using VelocityVector = typename Element<freeFlowMomentumIndex>::Geometry::GlobalCoordinate;
83
84 struct FreeFlowMomentumCouplingContext
85 {
86 Element<porousMediumIndex> element;
87 VolumeVariables<porousMediumIndex> volVars;
88 FVElementGeometry<porousMediumIndex> fvGeometry;
89 std::size_t freeFlowMomentumScvfIdx;
90 std::size_t porousMediumScvfIdx;
91 std::size_t porousMediumDofIdx;
92 VelocityVector gravity;
93 };
94
95 struct PorousMediumCouplingContext
96 {
97 Element<freeFlowMomentumIndex> element;
98 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
99 std::size_t porousMediumScvfIdx;
100 std::size_t freeFlowMomentumScvfIdx;
101 std::size_t freeFlowMomentumDofIdx;
102 VelocityVector faceVelocity;
103 };
104
105 using CouplingMapper = FreeFlowMomentumPorousMediumCouplingMapper<MDTraits, FreeFlowMomentumPorousMediumCouplingManager<MDTraits>>;
106
107public:
108
112 // \{
113
115 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> freeFlowMomentumProblem,
116 std::shared_ptr<Problem<porousMediumIndex>> porousMediumProblem,
118 {
119 this->setSubProblems(std::make_tuple(freeFlowMomentumProblem, porousMediumProblem));
120 this->attachSolution(curSol);
121
122 couplingMapper_.update(*this);
123 }
124
125 // \}
126
127
131 // \{
132
136 template<std::size_t i, class Assembler>
137 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<i>& element, const Assembler& assembler) const
138 {
139 bindCouplingContext_(domainI, element);
140 }
141
145 template<std::size_t i, std::size_t j, class LocalAssemblerI>
146 void updateCouplingContext(Dune::index_constant<i> domainI,
147 const LocalAssemblerI& localAssemblerI,
148 Dune::index_constant<j> domainJ,
149 std::size_t dofIdxGlobalJ,
150 const PrimaryVariables<j>& priVarsJ,
151 int pvIdxJ)
152 {
153 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
154
155 // we need to update all solution-depenent components of the coupling context
156 // the dof of domain J has been deflected
157
158 // update the faceVelocity in the PorousMediumCouplingContext
159 if constexpr (domainJ == freeFlowMomentumIndex)
160 {
161 // we only need to update if we are assembling the porous medium domain
162 // since the freeflow domain will not use the velocity from the context
163 if constexpr (domainI == porousMediumIndex)
164 {
165 auto& context = std::get<porousMediumIndex>(couplingContext_);
166 for (auto& c : context)
167 {
168 if (c.freeFlowMomentumDofIdx == dofIdxGlobalJ)
169 {
170 const auto& scvf = localAssemblerI.fvGeometry().scvf(c.porousMediumScvfIdx);
171 c.faceVelocity = faceVelocity(localAssemblerI.element(), scvf);
172 }
173 }
174 }
175 }
176
177 // update volVars in the FreeFlowMomentumCouplingContext
178 else if (domainJ == porousMediumIndex)
179 {
180 auto& context = std::get<freeFlowMomentumIndex>(couplingContext_);
181 for (auto& c : context)
182 {
183 if (c.porousMediumDofIdx == dofIdxGlobalJ)
184 {
185 const auto& ggJ = c.fvGeometry.gridGeometry();
186 const auto& scv = *scvs(c.fvGeometry).begin();
187 const auto elemSol = elementSolution(c.element, this->curSol(domainJ), ggJ);
188 c.volVars.update(elemSol, this->problem(domainJ), c.element, scv);
189 }
190 }
191 }
192 }
193
194 // \}
195
199 template<std::size_t i>
200 const auto& couplingContext(Dune::index_constant<i> domainI,
201 const FVElementGeometry<i>& fvGeometry,
202 const SubControlVolumeFace<i> scvf) const
203 {
204 auto& contexts = std::get<i>(couplingContext_);
205
206 if (contexts.empty() || couplingContextBoundForElement_[i] != scvf.insideScvIdx())
207 bindCouplingContext_(domainI, fvGeometry);
208
209 for (const auto& context : contexts)
210 {
211 const auto expectedScvfIdx = domainI == freeFlowMomentumIndex ? context.freeFlowMomentumScvfIdx : context.porousMediumScvfIdx;
212 if (scvf.index() == expectedScvfIdx)
213 return context;
214 }
215
216 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
217 }
218
222 // \{
223
227 const CouplingStencil& couplingStencil(Dune::index_constant<porousMediumIndex> domainI,
228 const Element<porousMediumIndex>& element,
229 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
230 {
231 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
232 return couplingMapper_.couplingStencil(domainI, eIdx, domainJ);
233 }
234
244 const CouplingStencil& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
245 const Element<freeFlowMomentumIndex>& elementI,
246 const SubControlVolume<freeFlowMomentumIndex>& scvI,
247 Dune::index_constant<porousMediumIndex> domainJ) const
248 {
249 return couplingMapper_.couplingStencil(domainI, elementI, scvI, domainJ);
250 }
251
253
258 template<class LocalAssemblerI, std::size_t j>
259 decltype(auto) evalCouplingResidual(Dune::index_constant<freeFlowMomentumIndex> domainI,
260 const LocalAssemblerI& localAssemblerI,
261 const SubControlVolume<freeFlowMomentumIndex>& scvI,
262 Dune::index_constant<j> domainJ,
263 std::size_t dofIdxGlobalJ) const
264 {
265 return localAssemblerI.evalLocalResidual();
266 }
267
268 // \}
269
273 bool isCoupledLateralScvf(Dune::index_constant<freeFlowMomentumIndex> domainI, const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
274 { return couplingMapper_.isCoupledLateralScvf(domainI, scvf); }
275
279 template<std::size_t i>
280 bool isCoupled(Dune::index_constant<i> domainI, const SubControlVolumeFace<i>& scvf) const
281 {
282 if constexpr (i == freeFlowMomentumIndex)
283 return couplingMapper_.isCoupled(domainI, scvf) || couplingMapper_.isCoupledLateralScvf(domainI, scvf);
284 else
285 return couplingMapper_.isCoupled(domainI, scvf);
286 }
287
293 bool isCoupled(Dune::index_constant<freeFlowMomentumIndex> domainI,
294 const SubControlVolume<freeFlowMomentumIndex>& scv) const
295 { return couplingMapper_.isCoupled(domainI, scv); }
296
300 auto faceVelocity(const Element<porousMediumIndex>& element,
301 const SubControlVolumeFace<porousMediumIndex>& scvf) const
302 {
303 // create a unit normal vector oriented in positive coordinate direction
304 auto velocity = scvf.unitOuterNormal();
305 using std::abs;
306 std::for_each(velocity.begin(), velocity.end(), [](auto& v){ v = abs(v); });
307
308 // create the actual velocity vector
309 velocity *= this->curSol(freeFlowMomentumIndex)[
310 couplingMapper_.outsideDofIndex(Dune::index_constant<porousMediumIndex>(), scvf)
311 ];
312
313 return velocity;
314 }
315
316private:
318 template<std::size_t i>
319 VolumeVariables<i> volVars_(Dune::index_constant<i> domainI,
320 const Element<i>& element,
321 const SubControlVolume<i>& scv) const
322 {
323 VolumeVariables<i> volVars;
324 const auto elemSol = elementSolution(
325 element, this->curSol(domainI), this->problem(domainI).gridGeometry()
326 );
327 volVars.update(elemSol, this->problem(domainI), element, scv);
328 return volVars;
329 }
330
334 template<std::size_t i>
335 bool isCoupledElement_(Dune::index_constant<i> domainI, std::size_t eIdx) const
336 { return couplingMapper_.isCoupledElement(domainI, eIdx); }
337
341 template<std::size_t i>
342 void bindCouplingContext_(Dune::index_constant<i> domainI, const Element<i>& element) const
343 {
344 const auto fvGeometry = localView(this->problem(domainI).gridGeometry()).bindElement(element);
345 bindCouplingContext_(domainI, fvGeometry);
346 }
347
351 template<std::size_t i>
352 void bindCouplingContext_(Dune::index_constant<i> domainI, const FVElementGeometry<i>& fvGeometry) const
353 {
354 auto& context = std::get<domainI>(couplingContext_);
355 context.clear();
356
357 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(fvGeometry.element());
358
359 // do nothing if the element is not coupled to the other domain
360 if (!isCoupledElement_(domainI, eIdx))
361 return;
362
363 couplingContextBoundForElement_[domainI] = eIdx;
364
365 for (const auto& scvf : scvfs(fvGeometry))
366 {
367 if (couplingMapper_.isCoupled(domainI, scvf))
368 {
369 if constexpr (domainI == freeFlowMomentumIndex)
370 {
371 const auto otherElementIdx = couplingMapper_.outsideElementIndex(domainI, scvf);
372 constexpr auto domainJ = Dune::index_constant<1-i>();
373 const auto& otherGridGeometry = this->problem(domainJ).gridGeometry();
374 const auto& otherElement = otherGridGeometry.element(otherElementIdx);
375 auto otherFvGeometry = localView(otherGridGeometry).bindElement(otherElement);
376
377 // there is only one scv for TPFA
378 context.push_back({
379 otherElement,
380 volVars_(domainJ, otherElement, *std::begin(scvs(otherFvGeometry))),
381 std::move(otherFvGeometry),
382 scvf.index(),
383 couplingMapper_.flipScvfIndex(domainI, scvf),
384 otherElementIdx,
385 this->problem(domainJ).spatialParams().gravity(scvf.center())
386 });
387 }
388
389 else if constexpr (domainI == porousMediumIndex)
390 {
391 const auto otherElementIdx = couplingMapper_.outsideElementIndex(domainI, scvf);
392 constexpr auto domainJ = Dune::index_constant<1-i>();
393 const auto& otherGridGeometry = this->problem(domainJ).gridGeometry();
394 const auto& otherElement = otherGridGeometry.element(otherElementIdx);
395 auto otherFvGeometry = localView(otherGridGeometry).bindElement(otherElement);
396
397 const auto otherScvfIdx = couplingMapper_.flipScvfIndex(domainI, scvf);
398 const auto& otherScvf = otherFvGeometry.scvf(otherScvfIdx);
399 const auto& otherScv = otherFvGeometry.scv(otherScvf.insideScvIdx());
400
401 context.push_back({
402 otherElement,
403 std::move(otherFvGeometry),
404 scvf.index(),
405 otherScvfIdx,
406 otherScv.dofIndex(),
407 faceVelocity(fvGeometry.element(), scvf)
408 });
409 }
410 }
411 }
412 }
413
414 mutable std::tuple<std::vector<FreeFlowMomentumCouplingContext>, std::vector<PorousMediumCouplingContext>> couplingContext_;
415 mutable std::array<std::size_t, 2> couplingContextBoundForElement_;
416
417 CouplingMapper couplingMapper_;
418};
419
420} // end namespace Dumux
421
422#endif
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
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:261
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
The type for a global container for the volume variables.
Definition: common/properties.hh:107
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:121
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:48
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:146
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:200
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:273
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:300
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:293
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:53
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:244
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:280
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:227
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:259
static constexpr auto porousMediumIndex
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:54
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:115
typename ParentType::SolutionVectorStorage SolutionVectorStorage
Definition: multidomain/boundary/freeflowporousmedium/ffmomentumpm/couplingmanager.hh:56
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:137
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:272
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:258
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:283
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:302
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:218
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:326
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:314
void update(const CouplingManager &couplingManager)
Main update routine.
Definition: boundary/freeflowporousmedium/ffmomentumpm/couplingmapper.hh:93
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:334
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
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:350
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:83
Declares all properties used in Dumux.
The local element solution class for staggered methods.
The interface of the coupling manager for multi domain problems.