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