3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/staggeredfreeflow/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 *****************************************************************************/
26#ifndef DUMUX_MULTIDOMAIN_STAGGERED_FREEFLOW_COUPLING_MANAGER_HH
27#define DUMUX_MULTIDOMAIN_STAGGERED_FREEFLOW_COUPLING_MANAGER_HH
28
29#include <memory>
30#include <tuple>
31#include <vector>
32
33#include <dune/common/exceptions.hh>
34#include <dune/common/indices.hh>
35#include <dune/common/float_cmp.hh>
36
42
43namespace Dumux {
44
50template<class Traits>
52: public CouplingManager<Traits>
53{
55public:
56 static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index();
57 static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index();
58
59 // this can be used if the coupling manager is used inside a meta-coupling manager (e.g. multi-binary)
60 // to manager the solution vector storage outside this class
62private:
63 template<std::size_t id> using SubDomainTypeTag = typename Traits::template SubDomain<id>::TypeTag;
64 template<std::size_t id> using PrimaryVariables = GetPropType<SubDomainTypeTag<id>, Properties::PrimaryVariables>;
65 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
66 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
67 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
68 template<std::size_t id> using FVElementGeometry = typename GridGeometry<id>::LocalView;
69 template<std::size_t id> using SubControlVolume = typename FVElementGeometry<id>::SubControlVolume;
70 template<std::size_t id> using SubControlVolumeFace = typename FVElementGeometry<id>::SubControlVolumeFace;
71 template<std::size_t id> using GridVariables = typename Traits::template SubDomain<id>::GridVariables;
72 template<std::size_t id> using ElementVolumeVariables = typename GridVariables<id>::GridVolumeVariables::LocalView;
73 template<std::size_t id> using GridFluxVariablesCache = typename GridVariables<id>::GridFluxVariablesCache;
74 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
75 template<std::size_t id> using VolumeVariables = GetPropType<SubDomainTypeTag<id>, Properties::VolumeVariables>;
76
77 using Scalar = typename Traits::Scalar;
78 using SolutionVector = typename Traits::SolutionVector;
79
80 using CouplingStencilType = std::vector<std::size_t>;
81
82 using GridVariablesTuple = typename Traits::template TupleOfSharedPtr<GridVariables>;
83
84 using FluidSystem = typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
85
86 using VelocityVector = typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
87 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
88
89 struct MomentumCouplingContext
90 {
91 FVElementGeometry<freeFlowMassIndex> fvGeometry;
92 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
93 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
94 std::size_t eIdx;
95 };
96
97 struct MassAndEnergyCouplingContext
98 {
99 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f, const std::size_t i)
100 : fvGeometry(std::move(f))
101 , eIdx(i)
102 {}
103
104 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
105 std::size_t eIdx;
106 };
107
108public:
109
110 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
111
115 // \{
116
118 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
119 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
120 GridVariablesTuple&& gridVariables,
121 const SolutionVector& curSol)
122 {
123 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
124 gridVariables_ = gridVariables;
125 this->updateSolution(curSol);
126
127 computeCouplingStencils_();
128 }
129
131 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
132 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
133 GridVariablesTuple&& gridVariables,
134 const SolutionVector& curSol,
135 const SolutionVector& prevSol)
136 {
137 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables), curSol);
138 prevSol_ = &prevSol;
139 isTransient_ = true;
140 }
141
143 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
144 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
145 GridVariablesTuple&& gridVariables,
147 {
148 this->setSubProblems(std::make_tuple(momentumProblem, massProblem));
149 gridVariables_ = gridVariables;
150 this->attachSolution(curSol);
151
152 computeCouplingStencils_();
153 }
154
155 // \}
156
158
178 template<std::size_t j, class LocalAssemblerI>
179 decltype(auto) evalCouplingResidual(Dune::index_constant<freeFlowMomentumIndex> domainI,
180 const LocalAssemblerI& localAssemblerI,
181 const SubControlVolume<freeFlowMomentumIndex>& scvI,
182 Dune::index_constant<j> domainJ,
183 std::size_t dofIdxGlobalJ) const
184 {
185 const auto& problem = localAssemblerI.problem();
186 const auto& element = localAssemblerI.element();
187 const auto& fvGeometry = localAssemblerI.fvGeometry();
188 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
189 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
190 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
191 const auto& localResidual = localAssemblerI.localResidual();
192
193 localResidual.evalSource(residual, problem, element, fvGeometry, curElemVolVars, scvI);
194
195 for (const auto& scvf : scvfs(fvGeometry, scvI))
196 localResidual.evalFlux(residual, problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
197
198 if (!localAssemblerI.assembler().isStationaryProblem())
199 {
200 assert(isTransient_);
201 localResidual.evalStorage(residual, problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
202 }
203
204 return residual;
205 }
206
210 // \{
211
215 Scalar pressure(const Element<freeFlowMomentumIndex>& element,
216 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
217 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
218 {
219 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
220 return this->curSol(freeFlowMassIndex)[fvGeometry.elementIndex()][pressureIdx];
221 }
222
229 Scalar cellPressure(const Element<freeFlowMassIndex>& element,
230 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
231 {
232 return this->curSol(freeFlowMassIndex)[scvf.insideScvIdx()][pressureIdx];
233 }
234
238 Scalar density(const Element<freeFlowMomentumIndex>& element,
239 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
240 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
241 const bool considerPreviousTimeStep = false) const
242 {
243 assert(!(considerPreviousTimeStep && !isTransient_));
244 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
245 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
246 const auto& insideMassScv = momentumCouplingContext_[0].fvGeometry.scv(insideMomentumScv.elementIndex());
247
248 const auto rho = [&](const auto& elemVolVars)
249 {
250 if (scvf.boundary())
251 return elemVolVars[insideMassScv].density();
252 else
253 {
254 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
255 const auto& outsideMassScv = momentumCouplingContext_[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
256 // TODO distance weighting
257 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
258 }
259 };
260
261 return considerPreviousTimeStep ? rho(momentumCouplingContext_[0].prevElemVolVars)
262 : rho(momentumCouplingContext_[0].curElemVolVars);
263 }
264
265 auto insideAndOutsideDensity(const Element<freeFlowMomentumIndex>& element,
266 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
267 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
268 const bool considerPreviousTimeStep = false) const
269 {
270 assert(!(considerPreviousTimeStep && !isTransient_));
271 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
272 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
273 const auto& insideMassScv = momentumCouplingContext_[0].fvGeometry.scv(insideMomentumScv.elementIndex());
274
275 const auto result = [&](const auto& elemVolVars)
276 {
277 if (scvf.boundary())
278 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[insideMassScv].density());
279 else
280 {
281 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
282 const auto& outsideMassScv = momentumCouplingContext_[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
283 return std::make_pair(elemVolVars[insideMassScv].density(), elemVolVars[outsideMassScv].density());
284 }
285 };
286
287 return considerPreviousTimeStep ? result(momentumCouplingContext_[0].prevElemVolVars)
288 : result(momentumCouplingContext_[0].curElemVolVars);
289 }
290
294 Scalar density(const Element<freeFlowMomentumIndex>& element,
295 const SubControlVolume<freeFlowMomentumIndex>& scv,
296 const bool considerPreviousTimeStep = false) const
297 {
298 assert(!(considerPreviousTimeStep && !isTransient_));
299 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
300 const auto& massScv = (*scvs(momentumCouplingContext_[0].fvGeometry).begin());
301
302 return considerPreviousTimeStep ? momentumCouplingContext_[0].prevElemVolVars[massScv].density()
303 : momentumCouplingContext_[0].curElemVolVars[massScv].density();
304 }
305
309 Scalar effectiveViscosity(const Element<freeFlowMomentumIndex>& element,
310 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
311 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf) const
312 {
313 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
314
315 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
316 const auto& insideMassScv = momentumCouplingContext_[0].fvGeometry.scv(insideMomentumScv.elementIndex());
317
318 if (scvf.boundary())
319 return momentumCouplingContext_[0].curElemVolVars[insideMassScv].viscosity();
320
321 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
322 const auto& outsideMassScv = momentumCouplingContext_[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
323
324 const auto mu = [&](const auto& elemVolVars)
325 {
326 // TODO distance weighting
327 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
328 };
329
330 return mu(momentumCouplingContext_[0].curElemVolVars);
331 }
332
336 VelocityVector faceVelocity(const Element<freeFlowMassIndex>& element,
337 const SubControlVolumeFace<freeFlowMassIndex>& scvf) const
338 {
339 // TODO: rethink this! Maybe we only need scvJ.dofIndex()
340 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx()/*eIdx*/);
341
342 // the TPFA scvf index corresponds to the staggered scv index (might need mapping)
343 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_[0].fvGeometry);
344 const auto& scvJ = massAndEnergyCouplingContext_[0].fvGeometry.scv(localMomentumScvIdx);
345
346 // create a unit normal vector oriented in positive coordinate direction
347 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
348 velocity[scvJ.dofAxis()] = 1.0;
349
350 // create the actual velocity vector
351 velocity *= this->curSol(freeFlowMomentumIndex)[scvJ.dofIndex()];
352
353 return velocity;
354 }
355
360 template<std::size_t j>
361 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
362 const Element<freeFlowMomentumIndex>& elementI,
363 const SubControlVolume<freeFlowMomentumIndex>& scvI,
364 Dune::index_constant<j> domainJ) const
365 { return emptyStencil_; }
366
381 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMassIndex> domainI,
382 const Element<freeFlowMassIndex>& elementI,
383 Dune::index_constant<freeFlowMomentumIndex> domainJ) const
384 {
385 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
386 return massAndEnergyToMomentumStencils_[eIdx];
387 }
388
398 const CouplingStencilType& couplingStencil(Dune::index_constant<freeFlowMomentumIndex> domainI,
399 const Element<freeFlowMomentumIndex>& elementI,
400 const SubControlVolume<freeFlowMomentumIndex>& scvI,
401 Dune::index_constant<freeFlowMassIndex> domainJ) const
402 {
403 return momentumToMassAndEnergyStencils_[scvI.index()];
404 }
405
406 // \}
407
411 // \{
412
433 template<std::size_t i, std::size_t j, class LocalAssemblerI>
434 void updateCouplingContext(Dune::index_constant<i> domainI,
435 const LocalAssemblerI& localAssemblerI,
436 Dune::index_constant<j> domainJ,
437 std::size_t dofIdxGlobalJ,
438 const PrimaryVariables<j>& priVarsJ,
439 int pvIdxJ)
440 {
441 this->curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
442
443 if constexpr (domainI == freeFlowMomentumIndex && domainJ == freeFlowMassIndex)
444 {
445 bindCouplingContext_(domainI, localAssemblerI.element());
446
447 const auto& problem = this->problem(domainJ);
448 const auto& deflectedElement = problem.gridGeometry().element(dofIdxGlobalJ);
449 const auto elemSol = elementSolution(deflectedElement, this->curSol(domainJ), problem.gridGeometry());
450 const auto& fvGeometry = momentumCouplingContext_[0].fvGeometry;
451 const auto scvIdxJ = dofIdxGlobalJ;
452 const auto& scv = fvGeometry.scv(scvIdxJ);
453
454 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
455 gridVars_(freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol), problem, deflectedElement, scv);
456 else
457 momentumCouplingContext_[0].curElemVolVars[scv].update(std::move(elemSol), problem, deflectedElement, scv);
458 }
459 }
460
461 // \}
462
463private:
464 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
465 const Element<freeFlowMomentumIndex>& elementI) const
466 {
467 const auto eIdx = this->problem(freeFlowMomentumIndex).gridGeometry().elementMapper().index(elementI);
468 bindCouplingContext_(domainI, elementI, eIdx);
469 }
470
471 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
472 const Element<freeFlowMomentumIndex>& elementI,
473 const std::size_t eIdx) const
474 {
475 if (momentumCouplingContext_.empty())
476 {
477 auto fvGeometry = localView(this->problem(freeFlowMassIndex).gridGeometry());
478 fvGeometry.bind(elementI);
479
480 auto curElemVolVars = localView(gridVars_(freeFlowMassIndex).curGridVolVars());
481 curElemVolVars.bind(elementI, fvGeometry, this->curSol(freeFlowMassIndex));
482
483 auto prevElemVolVars = isTransient_ ? localView(gridVars_(freeFlowMassIndex).prevGridVolVars())
484 : localView(gridVars_(freeFlowMassIndex).curGridVolVars());
485
486 if (isTransient_)
487 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[freeFlowMassIndex]);
488
489 momentumCouplingContext_.emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
490 }
491 else if (eIdx != momentumCouplingContext_[0].eIdx)
492 {
493 momentumCouplingContext_[0].eIdx = eIdx;
494 momentumCouplingContext_[0].fvGeometry.bind(elementI);
495 momentumCouplingContext_[0].curElemVolVars.bind(elementI, momentumCouplingContext_[0].fvGeometry, this->curSol(freeFlowMassIndex));
496
497 if (isTransient_)
498 momentumCouplingContext_[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_[0].fvGeometry, (*prevSol_)[freeFlowMassIndex]);
499 }
500 }
501
502 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
503 const Element<freeFlowMassIndex>& elementI) const
504 {
505 const auto eIdx = this->problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
506 bindCouplingContext_(domainI, elementI, eIdx);
507 }
508
509 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
510 const Element<freeFlowMassIndex>& elementI,
511 const std::size_t eIdx) const
512 {
513 if (massAndEnergyCouplingContext_.empty())
514 {
515 const auto& gridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
516 auto fvGeometry = localView(gridGeometry);
517 fvGeometry.bindElement(elementI);
518 massAndEnergyCouplingContext_.emplace_back(std::move(fvGeometry), eIdx);
519 }
520 else if (eIdx != massAndEnergyCouplingContext_[0].eIdx)
521 {
522 massAndEnergyCouplingContext_[0].eIdx = eIdx;
523 massAndEnergyCouplingContext_[0].fvGeometry.bindElement(elementI);
524 }
525 }
526
531 template<std::size_t i>
532 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx) const
533 {
534 if (std::get<i>(gridVariables_))
535 return *std::get<i>(gridVariables_);
536 else
537 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
538 }
539
544 template<std::size_t i>
545 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
546 {
547 if (std::get<i>(gridVariables_))
548 return *std::get<i>(gridVariables_);
549 else
550 DUNE_THROW(Dune::InvalidStateException, "The gridVariables pointer was not set. Use setGridVariables() before calling this function");
551 }
552
553
554 void computeCouplingStencils_()
555 {
556 // TODO higher order
557 const auto& momentumGridGeometry = this->problem(freeFlowMomentumIndex).gridGeometry();
558 auto momentumFvGeometry = localView(momentumGridGeometry);
559 massAndEnergyToMomentumStencils_.clear();
560 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
561
562 momentumToMassAndEnergyStencils_.clear();
563 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
564
565 for (const auto& element : elements(momentumGridGeometry.gridView()))
566 {
567 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
568 momentumFvGeometry.bind(element);
569 for (const auto& scv : scvs(momentumFvGeometry))
570 {
571 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
572 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
573
574 // extend the stencil for fluids with variable viscosity and density,
575 if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/))
576 // if constexpr (FluidSystem::isCompressible(0/*phaseIdx*/) || !FluidSystem::viscosityIsConstant(0/*phaseIdx*/)) // TODO fix on master
577 {
578 for (const auto& scvf : scvfs(momentumFvGeometry, scv))
579 {
580 if (scvf.isLateral() && !scvf.boundary())
581 {
582 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
583 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
584 }
585 }
586 }
587 }
588 }
589 }
590
591 std::size_t massScvfToMomentumScvIdx_(const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
592 [[maybe_unused]] const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry) const
593 {
594 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
595 return massScvf.index();
596 else
597 {
598 static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true);
599 if (!makeConsistentlyOriented)
600 return massScvf.index();
601
602 for (const auto& momentumScv : scvs(momentumFVGeometry))
603 {
604 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
605 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
606 if (Dune::FloatCmp::eq<typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
607 return momentumScv.index();
608 }
609 DUNE_THROW(Dune::InvalidStateException, "No Momentum SCV found");
610 }
611 }
612
613 CouplingStencilType emptyStencil_;
614 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
615 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
616
617 mutable std::vector<MomentumCouplingContext> momentumCouplingContext_;
618 mutable std::vector<MassAndEnergyCouplingContext> massAndEnergyCouplingContext_;
619
621 GridVariablesTuple gridVariables_;
622
623 const SolutionVector* prevSol_;
624 bool isTransient_;
625};
626
627} // end namespace Dumux
628
629#endif
An adapter class for local assemblers using numeric differentiation.
Type traits.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
decltype(auto) evalCouplingResidual(Dune::index_constant< freeFlowMomentumIndex > domainI, const LocalAssemblerI &localAssemblerI, const SubControlVolume< freeFlowMomentumIndex > &scvI, 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/staggeredfreeflow/couplingmanager.hh:179
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)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:434
decltype(auto) evalCouplingResidual(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluates the element residual of a coupled element of domain i which depends on the variables at the...
Definition: multidomain/couplingmanager.hh:261
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
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:57
Definition: common/properties.hh:102
The secondary variables within a sub-control volume.
Definition: common/properties.hh:107
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:334
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
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:76
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
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:83
The interface of the coupling manager for free flow systems.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:53
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given frontal sub control volume face.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:215
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
Returns the pressure at a given sub control volume face.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:309
static constexpr auto freeFlowMomentumIndex
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:56
Scalar density(const Element< freeFlowMomentumIndex > &element, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:294
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume face.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:238
Scalar cellPressure(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the pressure at the center of a sub control volume corresponding to a given sub control volum...
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:229
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMassIndex > domainI, const Element< freeFlowMassIndex > &elementI, Dune::index_constant< freeFlowMomentumIndex > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:381
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:336
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:265
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< freeFlowMassIndex > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:398
static constexpr auto freeFlowMassIndex
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:57
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, typename ParentType::SolutionVectorStorage &curSol)
use as binary coupling manager in multi model context
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:143
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, 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/staggeredfreeflow/couplingmanager.hh:361
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol)
Methods to be accessed by main.
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:118
static constexpr auto pressureIdx
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:110
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol, const SolutionVector &prevSol)
use as regular coupling manager in a transient setting
Definition: multidomain/staggeredfreeflow/couplingmanager.hh:131
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.