3.5-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>
37
38#include "couplingdata.hh"
39#include "couplingmapper.hh"
40
41namespace Dumux {
42
47template<class MDTraits>
49: public StaggeredCouplingManager<MDTraits>
50{
51 using Scalar = typename MDTraits::Scalar;
53
54public:
55 static constexpr auto stokesFaceIdx = typename MDTraits::template SubDomain<0>::Index();
56 static constexpr auto stokesCellCenterIdx = typename MDTraits::template SubDomain<1>::Index();
57 static constexpr auto stokesIdx = stokesFaceIdx;
58 static constexpr auto darcyIdx = typename MDTraits::template SubDomain<2>::Index();
59
60private:
61
62 using SolutionVector = typename MDTraits::SolutionVector;
63
64 // obtain the type tags of the sub problems
65 using StokesTypeTag = typename MDTraits::template SubDomain<0>::TypeTag;
66 using DarcyTypeTag = typename MDTraits::template SubDomain<2>::TypeTag;
67
68 using CouplingStencils = std::unordered_map<std::size_t, std::vector<std::size_t> >;
69 using CouplingStencil = CouplingStencils::mapped_type;
70
71 // the sub domain type tags
72 template<std::size_t id>
73 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
74
75 static constexpr bool isCompositional = GetPropType<SubDomainTypeTag<0>, Properties::ModelTraits>::numFluidComponents()> 1;
76
77 template<std::size_t id> using GridView = typename GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>::GridView;
78 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
79 template<std::size_t id> using PrimaryVariables = typename MDTraits::template SubDomain<id>::PrimaryVariables;
80 template<std::size_t id> using NumEqVector = Dumux::NumEqVector<PrimaryVariables<id>>;
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 SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
91
93
94 using VelocityVector = typename Element<stokesIdx>::Geometry::GlobalCoordinate;
95
97
98 struct StationaryStokesCouplingContext
99 {
100 Element<darcyIdx> element;
101 FVElementGeometry<darcyIdx> fvGeometry;
102 std::size_t darcyScvfIdx;
103 std::size_t stokesScvfIdx;
104 VolumeVariables<darcyIdx> volVars;
105 };
106
107 struct StationaryDarcyCouplingContext
108 {
109 Element<stokesIdx> element;
110 FVElementGeometry<stokesIdx> fvGeometry;
111 std::size_t stokesScvfIdx;
112 std::size_t darcyScvfIdx;
113 VelocityVector velocity;
114 VolumeVariables<stokesIdx> volVars;
115 };
116public:
117
121
123 StokesDarcyCouplingManager(std::shared_ptr<const GridGeometry<stokesIdx>> stokesFvGridGeometry,
124 std::shared_ptr<const GridGeometry<darcyIdx>> darcyFvGridGeometry)
125 { }
126
130 // \{
131
133 void init(std::shared_ptr<const Problem<stokesIdx>> stokesProblem,
134 std::shared_ptr<const Problem<darcyIdx>> darcyProblem,
135 const SolutionVector& curSol)
136 {
137 if (Dune::FloatCmp::ne(stokesProblem->gravity(), darcyProblem->spatialParams().gravity({})))
138 DUNE_THROW(Dune::InvalidStateException, "Both models must use the same gravity vector");
139
140 this->setSubProblems(std::make_tuple(stokesProblem, stokesProblem, darcyProblem));
141 this->updateSolution(curSol);
142 couplingData_ = std::make_shared<CouplingData>(*this);
144 }
145
146 // \}
147
150 {
151 couplingMapper_.computeCouplingMapsAndStencils(*this,
152 darcyToStokesCellCenterCouplingStencils_,
153 darcyToStokesFaceCouplingStencils_,
154 stokesCellCenterCouplingStencils_,
155 stokesFaceCouplingStencils_);
156
157 for(auto&& stencil : darcyToStokesCellCenterCouplingStencils_)
158 removeDuplicates_(stencil.second);
159 for(auto&& stencil : darcyToStokesFaceCouplingStencils_)
160 removeDuplicates_(stencil.second);
161 for(auto&& stencil : stokesCellCenterCouplingStencils_)
162 removeDuplicates_(stencil.second);
163 for(auto&& stencil : stokesFaceCouplingStencils_)
164 removeDuplicates_(stencil.second);
165 }
166
170 // \{
171
173
177 template<std::size_t i, class Assembler, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
178 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<stokesCellCenterIdx>& element, const Assembler& assembler) const
179 { bindCouplingContext(domainI, element); }
180
184 template<std::size_t i, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
185 void bindCouplingContext(Dune::index_constant<i> domainI, const Element<stokesCellCenterIdx>& element) const
186 {
187 stokesCouplingContext_.clear();
188
189 const auto stokesElementIdx = this->problem(stokesIdx).gridGeometry().elementMapper().index(element);
190 boundStokesElemIdx_ = stokesElementIdx;
191
192 // do nothing if the element is not coupled to the other domain
193 if(!couplingMapper_.stokesElementToDarcyElementMap().count(stokesElementIdx))
194 return;
195
196 // prepare the coupling context
197 const auto& darcyIndices = couplingMapper_.stokesElementToDarcyElementMap().at(stokesElementIdx);
198 auto darcyFvGeometry = localView(this->problem(darcyIdx).gridGeometry());
199 for(auto&& indices : darcyIndices)
200 {
201 const auto& darcyElement = this->problem(darcyIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
202 darcyFvGeometry.bindElement(darcyElement);
203 const auto& scv = (*scvs(darcyFvGeometry).begin());
204
205 const auto darcyElemSol = elementSolution(darcyElement, this->curSol(darcyIdx), this->problem(darcyIdx).gridGeometry());
206 VolumeVariables<darcyIdx> darcyVolVars;
207 darcyVolVars.update(darcyElemSol, this->problem(darcyIdx), darcyElement, scv);
208
209 // add the context
210 stokesCouplingContext_.push_back({darcyElement, darcyFvGeometry, indices.scvfIdx, indices.flipScvfIdx, darcyVolVars});
211 }
212 }
213
217 template<class Assembler>
218 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI, const Element<darcyIdx>& element, const Assembler& assembler) const
219 { bindCouplingContext(domainI, element); }
220
224 void bindCouplingContext(Dune::index_constant<darcyIdx> domainI, const Element<darcyIdx>& element) const
225 {
226 darcyCouplingContext_.clear();
227
228 const auto darcyElementIdx = this->problem(darcyIdx).gridGeometry().elementMapper().index(element);
229 boundDarcyElemIdx_ = darcyElementIdx;
230
231 // do nothing if the element is not coupled to the other domain
232 if(!couplingMapper_.darcyElementToStokesElementMap().count(darcyElementIdx))
233 return;
234
235 // prepare the coupling context
236 const auto& stokesElementIndices = couplingMapper_.darcyElementToStokesElementMap().at(darcyElementIdx);
237 auto stokesFvGeometry = localView(this->problem(stokesIdx).gridGeometry());
238
239 for(auto&& indices : stokesElementIndices)
240 {
241 const auto& stokesElement = this->problem(stokesIdx).gridGeometry().boundingBoxTree().entitySet().entity(indices.eIdx);
242 stokesFvGeometry.bindElement(stokesElement);
243
244 VelocityVector faceVelocity(0.0);
245
246 for(const auto& scvf : scvfs(stokesFvGeometry))
247 {
248 if(scvf.index() == indices.scvfIdx)
249 faceVelocity[scvf.directionIndex()] = this->curSol(stokesFaceIdx)[scvf.dofIndex()];
250 }
251
252 using PriVarsType = typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
253 const auto& cellCenterPriVars = this->curSol(stokesCellCenterIdx)[indices.eIdx];
254 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(cellCenterPriVars);
255
256 VolumeVariables<stokesIdx> stokesVolVars;
257 for(const auto& scv : scvs(stokesFvGeometry))
258 stokesVolVars.update(elemSol, this->problem(stokesIdx), stokesElement, scv);
259
260 // add the context
261 darcyCouplingContext_.push_back({stokesElement, stokesFvGeometry, indices.scvfIdx, indices.flipScvfIdx, faceVelocity, stokesVolVars});
262 }
263 }
264
268 template<class LocalAssemblerI>
269 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
270 const LocalAssemblerI& localAssemblerI,
271 Dune::index_constant<darcyIdx> domainJ,
272 std::size_t dofIdxGlobalJ,
273 const PrimaryVariables<darcyIdx>& priVarsJ,
274 int pvIdxJ)
275 {
276 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
277 }
278
282 template<class LocalAssemblerI>
283 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
284 const LocalAssemblerI& localAssemblerI,
285 Dune::index_constant<stokesCellCenterIdx> domainJ,
286 const std::size_t dofIdxGlobalJ,
287 const PrimaryVariables<stokesCellCenterIdx>& priVars,
288 int pvIdxJ)
289 {
290 this->curSol(domainJ)[dofIdxGlobalJ] = priVars;
291
292 for (auto& data : darcyCouplingContext_)
293 {
294 const auto stokesElemIdx = this->problem(stokesIdx).gridGeometry().elementMapper().index(data.element);
295
296 if(stokesElemIdx != dofIdxGlobalJ)
297 continue;
298
299 using PriVarsType = typename VolumeVariables<stokesCellCenterIdx>::PrimaryVariables;
300 const auto elemSol = makeElementSolutionFromCellCenterPrivars<PriVarsType>(priVars);
301
302 for(const auto& scv : scvs(data.fvGeometry))
303 data.volVars.update(elemSol, this->problem(stokesIdx), data.element, scv);
304 }
305 }
306
310 template<class LocalAssemblerI>
311 void updateCouplingContext(Dune::index_constant<darcyIdx> domainI,
312 const LocalAssemblerI& localAssemblerI,
313 Dune::index_constant<stokesFaceIdx> domainJ,
314 const std::size_t dofIdxGlobalJ,
315 const PrimaryVariables<stokesFaceIdx>& priVars,
316 int pvIdxJ)
317 {
318 this->curSol(domainJ)[dofIdxGlobalJ] = priVars;
319
320 for (auto& data : darcyCouplingContext_)
321 {
322 for(const auto& scvf : scvfs(data.fvGeometry))
323 {
324 if(scvf.dofIndex() == dofIdxGlobalJ)
325 data.velocity[scvf.directionIndex()] = priVars;
326 }
327 }
328 }
329
333 template<std::size_t i, class LocalAssemblerI, std::enable_if_t<(i == stokesCellCenterIdx || i == stokesFaceIdx), int> = 0>
334 void updateCouplingContext(Dune::index_constant<i> domainI,
335 const LocalAssemblerI& localAssemblerI,
336 Dune::index_constant<darcyIdx> domainJ,
337 const std::size_t dofIdxGlobalJ,
338 const PrimaryVariables<darcyIdx>& priVars,
339 int pvIdxJ)
340 {
341 this->curSol(domainJ)[dofIdxGlobalJ] = priVars;
342
343 for (auto& data : stokesCouplingContext_)
344 {
345 const auto darcyElemIdx = this->problem(darcyIdx).gridGeometry().elementMapper().index(data.element);
346
347 if(darcyElemIdx != dofIdxGlobalJ)
348 continue;
349
350 const auto darcyElemSol = elementSolution(data.element, this->curSol(darcyIdx), this->problem(darcyIdx).gridGeometry());
351
352 for(const auto& scv : scvs(data.fvGeometry))
353 data.volVars.update(darcyElemSol, this->problem(darcyIdx), data.element, scv);
354 }
355 }
356
357 // \}
358
362 const auto& couplingData() const
363 {
364 return *couplingData_;
365 }
366
370 const auto& stokesCouplingContext(const Element<stokesIdx>& element, const SubControlVolumeFace<stokesIdx>& scvf) const
371 {
372 if (stokesCouplingContext_.empty() || boundStokesElemIdx_ != scvf.insideScvIdx())
374
375 for(const auto& context : stokesCouplingContext_)
376 {
377 if(scvf.index() == context.stokesScvfIdx)
378 return context;
379 }
380
381 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
382 }
383
387 const auto& darcyCouplingContext(const Element<darcyIdx>& element, const SubControlVolumeFace<darcyIdx>& scvf) const
388 {
389 if (darcyCouplingContext_.empty() || boundDarcyElemIdx_ != scvf.insideScvIdx())
391
392 for(const auto& context : darcyCouplingContext_)
393 {
394 if(scvf.index() == context.darcyScvfIdx)
395 return context;
396 }
397
398 DUNE_THROW(Dune::InvalidStateException, "No coupling context found at scvf " << scvf.center());
399 }
400
404 // \{
405
409 const CouplingStencil& couplingStencil(Dune::index_constant<stokesCellCenterIdx> domainI,
410 const Element<stokesIdx>& element,
411 Dune::index_constant<darcyIdx> domainJ) const
412 {
413 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
414 if(stokesCellCenterCouplingStencils_.count(eIdx))
415 return stokesCellCenterCouplingStencils_.at(eIdx);
416 else
417 return emptyStencil_;
418 }
419
424 template<std::size_t i, std::size_t j>
425 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
426 const Element<i>& element,
427 Dune::index_constant<j> domainJ) const
428 { return emptyStencil_; }
429
434 const CouplingStencil& couplingStencil(Dune::index_constant<darcyIdx> domainI,
435 const Element<darcyIdx>& element,
436 Dune::index_constant<stokesCellCenterIdx> domainJ) const
437 {
438 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
439 if(darcyToStokesCellCenterCouplingStencils_.count(eIdx))
440 return darcyToStokesCellCenterCouplingStencils_.at(eIdx);
441 else
442 return emptyStencil_;
443 }
444
449 const CouplingStencil& couplingStencil(Dune::index_constant<darcyIdx> domainI,
450 const Element<darcyIdx>& element,
451 Dune::index_constant<stokesFaceIdx> domainJ) const
452 {
453 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
454 if (darcyToStokesFaceCouplingStencils_.count(eIdx))
455 return darcyToStokesFaceCouplingStencils_.at(eIdx);
456 else
457 return emptyStencil_;
458 }
459
464 template<std::size_t i, std::size_t j>
465 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
466 const SubControlVolumeFace<stokesIdx>& scvf,
467 Dune::index_constant<j> domainJ) const
468 { return emptyStencil_; }
469
473 const CouplingStencil& couplingStencil(Dune::index_constant<stokesFaceIdx> domainI,
474 const SubControlVolumeFace<stokesIdx>& scvf,
475 Dune::index_constant<darcyIdx> domainJ) const
476 {
477 const auto faceDofIdx = scvf.dofIndex();
478 if(stokesFaceCouplingStencils_.count(faceDofIdx))
479 return stokesFaceCouplingStencils_.at(faceDofIdx);
480 else
481 return emptyStencil_;
482 }
483
484 // \}
485
489 template<class IdType>
490 const std::vector<std::size_t>& getAdditionalDofDependencies(IdType id, std::size_t stokesElementIdx) const
491 { return emptyStencil_; }
492
496 template<class IdType>
497 const std::vector<std::size_t>& getAdditionalDofDependenciesInverse(IdType id, std::size_t darcyElementIdx) const
498 { return emptyStencil_; }
499
503 bool isCoupledEntity(Dune::index_constant<stokesIdx>, const SubControlVolumeFace<stokesFaceIdx>& scvf) const
504 {
505 return stokesFaceCouplingStencils_.count(scvf.dofIndex());
506 }
507
511 bool isCoupledEntity(Dune::index_constant<darcyIdx>, const SubControlVolumeFace<darcyIdx>& scvf) const
512 {
513 return couplingMapper_.isCoupledDarcyScvf(scvf.index());
514 }
515
516protected:
517
519 std::vector<std::size_t>& emptyStencil()
520 { return emptyStencil_; }
521
522 void removeDuplicates_(std::vector<std::size_t>& stencil)
523 {
524 std::sort(stencil.begin(), stencil.end());
525 stencil.erase(std::unique(stencil.begin(), stencil.end()), stencil.end());
526 }
527
528private:
529
530 std::vector<bool> isCoupledDarcyDof_;
531 std::shared_ptr<CouplingData> couplingData_;
532
533 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesCellCenterCouplingStencils_;
534 std::unordered_map<std::size_t, std::vector<std::size_t> > stokesFaceCouplingStencils_;
535 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesCellCenterCouplingStencils_;
536 std::unordered_map<std::size_t, std::vector<std::size_t> > darcyToStokesFaceCouplingStencils_;
537 std::vector<std::size_t> emptyStencil_;
538
542 mutable std::vector<StationaryStokesCouplingContext> stokesCouplingContext_;
543 mutable std::vector<StationaryDarcyCouplingContext> darcyCouplingContext_;
544
545 mutable std::size_t boundStokesElemIdx_;
546 mutable std::size_t boundDarcyElemIdx_;
547
548 CouplingMapper couplingMapper_;
549};
550
551} // end namespace Dumux
552
553#endif
A helper to deduce a vector with the same size as numbers of equations.
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==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
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:46
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Traits class encapsulating model specifications.
Definition: common/properties.hh:53
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:57
Stores the boundary types on an element.
Definition: common/properties.hh:99
Definition: common/properties.hh:102
The type for a global container for the volume variables.
Definition: common/properties.hh:109
The global vector of flux variable containers.
Definition: common/properties.hh:119
The grid variables object managing variable data on the grid (volvars/fluxvars cache)
Definition: common/properties.hh:123
Definition: couplingdata.hh:206
Coupling manager for Stokes and Darcy domains with equal dimension.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:50
void computeStencils()
Prepare the coupling stencils.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:149
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:185
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:425
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:269
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:473
std::vector< std::size_t > & emptyStencil()
Return a reference to an empty stencil.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:519
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:409
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:334
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:490
static constexpr auto darcyIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:58
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:370
StokesDarcyCouplingManager(std::shared_ptr< const GridGeometry< stokesIdx > > stokesFvGridGeometry, std::shared_ptr< const GridGeometry< darcyIdx > > darcyFvGridGeometry)
Constructor.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:123
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:224
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:387
void removeDuplicates_(std::vector< std::size_t > &stencil)
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:522
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:133
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:434
const auto & couplingData() const
Access the coupling data.
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:362
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:283
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:449
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:503
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:218
static constexpr auto stokesCellCenterIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:56
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:497
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:178
void updateCouplingContext(Dune::index_constant< darcyIdx > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< stokesFaceIdx > domainJ, const std::size_t dofIdxGlobalJ, const PrimaryVariables< stokesFaceIdx > &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:311
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:465
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:511
static constexpr auto stokesFaceIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:55
static constexpr auto stokesIdx
Definition: multidomain/boundary/stokesdarcy/couplingmanager.hh:57
Coupling mapper for Stokes and Darcy domains with equal dimension.
Definition: boundary/stokesdarcy/couplingmapper.hh:42
bool isCoupledDarcyScvf(std::size_t darcyScvfIdx) const
Returns whether a Darcy scvf is coupled to the other domain.
Definition: boundary/stokesdarcy/couplingmapper.hh:131
void computeCouplingMapsAndStencils(const CouplingManager &couplingManager, Stencils &darcyToStokesCellCenterStencils, Stencils &darcyToStokesFaceStencils, Stencils &stokesCellCenterToDarcyStencils, Stencils &stokesFaceToDarcyStencils)
Main update routine.
Definition: boundary/stokesdarcy/couplingmapper.hh:56
const auto & darcyElementToStokesElementMap() const
A map that returns all Stokes elements coupled to a Darcy element.
Definition: boundary/stokesdarcy/couplingmapper.hh:139
const auto & stokesElementToDarcyElementMap() const
A map that returns all Darcy elements coupled to a Stokes element.
Definition: boundary/stokesdarcy/couplingmapper.hh:147
decltype(auto) curSol()
the solution vector of the coupled problem
Definition: multidomain/couplingmanager.hh:370
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
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:231
Base coupling manager for the staggered discretization.
Definition: staggeredcouplingmanager.hh:43
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:150
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:94
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:69
Declares all properties used in Dumux.
The local element solution class for staggered methods.