12#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_HH
13#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_HH
20#include <dune/common/exceptions.hh>
21#include <dune/common/indices.hh>
22#include <dune/common/float_cmp.hh>
58 template<std::
size_t id>
using SubDomainTypeTag =
typename Traits::template SubDomain<id>::TypeTag;
61 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
62 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
63 template<std::
size_t id>
using ElementSeed =
typename GridView<id>::Grid::template Codim<0>::EntitySeed;
64 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
65 template<std::
size_t id>
using SubControlVolume =
typename FVElementGeometry<id>::SubControlVolume;
66 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
67 template<std::
size_t id>
using GridVariables =
typename Traits::template SubDomain<id>::GridVariables;
68 template<std::
size_t id>
using ElementVolumeVariables =
typename GridVariables<id>::GridVolumeVariables::LocalView;
69 template<std::
size_t id>
using GridFluxVariablesCache =
typename GridVariables<id>::GridFluxVariablesCache;
73 using Scalar =
typename Traits::Scalar;
76 template<std::
size_t id>
77 using SubSolutionVector
78 = std::decay_t<decltype(std::declval<SolutionVector>()[Dune::index_constant<id>()])>;
80 template<std::
size_t id>
81 using ConstSubSolutionVectorPtr =
const SubSolutionVector<id>*;
83 using PrevSolutionVectorStorage =
typename Traits::template Tuple<ConstSubSolutionVectorPtr>;
87 using GridVariablesTuple =
typename Traits::template TupleOfSharedPtr<GridVariables>;
89 using FluidSystem =
typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
91 using VelocityVector =
typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
92 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
94 struct MomentumCouplingContext
96 FVElementGeometry<freeFlowMassIndex> fvGeometry;
97 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
98 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
102 struct MassAndEnergyCouplingContext
104 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f,
const std::size_t i)
105 : fvGeometry(std::move(f))
109 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
115 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
123 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
124 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
125 GridVariablesTuple&& gridVariables,
128 this->momentumCouplingContext_().clear();
129 this->massAndEnergyCouplingContext_().clear();
131 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
132 gridVariables_ = gridVariables;
135 computeCouplingStencils_();
139 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
140 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
141 GridVariablesTuple&& gridVariables,
145 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables),
curSol);
147 Dune::Hybrid::forEach(std::make_index_sequence<Traits::numSubDomains>{}, [&](
auto i)
148 { std::get<i>(prevSolutions_) = &prevSol[i]; });
152 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
153 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
154 GridVariablesTuple&& gridVariables,
157 this->momentumCouplingContext_().clear();
158 this->massAndEnergyCouplingContext_().clear();
160 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
161 gridVariables_ = gridVariables;
164 computeCouplingStencils_();
168 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
169 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
170 GridVariablesTuple&& gridVariables,
172 const PrevSolutionVectorStorage& prevSol)
174 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables),
curSol);
175 prevSolutions_ = prevSol;
199 template<std::
size_t j,
class LocalAssemblerI>
201 const LocalAssemblerI& localAssemblerI,
202 const SubControlVolume<freeFlowMomentumIndex>& scvI,
203 Dune::index_constant<j> domainJ,
204 std::size_t dofIdxGlobalJ)
const
206 const auto&
problem = localAssemblerI.problem();
207 const auto& element = localAssemblerI.element();
208 const auto& fvGeometry = localAssemblerI.fvGeometry();
209 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
210 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
211 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
212 const auto& localResidual = localAssemblerI.localResidual();
214 localResidual.evalSource(residual,
problem, element, fvGeometry, curElemVolVars, scvI);
216 for (
const auto& scvf : scvfs(fvGeometry, scvI))
217 localResidual.evalFlux(residual,
problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
219 if (!localAssemblerI.assembler().isStationaryProblem())
221 assert(isTransient_());
222 localResidual.evalStorage(residual,
problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
236 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
237 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
238 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
240 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
251 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
259 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
260 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
261 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
262 const bool considerPreviousTimeStep =
false)
const
264 assert(!(considerPreviousTimeStep && !isTransient_()));
265 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
266 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
267 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
269 const auto rho = [&](
const auto& elemVolVars)
272 return elemVolVars[insideMassScv].
density();
275 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
276 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
278 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
282 return considerPreviousTimeStep ? rho(momentumCouplingContext_()[0].prevElemVolVars)
283 : rho(momentumCouplingContext_()[0].curElemVolVars);
287 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
288 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
289 const bool considerPreviousTimeStep =
false)
const
291 assert(!(considerPreviousTimeStep && !isTransient_()));
292 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
293 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
294 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
296 const auto result = [&](
const auto& elemVolVars)
299 return std::make_pair(elemVolVars[insideMassScv].
density(), elemVolVars[insideMassScv].density());
302 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
303 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
304 return std::make_pair(elemVolVars[insideMassScv].
density(), elemVolVars[outsideMassScv].
density());
308 return considerPreviousTimeStep ? result(momentumCouplingContext_()[0].prevElemVolVars)
309 : result(momentumCouplingContext_()[0].curElemVolVars);
315 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
316 const SubControlVolume<freeFlowMomentumIndex>& scv,
317 const bool considerPreviousTimeStep =
false)
const
319 assert(!(considerPreviousTimeStep && !isTransient_()));
320 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
321 const auto& massScv = (*
scvs(momentumCouplingContext_()[0].fvGeometry).begin());
323 return considerPreviousTimeStep ? momentumCouplingContext_()[0].prevElemVolVars[massScv].density()
324 : momentumCouplingContext_()[0].curElemVolVars[massScv].density();
331 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
332 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
334 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
336 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
337 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
340 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
342 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
343 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
345 const auto mu = [&](
const auto& elemVolVars)
348 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
351 return mu(momentumCouplingContext_()[0].curElemVolVars);
358 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
359 const SubControlVolume<freeFlowMomentumIndex>& scv)
const
361 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
362 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(scv.elementIndex());
363 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
370 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
373 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx());
376 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_()[0].fvGeometry);
377 const auto& scvJ = massAndEnergyCouplingContext_()[0].fvGeometry.scv(localMomentumScvIdx);
380 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
381 velocity[scvJ.dofAxis()] = 1.0;
398 template<std::
size_t j>
400 const Element<freeFlowMomentumIndex>& elementI,
401 const SubControlVolume<freeFlowMomentumIndex>& scvI,
402 Dune::index_constant<j> domainJ)
const
403 {
return emptyStencil_; }
420 const Element<freeFlowMassIndex>& elementI,
421 Dune::index_constant<freeFlowMomentumIndex> domainJ)
const
424 return massAndEnergyToMomentumStencils_[eIdx];
437 const Element<freeFlowMomentumIndex>& elementI,
438 const SubControlVolume<freeFlowMomentumIndex>& scvI,
439 Dune::index_constant<freeFlowMassIndex> domainJ)
const
441 return momentumToMassAndEnergyStencils_[scvI.index()];
452 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
454 const LocalAssemblerI& localAssemblerI,
455 Dune::index_constant<j> domainJ,
456 std::size_t dofIdxGlobalJ,
457 const PrimaryVariables<j>& priVarsJ,
460 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
464 bindCouplingContext_(domainI, localAssemblerI.element());
467 const auto& deflectedElement =
problem.gridGeometry().element(dofIdxGlobalJ);
469 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
470 const auto scvIdxJ = dofIdxGlobalJ;
471 const auto& scv = fvGeometry.scv(scvIdxJ);
473 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
474 gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
476 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
497 template<std::
size_t i,
class AssembleElementFunc>
500 if (elementSets_.empty())
501 DUNE_THROW(Dune::InvalidStateException,
"Call computeColorsForAssembly before assembling in parallel!");
508 for (
const auto& elements : elementSets_)
512 const auto element = grid.entity(elements[eIdx]);
513 assembleElement(element);
519 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
520 const Element<freeFlowMomentumIndex>& elementI)
const
523 bindCouplingContext_(domainI, elementI, eIdx);
526 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
527 const Element<freeFlowMomentumIndex>& elementI,
528 const std::size_t eIdx)
const
530 if (momentumCouplingContext_().empty())
533 fvGeometry.bind(elementI);
544 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
546 else if (eIdx != momentumCouplingContext_()[0].eIdx)
548 momentumCouplingContext_()[0].eIdx = eIdx;
549 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
550 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->
curSol(
freeFlowMassIndex));
553 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, prevSol_(
freeFlowMassIndex));
557 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
558 const Element<freeFlowMassIndex>& elementI)
const
561 bindCouplingContext_(domainI, elementI, eIdx);
564 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
565 const Element<freeFlowMassIndex>& elementI,
566 const std::size_t eIdx)
const
568 if (massAndEnergyCouplingContext_().empty())
571 auto fvGeometry =
localView(gridGeometry);
572 fvGeometry.bindElement(elementI);
573 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
575 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
577 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
578 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
586 template<std::
size_t i>
587 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
const
589 if (std::get<i>(gridVariables_))
590 return *std::get<i>(gridVariables_);
592 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
599 template<std::
size_t i>
600 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
602 if (std::get<i>(gridVariables_))
603 return *std::get<i>(gridVariables_);
605 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
609 void computeCouplingStencils_()
613 auto momentumFvGeometry =
localView(momentumGridGeometry);
614 massAndEnergyToMomentumStencils_.clear();
615 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
617 momentumToMassAndEnergyStencils_.clear();
618 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
620 for (
const auto& element : elements(momentumGridGeometry.gridView()))
622 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
623 momentumFvGeometry.bind(element);
624 for (
const auto& scv :
scvs(momentumFvGeometry))
626 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
627 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
630 if constexpr (FluidSystem::isCompressible(0))
633 for (
const auto& scvf : scvfs(momentumFvGeometry, scv))
635 if (scvf.isLateral() && !scvf.boundary())
637 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
638 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
646 std::size_t massScvfToMomentumScvIdx_(
const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
647 [[maybe_unused]]
const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry)
const
649 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
650 return massScvf.index();
653 static const bool makeConsistentlyOriented = getParam<bool>(
"Grid.MakeConsistentlyOriented",
true);
654 if (!makeConsistentlyOriented)
655 return massScvf.index();
657 for (
const auto& momentumScv :
scvs(momentumFVGeometry))
659 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
660 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
661 if (Dune::FloatCmp::eq<
typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
662 return momentumScv.index();
664 DUNE_THROW(Dune::InvalidStateException,
"No Momentum SCV found");
668 CouplingStencilType emptyStencil_;
669 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
670 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
672 std::vector<MomentumCouplingContext>& momentumCouplingContext_()
const
673 {
return momentumCouplingContextImpl_; }
675 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_()
const
676 {
return massAndEnergyCouplingContextImpl_; }
678 mutable std::vector<MassAndEnergyCouplingContext> massAndEnergyCouplingContextImpl_;
679 mutable std::vector<MomentumCouplingContext> momentumCouplingContextImpl_;
682 GridVariablesTuple gridVariables_;
684 bool isTransient_()
const
685 {
return std::get<0>(prevSolutions_) !=
nullptr; }
687 template<std::
size_t i>
688 const SubSolutionVector<i>& prevSol_(Dune::index_constant<i>)
const
689 {
return *std::get<i>(prevSolutions_); }
691 PrevSolutionVectorStorage prevSolutions_;
693 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
699:
public std::false_type {};
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
void attachSolution(const SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:311
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
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:53
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
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:60
typename Traits::SolutionVector SolutionVector
the type of the solution vector
Definition: multidomain/couplingmanager.hh:56
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_staggered.hh:48
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: couplingmanager_staggered.hh:236
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_staggered.hh:51
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: couplingmanager_staggered.hh:330
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: couplingmanager_staggered.hh:453
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: couplingmanager_staggered.hh:200
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_staggered.hh:52
void assembleMultithreaded(Dune::index_constant< i > domainI, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_staggered.hh:498
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: couplingmanager_staggered.hh:315
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: couplingmanager_staggered.hh:123
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 scv's residual depends ...
Definition: couplingmanager_staggered.hh:399
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: couplingmanager_staggered.hh:286
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const typename ParentType::SolutionVectorStorage &curSol)
use as binary coupling manager in multi model context
Definition: couplingmanager_staggered.hh:152
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< freeFlowMassIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmanager_staggered.hh:436
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: couplingmanager_staggered.hh:139
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv) const
Returns the pressure at a given sub control volume.
Definition: couplingmanager_staggered.hh:357
static constexpr auto pressureIdx
Definition: couplingmanager_staggered.hh:115
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMassIndex > domainI, const Element< freeFlowMassIndex > &elementI, Dune::index_constant< freeFlowMomentumIndex > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: couplingmanager_staggered.hh:419
void init(std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const typename ParentType::SolutionVectorStorage &curSol, const PrevSolutionVectorStorage &prevSol)
use as binary coupling manager in multi model context and for transient problems
Definition: couplingmanager_staggered.hh:168
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: couplingmanager_staggered.hh:259
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: couplingmanager_staggered.hh:250
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_staggered.hh:369
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_staggered.hh:485
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
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
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
The interface of the coupling manager for multi domain problems.
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:239
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
Parallel for loop (multithreading)
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78