24#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_HH
25#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_STAGGERED_HH
32#include <dune/common/exceptions.hh>
33#include <dune/common/indices.hh>
34#include <dune/common/float_cmp.hh>
70 template<std::
size_t id>
using SubDomainTypeTag =
typename Traits::template SubDomain<id>::TypeTag;
73 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
74 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
75 template<std::
size_t id>
using ElementSeed =
typename GridView<id>::Grid::template Codim<0>::EntitySeed;
76 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
77 template<std::
size_t id>
using SubControlVolume =
typename FVElementGeometry<id>::SubControlVolume;
78 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
79 template<std::
size_t id>
using GridVariables =
typename Traits::template SubDomain<id>::GridVariables;
80 template<std::
size_t id>
using ElementVolumeVariables =
typename GridVariables<id>::GridVolumeVariables::LocalView;
81 template<std::
size_t id>
using GridFluxVariablesCache =
typename GridVariables<id>::GridFluxVariablesCache;
85 using Scalar =
typename Traits::Scalar;
86 using SolutionVector =
typename Traits::SolutionVector;
90 using GridVariablesTuple =
typename Traits::template TupleOfSharedPtr<GridVariables>;
92 using FluidSystem =
typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
94 using VelocityVector =
typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
95 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
97 struct MomentumCouplingContext
99 FVElementGeometry<freeFlowMassIndex> fvGeometry;
100 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
101 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
105 struct MassAndEnergyCouplingContext
107 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f,
const std::size_t i)
108 : fvGeometry(std::move(f))
112 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
118 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
126 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
127 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
128 GridVariablesTuple&& gridVariables,
129 const SolutionVector&
curSol)
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,
142 const SolutionVector&
curSol,
143 const SolutionVector& prevSol)
145 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables),
curSol);
151 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
152 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
153 GridVariablesTuple&& gridVariables,
156 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
157 gridVariables_ = gridVariables;
160 computeCouplingStencils_();
185 template<std::
size_t j,
class LocalAssemblerI>
187 const LocalAssemblerI& localAssemblerI,
188 const SubControlVolume<freeFlowMomentumIndex>& scvI,
189 Dune::index_constant<j> domainJ,
190 std::size_t dofIdxGlobalJ)
const
192 const auto&
problem = localAssemblerI.problem();
193 const auto& element = localAssemblerI.element();
194 const auto& fvGeometry = localAssemblerI.fvGeometry();
195 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
196 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
197 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
198 const auto& localResidual = localAssemblerI.localResidual();
200 localResidual.evalSource(residual,
problem, element, fvGeometry, curElemVolVars, scvI);
202 for (
const auto& scvf : scvfs(fvGeometry, scvI))
203 localResidual.evalFlux(residual,
problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
205 if (!localAssemblerI.assembler().isStationaryProblem())
207 assert(isTransient_);
208 localResidual.evalStorage(residual,
problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
222 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
223 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
224 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
226 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
237 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
245 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
246 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
247 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
248 const bool considerPreviousTimeStep =
false)
const
250 assert(!(considerPreviousTimeStep && !isTransient_));
251 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
252 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
253 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
255 const auto rho = [&](
const auto& elemVolVars)
258 return elemVolVars[insideMassScv].
density();
261 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
262 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
264 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
268 return considerPreviousTimeStep ? rho(momentumCouplingContext_()[0].prevElemVolVars)
269 : rho(momentumCouplingContext_()[0].curElemVolVars);
273 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
274 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
275 const bool considerPreviousTimeStep =
false)
const
277 assert(!(considerPreviousTimeStep && !isTransient_));
278 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
279 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
280 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
282 const auto result = [&](
const auto& elemVolVars)
285 return std::make_pair(elemVolVars[insideMassScv].
density(), elemVolVars[insideMassScv].density());
288 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
289 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
290 return std::make_pair(elemVolVars[insideMassScv].
density(), elemVolVars[outsideMassScv].
density());
294 return considerPreviousTimeStep ? result(momentumCouplingContext_()[0].prevElemVolVars)
295 : result(momentumCouplingContext_()[0].curElemVolVars);
301 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
302 const SubControlVolume<freeFlowMomentumIndex>& scv,
303 const bool considerPreviousTimeStep =
false)
const
305 assert(!(considerPreviousTimeStep && !isTransient_));
306 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
307 const auto& massScv = (*scvs(momentumCouplingContext_()[0].fvGeometry).begin());
309 return considerPreviousTimeStep ? momentumCouplingContext_()[0].prevElemVolVars[massScv].density()
310 : momentumCouplingContext_()[0].curElemVolVars[massScv].density();
317 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
318 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
320 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
322 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
323 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
326 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
328 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
329 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
331 const auto mu = [&](
const auto& elemVolVars)
334 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
337 return mu(momentumCouplingContext_()[0].curElemVolVars);
344 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
347 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx());
350 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_()[0].fvGeometry);
351 const auto& scvJ = massAndEnergyCouplingContext_()[0].fvGeometry.scv(localMomentumScvIdx);
354 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
355 velocity[scvJ.dofAxis()] = 1.0;
367 template<std::
size_t j>
369 const Element<freeFlowMomentumIndex>& elementI,
370 const SubControlVolume<freeFlowMomentumIndex>& scvI,
371 Dune::index_constant<j> domainJ)
const
372 {
return emptyStencil_; }
389 const Element<freeFlowMassIndex>& elementI,
390 Dune::index_constant<freeFlowMomentumIndex> domainJ)
const
393 return massAndEnergyToMomentumStencils_[eIdx];
406 const Element<freeFlowMomentumIndex>& elementI,
407 const SubControlVolume<freeFlowMomentumIndex>& scvI,
408 Dune::index_constant<freeFlowMassIndex> domainJ)
const
410 return momentumToMassAndEnergyStencils_[scvI.index()];
438 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
440 const LocalAssemblerI& localAssemblerI,
441 Dune::index_constant<j> domainJ,
442 std::size_t dofIdxGlobalJ,
443 const PrimaryVariables<j>& priVarsJ,
446 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
450 bindCouplingContext_(domainI, localAssemblerI.element());
453 const auto& deflectedElement =
problem.gridGeometry().element(dofIdxGlobalJ);
455 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
456 const auto scvIdxJ = dofIdxGlobalJ;
457 const auto& scv = fvGeometry.scv(scvIdxJ);
459 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
460 gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
462 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
486 template<std::
size_t i,
class AssembleElementFunc>
489 if (elementSets_.empty())
490 DUNE_THROW(Dune::InvalidStateException,
"Call computeColorsForAssembly before assembling in parallel!");
497 for (
const auto& elements : elementSets_)
501 const auto element = grid.entity(elements[eIdx]);
502 assembleElement(element);
508 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
509 const Element<freeFlowMomentumIndex>& elementI)
const
512 bindCouplingContext_(domainI, elementI, eIdx);
515 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
516 const Element<freeFlowMomentumIndex>& elementI,
517 const std::size_t eIdx)
const
519 if (momentumCouplingContext_().empty())
522 fvGeometry.bind(elementI);
531 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
533 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
535 else if (eIdx != momentumCouplingContext_()[0].eIdx)
537 momentumCouplingContext_()[0].eIdx = eIdx;
538 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
539 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->
curSol(freeFlowMassIndex));
542 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
546 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
547 const Element<freeFlowMassIndex>& elementI)
const
549 const auto eIdx = this->
problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
550 bindCouplingContext_(domainI, elementI, eIdx);
553 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
554 const Element<freeFlowMassIndex>& elementI,
555 const std::size_t eIdx)
const
557 if (massAndEnergyCouplingContext_().empty())
560 auto fvGeometry =
localView(gridGeometry);
561 fvGeometry.bindElement(elementI);
562 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
564 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
566 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
567 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
575 template<std::
size_t i>
576 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
const
578 if (std::get<i>(gridVariables_))
579 return *std::get<i>(gridVariables_);
581 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
588 template<std::
size_t i>
589 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
591 if (std::get<i>(gridVariables_))
592 return *std::get<i>(gridVariables_);
594 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
598 void computeCouplingStencils_()
602 auto momentumFvGeometry =
localView(momentumGridGeometry);
603 massAndEnergyToMomentumStencils_.clear();
604 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
606 momentumToMassAndEnergyStencils_.clear();
607 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
609 for (
const auto& element : elements(momentumGridGeometry.gridView()))
611 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
612 momentumFvGeometry.bind(element);
613 for (
const auto& scv : scvs(momentumFvGeometry))
615 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
616 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
619 if constexpr (FluidSystem::isCompressible(0))
622 for (
const auto& scvf : scvfs(momentumFvGeometry, scv))
624 if (scvf.isLateral() && !scvf.boundary())
626 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
627 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
635 std::size_t massScvfToMomentumScvIdx_(
const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
636 [[maybe_unused]]
const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry)
const
638 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
639 return massScvf.index();
642 static const bool makeConsistentlyOriented = getParam<bool>(
"Grid.MakeConsistentlyOriented",
true);
643 if (!makeConsistentlyOriented)
644 return massScvf.index();
646 for (
const auto& momentumScv : scvs(momentumFVGeometry))
648 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
649 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
650 if (Dune::FloatCmp::eq<
typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
651 return momentumScv.index();
653 DUNE_THROW(Dune::InvalidStateException,
"No Momentum SCV found");
657 CouplingStencilType emptyStencil_;
658 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
659 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
664 std::vector<MomentumCouplingContext>& momentumCouplingContext_()
const
666 thread_local static std::vector<MomentumCouplingContext> c;
671 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_()
const
673 thread_local static std::vector<MassAndEnergyCouplingContext> c;
678 GridVariablesTuple gridVariables_;
680 const SolutionVector* prevSol_;
683 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
689:
public std::true_type {};
Coloring schemes for shared-memory-parallel assembly.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
The available discretization methods in Dumux.
Parallel for loop (multithreading)
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
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:172
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:186
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:251
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:55
Definition: common/properties.hh:100
The secondary variables within a sub-control volume.
Definition: common/properties.hh:105
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
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
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:350
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: couplingmanager_staggered.hh:60
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:222
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_staggered.hh:63
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:316
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:439
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_staggered.hh:64
void assembleMultithreaded(Dune::index_constant< i > domainI, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_staggered.hh:487
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:301
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:126
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: couplingmanager_staggered.hh:368
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: couplingmanager_staggered.hh:272
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:405
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
static constexpr auto pressureIdx
Definition: couplingmanager_staggered.hh:118
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:388
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:245
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: couplingmanager_staggered.hh:151
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:236
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:343
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_staggered.hh:474
trait that is specialized for coupling manager supporting multithreaded assembly
Definition: multidomain/fvassembler.hh:85
A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
Declares all properties used in Dumux.
The interface of the coupling manager for multi domain problems.