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;
74 using SolutionVector =
typename Traits::SolutionVector;
78 using GridVariablesTuple =
typename Traits::template TupleOfSharedPtr<GridVariables>;
80 using FluidSystem =
typename VolumeVariables<freeFlowMassIndex>::FluidSystem;
82 using VelocityVector =
typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
83 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
85 struct MomentumCouplingContext
87 FVElementGeometry<freeFlowMassIndex> fvGeometry;
88 ElementVolumeVariables<freeFlowMassIndex> curElemVolVars;
89 ElementVolumeVariables<freeFlowMassIndex> prevElemVolVars;
93 struct MassAndEnergyCouplingContext
95 MassAndEnergyCouplingContext(FVElementGeometry<freeFlowMomentumIndex>&& f,
const std::size_t i)
96 : fvGeometry(std::move(f))
100 FVElementGeometry<freeFlowMomentumIndex> fvGeometry;
106 static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx;
114 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
115 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
116 GridVariablesTuple&& gridVariables,
117 const SolutionVector&
curSol)
119 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
120 gridVariables_ = gridVariables;
123 computeCouplingStencils_();
127 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
128 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
129 GridVariablesTuple&& gridVariables,
130 const SolutionVector&
curSol,
131 const SolutionVector& prevSol)
133 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables),
curSol);
139 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
140 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
141 GridVariablesTuple&& gridVariables,
144 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
145 gridVariables_ = gridVariables;
148 computeCouplingStencils_();
172 template<std::
size_t j,
class LocalAssemblerI>
174 const LocalAssemblerI& localAssemblerI,
175 const SubControlVolume<freeFlowMomentumIndex>& scvI,
176 Dune::index_constant<j> domainJ,
177 std::size_t dofIdxGlobalJ)
const
179 const auto&
problem = localAssemblerI.problem();
180 const auto& element = localAssemblerI.element();
181 const auto& fvGeometry = localAssemblerI.fvGeometry();
182 const auto& curElemVolVars = localAssemblerI.curElemVolVars();
183 const auto& prevElemVolVars = localAssemblerI.prevElemVolVars();
184 typename LocalAssemblerI::ElementResidualVector residual(localAssemblerI.element().subEntities(1));
185 const auto& localResidual = localAssemblerI.localResidual();
187 localResidual.evalSource(residual,
problem, element, fvGeometry, curElemVolVars, scvI);
189 for (
const auto& scvf : scvfs(fvGeometry, scvI))
190 localResidual.evalFlux(residual,
problem, element, fvGeometry, curElemVolVars, localAssemblerI.elemBcTypes(), localAssemblerI.elemFluxVarsCache(), scvf);
192 if (!localAssemblerI.assembler().isStationaryProblem())
194 assert(isTransient_);
195 localResidual.evalStorage(residual,
problem, element, fvGeometry, prevElemVolVars, curElemVolVars, scvI);
209 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
210 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
211 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
213 assert(scvf.isFrontal() && !scvf.isLateral() && !scvf.boundary());
224 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
232 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
233 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
234 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
235 const bool considerPreviousTimeStep =
false)
const
237 assert(!(considerPreviousTimeStep && !isTransient_));
238 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
239 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
240 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
242 const auto rho = [&](
const auto& elemVolVars)
245 return elemVolVars[insideMassScv].
density();
248 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
249 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
251 return 0.5*(elemVolVars[insideMassScv].density() + elemVolVars[outsideMassScv].density());
255 return considerPreviousTimeStep ? rho(momentumCouplingContext_()[0].prevElemVolVars)
256 : rho(momentumCouplingContext_()[0].curElemVolVars);
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 result = [&](
const auto& elemVolVars)
272 return std::make_pair(elemVolVars[insideMassScv].
density(), elemVolVars[insideMassScv].density());
275 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
276 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
277 return std::make_pair(elemVolVars[insideMassScv].
density(), elemVolVars[outsideMassScv].
density());
281 return considerPreviousTimeStep ? result(momentumCouplingContext_()[0].prevElemVolVars)
282 : result(momentumCouplingContext_()[0].curElemVolVars);
288 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
289 const SubControlVolume<freeFlowMomentumIndex>& scv,
290 const bool considerPreviousTimeStep =
false)
const
292 assert(!(considerPreviousTimeStep && !isTransient_));
293 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, scv.elementIndex());
294 const auto& massScv = (*scvs(momentumCouplingContext_()[0].fvGeometry).begin());
296 return considerPreviousTimeStep ? momentumCouplingContext_()[0].prevElemVolVars[massScv].density()
297 : momentumCouplingContext_()[0].curElemVolVars[massScv].density();
304 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
305 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf)
const
307 bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex>(), element, fvGeometry.elementIndex());
309 const auto& insideMomentumScv = fvGeometry.scv(scvf.insideScvIdx());
310 const auto& insideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(insideMomentumScv.elementIndex());
313 return momentumCouplingContext_()[0].curElemVolVars[insideMassScv].viscosity();
315 const auto& outsideMomentumScv = fvGeometry.scv(scvf.outsideScvIdx());
316 const auto& outsideMassScv = momentumCouplingContext_()[0].fvGeometry.scv(outsideMomentumScv.elementIndex());
318 const auto mu = [&](
const auto& elemVolVars)
321 return 0.5*(elemVolVars[insideMassScv].viscosity() + elemVolVars[outsideMassScv].viscosity());
324 return mu(momentumCouplingContext_()[0].curElemVolVars);
331 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
334 bindCouplingContext_(Dune::index_constant<freeFlowMassIndex>(), element, scvf.insideScvIdx());
337 const auto localMomentumScvIdx = massScvfToMomentumScvIdx_(scvf, massAndEnergyCouplingContext_()[0].fvGeometry);
338 const auto& scvJ = massAndEnergyCouplingContext_()[0].fvGeometry.scv(localMomentumScvIdx);
341 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition velocity;
342 velocity[scvJ.dofAxis()] = 1.0;
359 template<std::
size_t j>
361 const Element<freeFlowMomentumIndex>& elementI,
362 const SubControlVolume<freeFlowMomentumIndex>& scvI,
363 Dune::index_constant<j> domainJ)
const
364 {
return emptyStencil_; }
381 const Element<freeFlowMassIndex>& elementI,
382 Dune::index_constant<freeFlowMomentumIndex> domainJ)
const
385 return massAndEnergyToMomentumStencils_[eIdx];
398 const Element<freeFlowMomentumIndex>& elementI,
399 const SubControlVolume<freeFlowMomentumIndex>& scvI,
400 Dune::index_constant<freeFlowMassIndex> domainJ)
const
402 return momentumToMassAndEnergyStencils_[scvI.index()];
413 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
415 const LocalAssemblerI& localAssemblerI,
416 Dune::index_constant<j> domainJ,
417 std::size_t dofIdxGlobalJ,
418 const PrimaryVariables<j>& priVarsJ,
421 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
425 bindCouplingContext_(domainI, localAssemblerI.element());
428 const auto& deflectedElement =
problem.gridGeometry().element(dofIdxGlobalJ);
430 const auto& fvGeometry = momentumCouplingContext_()[0].fvGeometry;
431 const auto scvIdxJ = dofIdxGlobalJ;
432 const auto& scv = fvGeometry.scv(scvIdxJ);
434 if constexpr (ElementVolumeVariables<freeFlowMassIndex>::GridVolumeVariables::cachingEnabled)
435 gridVars_(
freeFlowMassIndex).curGridVolVars().volVars(scv).update(std::move(elemSol),
problem, deflectedElement, scv);
437 momentumCouplingContext_()[0].curElemVolVars[scv].update(std::move(elemSol),
problem, deflectedElement, scv);
458 template<std::
size_t i,
class AssembleElementFunc>
461 if (elementSets_.empty())
462 DUNE_THROW(Dune::InvalidStateException,
"Call computeColorsForAssembly before assembling in parallel!");
469 for (
const auto& elements : elementSets_)
473 const auto element = grid.entity(elements[eIdx]);
474 assembleElement(element);
480 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
481 const Element<freeFlowMomentumIndex>& elementI)
const
484 bindCouplingContext_(domainI, elementI, eIdx);
487 void bindCouplingContext_(Dune::index_constant<freeFlowMomentumIndex> domainI,
488 const Element<freeFlowMomentumIndex>& elementI,
489 const std::size_t eIdx)
const
491 if (momentumCouplingContext_().empty())
494 fvGeometry.bind(elementI);
503 prevElemVolVars.bindElement(elementI, fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
505 momentumCouplingContext_().emplace_back(MomentumCouplingContext{std::move(fvGeometry), std::move(curElemVolVars), std::move(prevElemVolVars), eIdx});
507 else if (eIdx != momentumCouplingContext_()[0].eIdx)
509 momentumCouplingContext_()[0].eIdx = eIdx;
510 momentumCouplingContext_()[0].fvGeometry.bind(elementI);
511 momentumCouplingContext_()[0].curElemVolVars.bind(elementI, momentumCouplingContext_()[0].fvGeometry, this->
curSol(freeFlowMassIndex));
514 momentumCouplingContext_()[0].prevElemVolVars.bindElement(elementI, momentumCouplingContext_()[0].fvGeometry, (*prevSol_)[
freeFlowMassIndex]);
518 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
519 const Element<freeFlowMassIndex>& elementI)
const
521 const auto eIdx = this->
problem(freeFlowMassIndex).gridGeometry().elementMapper().index(elementI);
522 bindCouplingContext_(domainI, elementI, eIdx);
525 void bindCouplingContext_(Dune::index_constant<freeFlowMassIndex> domainI,
526 const Element<freeFlowMassIndex>& elementI,
527 const std::size_t eIdx)
const
529 if (massAndEnergyCouplingContext_().empty())
532 auto fvGeometry =
localView(gridGeometry);
533 fvGeometry.bindElement(elementI);
534 massAndEnergyCouplingContext_().emplace_back(std::move(fvGeometry), eIdx);
536 else if (eIdx != massAndEnergyCouplingContext_()[0].eIdx)
538 massAndEnergyCouplingContext_()[0].eIdx = eIdx;
539 massAndEnergyCouplingContext_()[0].fvGeometry.bindElement(elementI);
547 template<std::
size_t i>
548 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
const
550 if (std::get<i>(gridVariables_))
551 return *std::get<i>(gridVariables_);
553 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
560 template<std::
size_t i>
561 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
563 if (std::get<i>(gridVariables_))
564 return *std::get<i>(gridVariables_);
566 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
570 void computeCouplingStencils_()
574 auto momentumFvGeometry =
localView(momentumGridGeometry);
575 massAndEnergyToMomentumStencils_.clear();
576 massAndEnergyToMomentumStencils_.resize(momentumGridGeometry.gridView().size(0));
578 momentumToMassAndEnergyStencils_.clear();
579 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.numScv());
581 for (
const auto& element : elements(momentumGridGeometry.gridView()))
583 const auto eIdx = momentumGridGeometry.elementMapper().index(element);
584 momentumFvGeometry.bind(element);
585 for (
const auto& scv : scvs(momentumFvGeometry))
587 massAndEnergyToMomentumStencils_[eIdx].push_back(scv.dofIndex());
588 momentumToMassAndEnergyStencils_[scv.index()].push_back(eIdx);
591 if constexpr (FluidSystem::isCompressible(0))
594 for (
const auto& scvf : scvfs(momentumFvGeometry, scv))
596 if (scvf.isLateral() && !scvf.boundary())
598 const auto& outsideScv = momentumFvGeometry.scv(scvf.outsideScvIdx());
599 momentumToMassAndEnergyStencils_[scv.index()].push_back(outsideScv.elementIndex());
607 std::size_t massScvfToMomentumScvIdx_(
const SubControlVolumeFace<freeFlowMassIndex>& massScvf,
608 [[maybe_unused]]
const FVElementGeometry<freeFlowMomentumIndex>& momentumFVGeometry)
const
610 if constexpr (ConsistentlyOrientedGrid<typename GridView<freeFlowMomentumIndex>::Grid>{})
611 return massScvf.index();
614 static const bool makeConsistentlyOriented = getParam<bool>(
"Grid.MakeConsistentlyOriented",
true);
615 if (!makeConsistentlyOriented)
616 return massScvf.index();
618 for (
const auto& momentumScv : scvs(momentumFVGeometry))
620 typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition momentumUnitOuterNormal(0.0);
621 momentumUnitOuterNormal[momentumScv.dofAxis()] = momentumScv.directionSign();
622 if (Dune::FloatCmp::eq<
typename GridView<freeFlowMomentumIndex>::ctype>(massScvf.unitOuterNormal()*momentumUnitOuterNormal, 1.0))
623 return momentumScv.index();
625 DUNE_THROW(Dune::InvalidStateException,
"No Momentum SCV found");
629 CouplingStencilType emptyStencil_;
630 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
631 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
636 std::vector<MomentumCouplingContext>& momentumCouplingContext_()
const
638 thread_local static std::vector<MomentumCouplingContext> c;
643 std::vector<MassAndEnergyCouplingContext>& massAndEnergyCouplingContext_()
const
645 thread_local static std::vector<MassAndEnergyCouplingContext> c;
650 GridVariablesTuple gridVariables_;
652 const SolutionVector* prevSol_;
655 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
661:
public std::false_type {};
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
void attachSolution(SolutionVectorStorage &curSol)
Attach a solution vector stored outside of this class.
Definition: multidomain/couplingmanager.hh:322
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:287
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:309
std::vector< std::size_t > CouplingStencilType
default type used for coupling element stencils
Definition: multidomain/couplingmanager.hh:64
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:338
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:219
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:71
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:209
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:303
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:414
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:173
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:459
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:288
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:114
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:360
auto insideAndOutsideDensity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Definition: couplingmanager_staggered.hh:259
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:397
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:127
static constexpr auto pressureIdx
Definition: couplingmanager_staggered.hh:106
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:380
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:232
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:139
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:223
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:330
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_staggered.hh:446
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
Parallel for loop (multithreading)
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78