3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_STOKES_DARCY_COUPLINGMANAGER_HH
26#define DUMUX_STOKES_DARCY_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 "couplingdata.hh"
38#include "couplingmapper.hh"
39
40namespace Dumux {
41
46template<class MDTraits>
48: public StaggeredCouplingManager<MDTraits>
49{
50 using Scalar = typename MDTraits::Scalar;
52
53public:
54 static constexpr auto stokesCellCenterIdx = typename MDTraits::template SubDomain<0>::Index();
55 static constexpr auto stokesFaceIdx = typename MDTraits::template SubDomain<1>::Index();
56 static constexpr auto cellCenterIdx = typename MDTraits::template SubDomain<0>::Index();
57 static constexpr auto faceIdx = typename MDTraits::template SubDomain<1>::Index();
58 static constexpr auto stokesIdx = stokesCellCenterIdx;
59 static constexpr auto darcyIdx = typename MDTraits::template SubDomain<2>::Index();
60
61private:
62
63 using SolutionVector = typename MDTraits::SolutionVector;
64
65 // obtain the type tags of the sub problems
66 using StokesTypeTag = typename MDTraits::template SubDomain<0>::TypeTag;
67 using DarcyTypeTag = typename MDTraits::template SubDomain<2>::TypeTag;
68
69 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
70 using CouplingStencil = CouplingStencils::mapped_type;
71
72 // the sub domain type tags
73 template<std::size_t id>
74 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
75
76 static constexpr bool isCompositional = GetPropType<SubDomainTypeTag<0>, Properties::ModelTraits>::numFluidComponents()> 1;
77
78 template<std::size_t id> using GridView = GetPropType<SubDomainTypeTag<id>, Properties::GridView>;
79 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
80 template<std::size_t id> using NumEqVector = GetPropType<SubDomainTypeTag<id>, Properties::NumEqVector>;
81 template<std::size_t id> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::LocalView;
82 template<std::size_t id> using GridVolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>;
83 template<std::size_t id> using VolumeVariables = typename GetPropType<SubDomainTypeTag<id>, Properties::GridVolumeVariables>::VolumeVariables;
84 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
85 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
86 template<std::size_t id> using ElementBoundaryTypes = GetPropType<SubDomainTypeTag<id>, Properties::ElementBoundaryTypes>;
87 template<std::size_t id> using ElementFluxVariablesCache = typename GetPropType<SubDomainTypeTag<id>, Properties::GridFluxVariablesCache>::LocalView;
88 template<std::size_t id> using GridVariables = GetPropType<SubDomainTypeTag<id>, Properties::GridVariables>;
89 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
90 template<std::size_t id> using PrimaryVariables = typename MDTraits::template SubDomain<id>::PrimaryVariables;
91 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
92
94
95 using VelocityVector = typename Element<stokesIdx>::Geometry::GlobalCoordinate;
96
98
99 struct StationaryStokesCouplingContext
100 {
101 Element<darcyIdx> element;
102 FVElementGeometry<darcyIdx> fvGeometry;
103 std::size_t darcyScvfIdx;
104 std::size_t stokesScvfIdx;
105 VolumeVariables<darcyIdx> volVars;
106 };
107
108 struct StationaryDarcyCouplingContext
109 {
110 Element<stokesIdx> element;
111 FVElementGeometry<stokesIdx> fvGeometry;
112 std::size_t stokesScvfIdx;
113 std::size_t darcyScvfIdx;
114 VelocityVector velocity;
115 VolumeVariables<stokesIdx> volVars;
116 };
117public:
118
122
124 StokesDarcyCouplingManager(std::shared_ptr<const GridGeometry<stokesIdx>> stokesFvGridGeometry,
125 std::shared_ptr<const GridGeometry<darcyIdx>> darcyFvGridGeometry) : couplingMapper_(*this)
126 { }
127
131 // \{
132
134 void init(std::shared_ptr<const Problem<stokesIdx>> stokesProblem,
135 std::shared_ptr<const Problem<darcyIdx>> darcyProblem,
136 const SolutionVector& curSol)
137 {
138 if (Dune::FloatCmp::ne(stokesProblem->gravity(), darcyProblem->spatialParams().gravity({})))
139 DUNE_THROW(Dune::InvalidStateException, "Both models must use the same gravity vector");
140
141 this->setSubProblems(std::make_tuple(stokesProblem, stokesProblem, darcyProblem));
142 this->curSol() = curSol;
143 couplingData_ = std::make_shared<CouplingData>(*this);
145 }
146
148 void update()
149 { }
150
151 // \}
152
154 void updateSolution(const SolutionVector& curSol)
155 { this->curSol() = curSol; }
156
159 {
160 couplingMapper_.computeCouplingMapsAndStencils(darcyToStokesCellCenterCouplingStencils_,
161 darcyToStokesFaceCouplingStencils_,
162 stokesCellCenterCouplingStencils_,
163 stokesFaceCouplingStencils_);
164
165 for(auto&& stencil : darcyToStokesCellCenterCouplingStencils_)
166 removeDuplicates_(stencil.second);
167 for(auto&& stencil : darcyToStokesFaceCouplingStencils_)
168 removeDuplicates_(stencil.second);
169 for(auto&& stencil : stokesCellCenterCouplingStencils_)
170 removeDuplicates_(stencil.second);
171 for(auto&& stencil : stokesFaceCouplingStencils_)
172 removeDuplicates_(stencil.second);
173 }
174
178 // \{
179
181
185 template<std::size_t i, class Assembler, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
186 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<stokesCellCenterIdx>& element, const Assembler& assembler) const
187 { bindCouplingContext(domainI, element); }
188
192 template<std::size_t i, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
193 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<stokesCellCenterIdx>& element) const
194 {
195 stokesCouplingContext_.clear();
196
197 const auto stokesElementIdx = this->problem(stokesIdx).gridGeometry().elementMapper().index(element);
198 boundStokesElemIdx_ = stokesElementIdx;
199
200 // do nothing if the element is not coupled to the other domain
201 if(!couplingMapper_.stokesElementToDarcyElementMap().count(stokesElementIdx))
202 return;
203
204 // prepare the coupling context
205 const auto& darcyIndices = couplingMapper_.stokesElementToDarcyElementMap().at(stokesElementIdx);
206 auto darcyFvGeometry = localView(this->problem(darcyIdx).gridGeometry());
207
208 for(auto&& indices : darcyIndices)
209 {
210 const auto& darcyElement = this->problem(darcyIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
211 darcyFvGeometry.bindElement(darcyElement);
212 const auto& scv = (*scvs(darcyFvGeometry).begin());
213
214 const auto darcyElemSol = elementSolution(darcyElement, this->curSol()[darcyIdx], this->problem(darcyIdx).gridGeometry());
215 VolumeVariables<darcyIdx> darcyVolVars;
216 darcyVolVars.update(darcyElemSol, this->problem(darcyIdx), darcyElement, scv);
217
218 // add the context
219 stokesCouplingContext_.push_back({darcyElement, darcyFvGeometry, indices.scvfIdx, indices.flipScvfIdx, darcyVolVars});
220 }
221 }
222
226 template<class Assembler>
227 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI, const Element<darcyIdx>& element, const Assembler& assembler) const
228 { bindCouplingContext(domainI, element); }
229
233 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI, const Element<darcyIdx>& element) const
234 {
235 darcyCouplingContext_.clear();
236
237 const auto darcyElementIdx = this->problem(darcyIdx).gridGeometry().elementMapper().index(element);
238 boundDarcyElemIdx_ = darcyElementIdx;
239
240 // do nothing if the element is not coupled to the other domain
241 if(!couplingMapper_.darcyElementToStokesElementMap().count(darcyElementIdx))
242 return;
243
244 // prepare the coupling context
245 const auto& stokesElementIndices = couplingMapper_.darcyElementToStokesElementMap().at(darcyElementIdx);
246 auto stokesFvGeometry = localView(this->problem(stokesIdx).gridGeometry());
247
248 for(auto&& indices : stokesElementIndices)
249 {
250 const auto& stokesElement = this->problem(stokesIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
251 stokesFvGeometry.bindElement(stokesElement);
252
253 VelocityVector faceVelocity(0.0);
254
255 for(const auto& scvf : scvfs(stokesFvGeometry))
256 {
257 if(scvf.index() == indices.scvfIdx)
258 faceVelocity[scvf.directionIndex()] = this->curSol()[stokesFaceIdx][scvf.dofIndex()];
259 }
260
261 using PriVarsType = typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
262 const auto& cellCenterPriVars = this->curSol()[stokesCellCenterIdx][indices.eIdx];
263 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(cellCenterPriVars);
264
265 VolumeVariables<stokesIdx> stokesVolVars;
266 for(const auto& scv : scvs(stokesFvGeometry))
267 stokesVolVars.update(elemSol, this->problem(stokesIdx), stokesElement, scv);
268
269 // add the context
270 darcyCouplingContext_.push_back({stokesElement, stokesFvGeometry, indices.scvfIdx, indices.flipScvfIdx, faceVelocity, stokesVolVars});
271 }
272 }
273
277 template<class LocalAssemblerI>
278 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
279 const LocalAssemblerI& localAssemblerI,
280 Dune::index_constant<darcyIdx> domainJ,
281 std::size_t dofIdxGlobalJ,
282 const PrimaryVariables<darcyIdx>& priVarsJ,
283 int pvIdxJ)
284 {
285 this->curSol()[domainJ][dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
286 }
287
291 template<class LocalAssemblerI>
292 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
293 const LocalAssemblerI& localAssemblerI,
294 Dune::index_constant<stokesCellCenterIdx> domainJ,
295 const std::size_t dofIdxGlobalJ,
296 const PrimaryVariables<stokesCellCenterIdx>& priVars,
297 int pvIdxJ)
298 {
299 this->curSol()[domainJ][dofIdxGlobalJ] = priVars;
300
301 for (auto& data : darcyCouplingContext_)
302 {
303 const auto stokesElemIdx = this->problem(stokesIdx).gridGeometry().elementMapper().index(data.element);
304
305 if(stokesElemIdx != dofIdxGlobalJ)
306 continue;
307
308 using PriVarsType = typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
309 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(priVars);
310
311 for(const auto& scv : scvs(data.fvGeometry))
312 data.volVars.update(elemSol, this->problem(stokesIdx), data.element, scv);
313 }
314 }
315
319 template<class LocalAssemblerI>
320 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
321 const LocalAssemblerI& localAssemblerI,
322 Dune::index_constant<stokesFaceIdx> domainJ,
323 const std::size_t dofIdxGlobalJ,
324 const PrimaryVariables<1>& priVars,
325 int pvIdxJ)
326 {
327 this->curSol()[domainJ][dofIdxGlobalJ] = priVars;
328
329 for (auto& data : darcyCouplingContext_)
330 {
331 for(const auto& scvf : scvfs(data.fvGeometry))
332 {
333 if(scvf.dofIndex() == dofIdxGlobalJ)
334 data.velocity[scvf.directionIndex()] = priVars;
335 }
336 }
337 }
338
342 template<std::size_t i, class LocalAssemblerI, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
343 void updateCouplingContext(Dune::index_constant<i> domainI,
344 const LocalAssemblerI& localAssemblerI,
345 Dune::index_constant<darcyIdx> domainJ,
346 const std::size_t dofIdxGlobalJ,
347 const PrimaryVariables<darcyIdx>& priVars,
348 int pvIdxJ)
349 {
350 this->curSol()[domainJ][dofIdxGlobalJ] = priVars;
351
352 for (auto& data : stokesCouplingContext_)
353 {
354 const auto darcyElemIdx = this->problem(darcyIdx).gridGeometry().elementMapper().index(data.element);
355
356 if(darcyElemIdx != dofIdxGlobalJ)
357 continue;
358
359 const auto darcyElemSol = elementSolution(data.element, this->curSol()[darcyIdx], this->problem(darcyIdx).gridGeometry());
360
361 for(const auto& scv : scvs(data.fvGeometry))
362 data.volVars.update(darcyElemSol, this->problem(darcyIdx), data.element, scv);
363 }
364 }
365
366 // \}
367
371 const auto& couplingData() const
372 {
373 return *couplingData_;
374 }
375
379 const auto& stokesCouplingContext(const Element<stokesIdx>& element, const SubControlVolumeFace<stokesIdx>& scvf) const
380 {
381 if (stokesCouplingContext_.empty() || boundStokesElemIdx_ != scvf.insideScvIdx())
383
384 for(const auto& context : stokesCouplingContext_)
385 {
386 if(scvf.index() == context.stokesScvfIdx)
387 return context;
388 }
389
390 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
391 }
392
396 const auto& darcyCouplingContext(const Element<darcyIdx>& element, const SubControlVolumeFace<darcyIdx>& scvf) const
397 {
398 if (darcyCouplingContext_.empty() || boundDarcyElemIdx_ != scvf.insideScvIdx())
400
401 for(const auto& context : darcyCouplingContext_)
402 {
403 if(scvf.index() == context.darcyScvfIdx)
404 return context;
405 }
406
407 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
408 }
409
413 // \{
414
418 const CouplingStencil& couplingStencil(Dune::index_constant<stokesCellCenterIdx> domainI,
419 const Element<stokesIdx>& element,
420 Dune::index_constant<darcyIdx> domainJ) const
421 {
422 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
423 if(stokesCellCenterCouplingStencils_.count(eIdx))
424 return stokesCellCenterCouplingStencils_.at(eIdx);
425 else
426 return emptyStencil_;
427 }
428
433 template<std::size_t i, std::size_t j>
434 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
435 const Element<i>& element,
436 Dune::index_constant<j> domainJ) const
437 { return emptyStencil_; }
438
443 const CouplingStencil& couplingStencil(Dune::index_constant<darcyIdx> domainI,
444 const Element<darcyIdx>& element,
445 Dune::index_constant<stokesCellCenterIdx> domainJ) const
446 {
447 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
448 if(darcyToStokesCellCenterCouplingStencils_.count(eIdx))
449 return darcyToStokesCellCenterCouplingStencils_.at(eIdx);
450 else
451 return emptyStencil_;
452 }
453
458 const CouplingStencil& couplingStencil(Dune::index_constant<darcyIdx> domainI,
459 const Element<darcyIdx>& element,
460 Dune::index_constant<stokesFaceIdx> domainJ) const
461 {
462 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
463 if (darcyToStokesFaceCouplingStencils_.count(eIdx))
464 return darcyToStokesFaceCouplingStencils_.at(eIdx);
465 else
466 return emptyStencil_;
467 }
468
473 template<std::size_t i, std::size_t j>
474 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
475 const SubControlVolumeFace<stokesIdx>& scvf,
476 Dune::index_constant<j> domainJ) const
477 { return emptyStencil_; }
478
482 const CouplingStencil& couplingStencil(Dune::index_constant<stokesFaceIdx> domainI,
483 const SubControlVolumeFace<stokesIdx>& scvf,
484 Dune::index_constant<darcyIdx> domainJ) const
485 {
486 const auto faceDofIdx = scvf.dofIndex();
487 if(stokesFaceCouplingStencils_.count(faceDofIdx))
488 return stokesFaceCouplingStencils_.at(faceDofIdx);
489 else
490 return emptyStencil_;
491 }
492
493 // \}
494
498 template<class IdType>
499 const std::vector<std::size_t>& getAdditionalDofDependencies(IdType id, std::size_t stokesElementIdx) const
500 { return emptyStencil_; }
501
505 template<class IdType>
506 const std::vector<std::size_t>& getAdditionalDofDependenciesInverse(IdType id, std::size_t darcyElementIdx) const
507 { return emptyStencil_; }
508
512 bool isCoupledEntity(Dune::index_constant<stokesIdx>, const SubControlVolumeFace<stokesFaceIdx>& scvf) const
513 {
514 return stokesFaceCouplingStencils_.count(scvf.dofIndex());
515 }
516
520 bool isCoupledEntity(Dune::index_constant<darcyIdx>, const SubControlVolumeFace<darcyIdx>& scvf) const
521 {
522 return couplingMapper_.isCoupledDarcyScvf(scvf.index());
523 }
524
525protected:
526
528 std::vector<std::size_t>& emptyStencil()
529 { return emptyStencil_; }
530
531 void removeDuplicates_(std::vector<std::size_t>& stencil)
532 {
533 std::sort(stencil.begin(), stencil.end());
534 stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
535 }
536
537private:
538
539 std::vector<bool> isCoupledDarcyDof_;
540 std::shared_ptr<CouplingData> couplingData_;
541
542 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesCellCenterCouplingStencils_;
543 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesFaceCouplingStencils_;
544 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesCellCenterCouplingStencils_;
545 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesFaceCouplingStencils_;
546 std::vector<std::size_t> emptyStencil_;
547
551 mutable std::vector<StationaryStokesCouplingContext> stokesCouplingContext_;
552 mutable std::vector<StationaryDarcyCouplingContext> darcyCouplingContext_;
553
554 mutable std::size_t boundStokesElemIdx_;
555 mutable std::size_t boundDarcyElemIdx_;
556
557 CouplingMapper couplingMapper_;
558};
559
560} // end namespace Dumux
561
562#endif
The interface of the coupling manager for multi domain problems.
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==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:115
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
A vector of size number equations that can be used for Neumann fluxes, sources, residuals,...
Definition: common/properties.hh:61
The type of the grid view according to the grid type.
Definition: common/properties.hh:63
Traits class encapsulating model specifications.
Definition: common/properties.hh:65
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Stores the boundary types on an element.
Definition: common/properties.hh:112
Definition: common/properties.hh:150
The type for a global container for the volume variables.
Definition: common/properties.hh:176
The global vector of flux variable containers.
Definition: common/properties.hh:186
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:190
Definition: couplingdata.hh:206
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:49
void computeStencils()
Prepare the coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:158
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:193
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &element, 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:434
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)
Update the coupling context for the Darcy residual w.r.t. Darcy DOFs.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:278
const CouplingStencil & couplingStencil(Dune::index_constant< stokesFaceIdx > domainI, const SubControlVolumeFace< stokesIdx > &scvf, Dune::index_constant< darcyIdx > domainJ) const
The coupling stencil of a Stokes face w.r.t. Darcy DOFs.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:482
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:528
const CouplingStencil & couplingStencil(Dune::index_constant< stokesCellCenterIdx > domainI, const Element< stokesIdx > &element, Dune::index_constant< darcyIdx > domainJ) const
The coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:418
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< darcyIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< darcyIdx > &priVars, int pvIdxJ)
Update the coupling context for the Stokes cc residual w.r.t. the Darcy DOFs (FaceToDarcy)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:343
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:499
static constexpr auto darcyIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:59
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:379
StokesDarcyCouplingManager(std::shared_ptr< const GridGeometry< stokesIdx > > stokesFvGridGeometry, std::shared_ptr< const GridGeometry< darcyIdx > > darcyFvGridGeometry)
Constructor.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:124
void update()
Update after the grid has changed.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:148
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:233
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:396
void removeDuplicates_(std::vector< std::size_t > &stencil)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:531
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:134
const CouplingStencil & couplingStencil(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &element, Dune::index_constant< stokesCellCenterIdx > 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:443
const auto & couplingData() const
Access the coupling data.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:371
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesFaceIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< 1 > &priVars, int pvIdxJ)
Update the coupling context for the Darcy residual w.r.t. the Stokes face DOFs (DarcyToFace)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:320
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesCellCenterIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< stokesCellCenterIdx > &priVars, int pvIdxJ)
Update the coupling context for the Darcy residual w.r.t. the Stokes cell-center DOFs (DarcyToCC)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:292
static constexpr auto faceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:57
const CouplingStencil & couplingStencil(Dune::index_constant< darcyIdx > domainI, const Element< darcyIdx > &element, Dune::index_constant< stokesFaceIdx > 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:458
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:512
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:227
static constexpr auto stokesCellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:54
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:506
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:186
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const SubControlVolumeFace< stokesIdx > &scvf, 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:474
static constexpr auto cellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:56
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:520
void updateSolution(const SolutionVector &curSol)
Update the solution vector before assembly.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:154
static constexpr auto stokesFaceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:55
static constexpr auto stokesIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:58
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:49
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:167
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:175
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:183
void computeCouplingMapsAndStencils(Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:94
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:247
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264
SolutionVector & curSol()
the solution vector of the coupled problem
Definition: multidomain/couplingmanager.hh:278
Base coupling manager for the staggered discretization.
Definition: staggeredcouplingmanager.hh:42
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:149
const CouplingStencil & couplingStencil(Dune::index_constant< cellCenterIdx > domainI, const Element &elementI, Dune::index_constant< faceIdx > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: staggeredcouplingmanager.hh:93
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, const std::size_t dofIdxGlobalJ, const PriVarsJ &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: staggeredcouplingmanager.hh:68
Declares all properties used in Dumux.
The local element solution class for staggered methods.