version 3.10-dev
multidomain/boundary/stokesdarcy/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_STOKES_DARCY_COUPLINGMANAGER_HH
14#define DUMUX_STOKES_DARCY_COUPLINGMANAGER_HH
15
16#include <utility>
17#include <memory>
18
19#include <dune/common/float_cmp.hh>
20#include <dune/common/exceptions.hh>
25
26#include "couplingdata.hh"
27#include "couplingmapper.hh"
28
29namespace Dumux {
30
35template<class MDTraits>
37: public StaggeredCouplingManager<MDTraits>
38{
39 using Scalar = typename MDTraits::Scalar;
41
42public:
43 static constexpr auto stokesFaceIdx = typename MDTraits::template SubDomain<0>::Index();
44 static constexpr auto stokesCellCenterIdx = typename MDTraits::template SubDomain<1>::Index();
45 static constexpr auto stokesIdx = stokesFaceIdx;
46 static constexpr auto darcyIdx = typename MDTraits::template SubDomain<2>::Index();
47
48private:
49
50 using SolutionVector = typename MDTraits::SolutionVector;
51
52 // obtain the type tags of the sub problems
53 using StokesTypeTag = typename MDTraits::template SubDomain<0>::TypeTag;
54 using DarcyTypeTag = typename MDTraits::template SubDomain<2>::TypeTag;
55
56 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
57 using CouplingStencil = CouplingStencils::mapped_type;
58
59 // the sub domain type tags
60 template<std::size_t id>
61 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
62
63 static constexpr bool isCompositional = GetPropType<SubDomainTypeTag<0>, Properties::ModelTraits>::numFluidComponents()> 1;
64
65 template<std::size_t id> using GridView = typename GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>::GridView;
66 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
67 template<std::size_t id> using PrimaryVariables = typename MDTraits::template SubDomain<id>::PrimaryVariables;
68 template<std::size_t id> using NumEqVector = Dumux::NumEqVector<PrimaryVariables<id>>;
69 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
70 template<std::size_t id> using GridVolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>;
71 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
72 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
73 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
74 template<std::size_t id> using ElementBoundaryTypes = GetPropType<SubDomainTypeTag<id>, Properties::ElementBoundaryTypes>;
75 template<std::size_t id> using ElementFluxVariablesCache = typename GetPropType<SubDomainTypeTag<id>, Properties::GridFluxVariablesCache>::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 SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
79
81
82 using VelocityVector = typename Element<stokesIdx>::Geometry::GlobalCoordinate;
83
85
86 struct StationaryStokesCouplingContext
87 {
88 Element<darcyIdx> element;
89 FVElementGeometry<darcyIdx> fvGeometry;
90 std::size_t darcyScvfIdx;
91 std::size_t stokesScvfIdx;
92 VolumeVariables<darcyIdx> volVars;
93 };
94
95 struct StationaryDarcyCouplingContext
96 {
97 Element<stokesIdx> element;
98 FVElementGeometry<stokesIdx> fvGeometry;
99 std::size_t stokesScvfIdx;
100 std::size_t darcyScvfIdx;
101 VelocityVector velocity;
102 VolumeVariables<stokesIdx> volVars;
103 };
104public:
105
107
109 StokesDarcyCouplingManager(std::shared_ptr<const GridGeometry<stokesIdx>> stokesFvGridGeometry,
110 std::shared_ptr<const GridGeometry<darcyIdx>> darcyFvGridGeometry)
111 { }
112
116 // \{
117
119 void init(std::shared_ptr<const Problem<stokesIdx>> stokesProblem,
120 std::shared_ptr<const Problem<darcyIdx>> darcyProblem,
121 const SolutionVector& curSol)
122 {
123 if (Dune::FloatCmp::ne(stokesProblem->gravity(), darcyProblem->spatialParams().gravity({})))
124 DUNE_THROW(Dune::InvalidStateException, "Both models must use the same gravity vector");
125
126 this->setSubProblems(std::make_tuple(stokesProblem, stokesProblem, darcyProblem));
127 this->updateSolution(curSol);
128 couplingData_ = std::make_shared<CouplingData>(*this);
130 }
131
132 // \}
133
136 {
137 couplingMapper_.computeCouplingMapsAndStencils(*this,
138 darcyToStokesCellCenterCouplingStencils_,
139 darcyToStokesFaceCouplingStencils_,
140 stokesCellCenterCouplingStencils_,
141 stokesFaceCouplingStencils_);
142
143 for(auto&& stencil : darcyToStokesCellCenterCouplingStencils_)
144 removeDuplicates_(stencil.second);
145 for(auto&& stencil : darcyToStokesFaceCouplingStencils_)
146 removeDuplicates_(stencil.second);
147 for(auto&& stencil : stokesCellCenterCouplingStencils_)
148 removeDuplicates_(stencil.second);
149 for(auto&& stencil : stokesFaceCouplingStencils_)
150 removeDuplicates_(stencil.second);
151 }
152
156 // \{
157
159
163 template<std::size_t i, class Assembler, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
164 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<stokesCellCenterIdx>& element, const Assembler& assembler) const
165 { bindCouplingContext(domainI, element); }
166
170 template<std::size_t i, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
171 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<stokesCellCenterIdx>& element) const
172 {
173 stokesCouplingContext_.clear();
174
175 const auto stokesElementIdx = this->problem(stokesIdx).gridGeometry().elementMapper().index(element);
176 boundStokesElemIdx_ = stokesElementIdx;
177
178 // do nothing if the element is not coupled to the other domain
179 if(!couplingMapper_.stokesElementToDarcyElementMap().count(stokesElementIdx))
180 return;
181
182 // prepare the coupling context
183 const auto& darcyIndices = couplingMapper_.stokesElementToDarcyElementMap().at(stokesElementIdx);
184 auto darcyFvGeometry = localView(this->problem(darcyIdx).gridGeometry());
185 for(auto&& indices : darcyIndices)
186 {
187 const auto& darcyElement = this->problem(darcyIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
188 darcyFvGeometry.bindElement(darcyElement);
189 const auto& scv = (*scvs(darcyFvGeometry).begin());
190
191 const auto darcyElemSol = elementSolution(darcyElement, this->curSol(darcyIdx), this->problem(darcyIdx).gridGeometry());
192 VolumeVariables<darcyIdx> darcyVolVars;
193 darcyVolVars.update(darcyElemSol, this->problem(darcyIdx), darcyElement, scv);
194
195 // add the context
196 stokesCouplingContext_.push_back({darcyElement, darcyFvGeometry, indices.scvfIdx, indices.flipScvfIdx, darcyVolVars});
197 }
198 }
199
203 template<class Assembler>
204 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI, const Element<darcyIdx>& element, const Assembler& assembler) const
205 { bindCouplingContext(domainI, element); }
206
210 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI, const Element<darcyIdx>& element) const
211 {
212 darcyCouplingContext_.clear();
213
214 const auto darcyElementIdx = this->problem(darcyIdx).gridGeometry().elementMapper().index(element);
215 boundDarcyElemIdx_ = darcyElementIdx;
216
217 // do nothing if the element is not coupled to the other domain
218 if(!couplingMapper_.darcyElementToStokesElementMap().count(darcyElementIdx))
219 return;
220
221 // prepare the coupling context
222 const auto& stokesElementIndices = couplingMapper_.darcyElementToStokesElementMap().at(darcyElementIdx);
223 auto stokesFvGeometry = localView(this->problem(stokesIdx).gridGeometry());
224
225 for(auto&& indices : stokesElementIndices)
226 {
227 const auto& stokesElement = this->problem(stokesIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
228 stokesFvGeometry.bindElement(stokesElement);
229
230 VelocityVector faceVelocity(0.0);
231
232 for(const auto& scvf : scvfs(stokesFvGeometry))
233 {
234 if(scvf.index() == indices.scvfIdx)
235 faceVelocity[scvf.directionIndex()] = this->curSol(stokesFaceIdx)[scvf.dofIndex()];
236 }
237
238 using PriVarsType = typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
239 const auto& cellCenterPriVars = this->curSol(stokesCellCenterIdx)[indices.eIdx];
240 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(cellCenterPriVars);
241
242 VolumeVariables<stokesIdx> stokesVolVars;
243 for(const auto& scv : scvs(stokesFvGeometry))
244 stokesVolVars.update(elemSol, this->problem(stokesIdx), stokesElement, scv);
245
246 // add the context
247 darcyCouplingContext_.push_back({stokesElement, stokesFvGeometry, indices.scvfIdx, indices.flipScvfIdx, faceVelocity, stokesVolVars});
248 }
249 }
250
252
254 template<class LocalAssemblerI>
255 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
256 const LocalAssemblerI& localAssemblerI,
257 Dune::index_constant<darcyIdx> domainJ,
258 std::size_t dofIdxGlobalJ,
259 const PrimaryVariables<darcyIdx>& priVarsJ,
260 int pvIdxJ)
261 {
262 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
263 }
264
266 template<class LocalAssemblerI>
267 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
268 const LocalAssemblerI& localAssemblerI,
269 Dune::index_constant<stokesCellCenterIdx> domainJ,
270 const std::size_t dofIdxGlobalJ,
271 const PrimaryVariables<stokesCellCenterIdx>& priVarsJ,
272 int pvIdxJ)
273 {
274 this->curSol(domainJ)[dofIdxGlobalJ] = priVarsJ;
275
276 for (auto& data : darcyCouplingContext_)
277 {
278 const auto stokesElemIdx = this->problem(stokesIdx).gridGeometry().elementMapper().index(data.element);
279
280 if(stokesElemIdx != dofIdxGlobalJ)
281 continue;
282
283 using PriVarsType = typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
284 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(priVarsJ);
285
286 for(const auto& scv : scvs(data.fvGeometry))
287 data.volVars.update(elemSol, this->problem(stokesIdx), data.element, scv);
288 }
289 }
290
292 template<class LocalAssemblerI>
293 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
294 const LocalAssemblerI& localAssemblerI,
295 Dune::index_constant<stokesFaceIdx> domainJ,
296 const std::size_t dofIdxGlobalJ,
297 const PrimaryVariables<stokesFaceIdx>& priVarsJ,
298 int pvIdxJ)
299 {
300 this->curSol(domainJ)[dofIdxGlobalJ] = priVarsJ;
301
302 for (auto& data : darcyCouplingContext_)
303 {
304 for(const auto& scvf : scvfs(data.fvGeometry))
305 {
306 if(scvf.dofIndex() == dofIdxGlobalJ)
307 data.velocity[scvf.directionIndex()] = priVarsJ;
308 }
309 }
310 }
311
313 template<std::size_t i, class LocalAssemblerI, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
314 void updateCouplingContext(Dune::index_constant<i> domainI,
315 const LocalAssemblerI& localAssemblerI,
316 Dune::index_constant<darcyIdx> domainJ,
317 const std::size_t dofIdxGlobalJ,
318 const PrimaryVariables<darcyIdx>& priVarsJ,
319 int pvIdxJ)
320 {
321 this->curSol(domainJ)[dofIdxGlobalJ] = priVarsJ;
322
323 for (auto& data : stokesCouplingContext_)
324 {
325 const auto darcyElemIdx = this->problem(darcyIdx).gridGeometry().elementMapper().index(data.element);
326
327 if(darcyElemIdx != dofIdxGlobalJ)
328 continue;
329
330 const auto darcyElemSol = elementSolution(data.element, this->curSol(darcyIdx), this->problem(darcyIdx).gridGeometry());
331
332 for(const auto& scv : scvs(data.fvGeometry))
333 data.volVars.update(darcyElemSol, this->problem(darcyIdx), data.element, scv);
334 }
335 }
336
337 // \}
338
342 const auto& couplingData() const
343 {
344 return *couplingData_;
345 }
346
350 const auto& stokesCouplingContext(const Element<stokesIdx>& element, const SubControlVolumeFace<stokesIdx>& scvf) const
351 {
352 if (stokesCouplingContext_.empty() || boundStokesElemIdx_ != scvf.insideScvIdx())
354
355 for(const auto& context : stokesCouplingContext_)
356 {
357 if(scvf.index() == context.stokesScvfIdx)
358 return context;
359 }
360
361 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
362 }
363
367 const auto& darcyCouplingContext(const Element<darcyIdx>& element, const SubControlVolumeFace<darcyIdx>& scvf) const
368 {
369 if (darcyCouplingContext_.empty() || boundDarcyElemIdx_ != scvf.insideScvIdx())
371
372 for(const auto& context : darcyCouplingContext_)
373 {
374 if(scvf.index() == context.darcyScvfIdx)
375 return context;
376 }
377
378 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
379 }
380
384 // \{
385
387
395 const CouplingStencil& couplingStencil(Dune::index_constant<stokesCellCenterIdx> domainI,
396 const Element<stokesIdx>& elementI,
397 Dune::index_constant<darcyIdx> domainJ) const
398 {
399 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(elementI);
400 if(stokesCellCenterCouplingStencils_.count(eIdx))
401 return stokesCellCenterCouplingStencils_.at(eIdx);
402 else
403 return emptyStencil_;
404 }
405
414 template<std::size_t i, std::size_t j>
415 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
416 const Element<i>& elementI,
417 Dune::index_constant<j> domainJ) const
418 { return emptyStencil_; }
419
428 const CouplingStencil& couplingStencil(Dune::index_constant<darcyIdx> domainI,
429 const Element<darcyIdx>& elementI,
430 Dune::index_constant<stokesCellCenterIdx> domainJ) const
431 {
432 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(elementI);
433 if(darcyToStokesCellCenterCouplingStencils_.count(eIdx))
434 return darcyToStokesCellCenterCouplingStencils_.at(eIdx);
435 else
436 return emptyStencil_;
437 }
438
447 const CouplingStencil& couplingStencil(Dune::index_constant<darcyIdx> domainI,
448 const Element<darcyIdx>& elementI,
449 Dune::index_constant<stokesFaceIdx> domainJ) const
450 {
451 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(elementI);
452 if (darcyToStokesFaceCouplingStencils_.count(eIdx))
453 return darcyToStokesFaceCouplingStencils_.at(eIdx);
454 else
455 return emptyStencil_;
456 }
457
466 template<std::size_t i, std::size_t j>
467 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
468 const SubControlVolumeFace<stokesIdx>& scvfI,
469 Dune::index_constant<j> domainJ) const
470 { return emptyStencil_; }
471
479 const CouplingStencil& couplingStencil(Dune::index_constant<stokesFaceIdx> domainI,
480 const SubControlVolumeFace<stokesIdx>& scvfI,
481 Dune::index_constant<darcyIdx> domainJ) const
482 {
483 const auto faceDofIdx = scvfI.dofIndex();
484 if(stokesFaceCouplingStencils_.count(faceDofIdx))
485 return stokesFaceCouplingStencils_.at(faceDofIdx);
486 else
487 return emptyStencil_;
488 }
489
490 // \}
491
495 template<class IdType>
496 const std::vector<std::size_t>& getAdditionalDofDependencies(IdType id, std::size_t stokesElementIdx) const
497 { return emptyStencil_; }
498
502 template<class IdType>
503 const std::vector<std::size_t>& getAdditionalDofDependenciesInverse(IdType id, std::size_t darcyElementIdx) const
504 { return emptyStencil_; }
505
509 bool isCoupledEntity(Dune::index_constant<stokesIdx>, const SubControlVolumeFace<stokesFaceIdx>& scvf) const
510 {
511 return stokesFaceCouplingStencils_.count(scvf.dofIndex());
512 }
513
517 bool isCoupledEntity(Dune::index_constant<darcyIdx>, const SubControlVolumeFace<darcyIdx>& scvf) const
518 {
519 return couplingMapper_.isCoupledDarcyScvf(scvf.index());
520 }
521
522protected:
523
525 std::vector<std::size_t>& emptyStencil()
526 { return emptyStencil_; }
527
528 void removeDuplicates_(std::vector<std::size_t>& stencil)
529 {
530 std::sort(stencil.begin(), stencil.end());
531 stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
532 }
533
534private:
535
536 std::vector<bool> isCoupledDarcyDof_;
537 std::shared_ptr<CouplingData> couplingData_;
538
539 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesCellCenterCouplingStencils_;
540 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesFaceCouplingStencils_;
541 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesCellCenterCouplingStencils_;
542 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesFaceCouplingStencils_;
543 std::vector<std::size_t> emptyStencil_;
544
548 mutable std::vector<StationaryStokesCouplingContext> stokesCouplingContext_;
549 mutable std::vector<StationaryDarcyCouplingContext> darcyCouplingContext_;
550
551 mutable std::size_t boundStokesElemIdx_;
552 mutable std::size_t boundDarcyElemIdx_;
553
554 CouplingMapper couplingMapper_;
555};
556
557} // end namespace Dumux
558
559#endif
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:276
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:298
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:327
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:208
Base coupling manager for the staggered discretization.
Definition: staggeredcouplingmanager.hh:31
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:138
Definition: couplingdata.hh:194
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:38
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const SubControlVolumeFace< stokesIdx > &scvfI, Dune::index_constant< j > domainJ) const
Return the dof indices of a subdomain that influence the residual of a sub-control volume face of the...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:467
void computeStencils()
Prepare the coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:135
const CouplingStencil & couplingStencil(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &elementI, Dune::index_constant< stokesCellCenterIdx > domainJ) const
Return the Stokes cell indices that influence the residual of an element in the Darcy domain.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:428
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:171
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesFaceIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< stokesFaceIdx > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:293
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)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:255
const CouplingStencil & couplingStencil(Dune::index_constant< stokesCellCenterIdx > domainI, const Element< stokesIdx > &elementI, Dune::index_constant< darcyIdx > domainJ) const
The Stokes cell center coupling stencil w.r.t. Darcy DOFs.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:395
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:525
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:496
static constexpr auto darcyIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:46
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:350
StokesDarcyCouplingManager(std::shared_ptr< const GridGeometry< stokesIdx > > stokesFvGridGeometry, std::shared_ptr< const GridGeometry< darcyIdx > > darcyFvGridGeometry)
Constructor.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:109
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:210
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< darcyIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< darcyIdx > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:314
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:367
void removeDuplicates_(std::vector< std::size_t > &stencil)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:528
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:119
const auto & couplingData() const
Access the coupling data.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:342
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesCellCenterIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< stokesCellCenterIdx > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:267
const CouplingStencil & couplingStencil(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &elementI, Dune::index_constant< stokesFaceIdx > domainJ) const
Return the Stokes face indices that influence the residual of an element in the Darcy domain.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:447
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:509
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:204
static constexpr auto stokesCellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:44
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:503
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:164
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:517
static constexpr auto stokesFaceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:43
static constexpr auto stokesIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:45
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &elementI, 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:415
const CouplingStencil & couplingStencil(Dune::index_constant< stokesFaceIdx > domainI, const SubControlVolumeFace< stokesIdx > &scvfI, Dune::index_constant< darcyIdx > domainJ) const
The coupling stencil of a Stokes face w.r.t. Darcy DOFs.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:479
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:30
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:119
void computeCouplingMapsAndStencils(const CouplingManager &couplingManager, Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:44
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:127
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:135
Defines all properties used in Dumux.
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:34
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:296
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
The local element solution class for staggered methods.
The interface of the coupling manager for multi domain problems.