12#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_CVFE_HH
13#define DUMUX_MULTIDOMAIN_FREEFLOW_COUPLING_MANAGER_CVFE_HH
21#include <dune/common/exceptions.hh>
22#include <dune/common/indices.hh>
23#include <dune/common/float_cmp.hh>
24#include <dune/geometry/referenceelements.hh>
27#include <dumux/common/concepts/variables_.hh>
64 template<std::
size_t id>
using SubDomainTypeTag =
typename Traits::template SubDomain<id>::TypeTag;
67 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
68 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
69 template<std::
size_t id>
using ElementSeed =
typename GridView<id>::Grid::template Codim<0>::EntitySeed;
70 template<std::
size_t id>
using FVElementGeometry =
typename GridGeometry<id>::LocalView;
71 template<std::
size_t id>
using SubControlVolume =
typename FVElementGeometry<id>::SubControlVolume;
72 template<std::
size_t id>
using SubControlVolumeFace =
typename FVElementGeometry<id>::SubControlVolumeFace;
73 template<std::
size_t id>
using GridVariables =
typename Traits::template SubDomain<id>::GridVariables;
74 template<std::
size_t id>
using GridVariablesCache = Concept::GridVariablesCache_t<GridVariables<id>>;
75 template<std::
size_t id>
using ElementVariables =
typename GridVariablesCache<id>::LocalView;
77 template<std::
size_t id>
using Variables = Concept::Variables_t<GridVariables<id>>;
79 using Scalar =
typename Traits::Scalar;
80 using SolutionVector =
typename Traits::SolutionVector;
84 using GridVariablesTuple =
typename Traits::template TupleOfSharedPtr<GridVariables>;
86 using FluidSystem =
typename Variables<freeFlowMassIndex>::FluidSystem;
88 using GlobalPosition =
typename SubControlVolumeFace<freeFlowMassIndex>::GlobalPosition;
89 using VelocityVector = GlobalPosition;
90 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
92 static_assert(std::is_same_v<VelocityVector, typename SubControlVolumeFace<freeFlowMomentumIndex>::GlobalPosition>);
94 template<
class ElementSolution>
95 struct MomentumCouplingContextNoCaching
97 MomentumCouplingContextNoCaching(ElementSolution&& elemSol)
98 : elemSol_(std::move(elemSol)) {}
100 template<
class Gr
idVarsCache,
class FvElementGeometry,
class SubControlVolume>
101 auto vars(
const GridVarsCache& gridVarsCache,
const FvElementGeometry& fvGeometry,
const SubControlVolume& scv)
const
103 const auto&
problem = gridVarsCache.problem();
104 Variables<freeFlowMassIndex> variables;
105 if constexpr (Concept::FVGridVariables<GridVariables<freeFlowMassIndex>>)
106 variables.update(elemSol_,
problem, fvGeometry.element(), scv);
108 variables.update(elemSol_,
problem, fvGeometry, ipData(fvGeometry, scv));
112 ElementSolution elemSol_;
115 struct MomentumCouplingContextGlobalCaching
117 template<
class Gr
idVarsCache,
class FvElementGeometry,
class SubControlVolume>
118 const auto& vars(
const GridVarsCache& gridVarsCache,
const FvElementGeometry& fvGeometry,
const SubControlVolume& scv)
const
120 if constexpr (
requires { gridVarsCache.volVars(scv); })
121 return gridVarsCache.volVars(scv);
123 return gridVarsCache.variables(scv);
127 using MomentumDiscretizationMethod =
typename GridGeometry<freeFlowMomentumIndex>::DiscretizationMethod;
128 using MassDiscretizationMethod =
typename GridGeometry<freeFlowMassIndex>::DiscretizationMethod;
130 template<std::
size_t id>
using IpData
132 typename GridView<id>::template Codim<0>::Entity::Geometry::GlobalCoordinate>;
136 static constexpr auto pressureIdx = Variables<freeFlowMassIndex>::Indices::pressureIdx;
144 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
145 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
146 GridVariablesTuple&& gridVariables,
147 const SolutionVector&
curSol)
149 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
150 gridVariables_ = gridVariables;
153 computeCouplingStencils_();
157 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
158 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
159 GridVariablesTuple&& gridVariables,
160 const SolutionVector&
curSol,
161 const SolutionVector& prevSol)
163 init(momentumProblem, massProblem, std::forward<GridVariablesTuple>(gridVariables),
curSol);
169 void init(std::shared_ptr<Problem<freeFlowMomentumIndex>> momentumProblem,
170 std::shared_ptr<Problem<freeFlowMassIndex>> massProblem,
171 GridVariablesTuple&& gridVariables,
174 this->
setSubProblems(std::make_tuple(momentumProblem, massProblem));
175 gridVariables_ = gridVariables;
178 computeCouplingStencils_();
191 [[deprecated(
"This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
192 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
193 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
194 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
195 const bool considerPreviousTimeStep =
false)
const
197 const auto& globalPos = scvf.ipGlobal();
198 const auto& localPos = element.geometry().local(globalPos);
205 [[deprecated(
"This method will be removed after release (3.11). Use pressure(..., ipData) instead!")]]
206 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
207 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
208 const SubControlVolume<freeFlowMomentumIndex>& scv,
209 const bool considerPreviousTimeStep =
false)
const
211 return this->
pressure(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
217 template <
class IpData>
218 Scalar
pressure(
const Element<freeFlowMomentumIndex>& element,
219 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
221 const bool considerPreviousTimeStep =
false)
const
223 assert(!(considerPreviousTimeStep && !this->isTransient_));
234 [[deprecated(
"This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
235 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
236 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
237 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
238 const bool considerPreviousTimeStep =
false)
const
240 const auto& globalPos = scvf.ipGlobal();
241 const auto& localPos = element.geometry().local(globalPos);
248 [[deprecated(
"This method will be removed after release (3.11). Use density(..., ipData) instead!")]]
249 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
250 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
251 const SubControlVolume<freeFlowMomentumIndex>& scv,
252 const bool considerPreviousTimeStep =
false)
const
254 return this->
density(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
260 template <
class IpData>
261 Scalar
density(
const Element<freeFlowMomentumIndex>& element,
262 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
264 const bool considerPreviousTimeStep =
false)
const
266 assert(!(considerPreviousTimeStep && !this->isTransient_));
272 massFvGeometry.bind(element);
273 const auto context = makeMomentumCouplingContext_(massFvGeometry, sol);
274 const auto& gridVarsCache = subDomainGridVars_(Dune::index_constant<freeFlowMassIndex>{}, considerPreviousTimeStep);
278 const auto eIdx = fvGeometry.elementIndex();
279 const auto& scv = massFvGeometry.scv(eIdx);
281 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
283 return variables.density();
289 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
290 const auto& localBasis = massFvGeometry.feLocalBasis();
291 std::vector<ShapeValue> shapeValues;
292 localBasis.evaluateFunction(ipData.
local(), shapeValues);
295 for (
const auto& scv :
scvs(massFvGeometry))
297 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
298 rho += variables.density()*shapeValues[scv.localDofIndex()][0];
304 DUNE_THROW(Dune::NotImplemented,
305 "Density interpolation for discretization scheme " << MassDiscretizationMethod{}
312 [[deprecated(
"This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
314 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
315 const SubControlVolumeFace<freeFlowMomentumIndex>& scvf,
316 const bool considerPreviousTimeStep =
false)
const
318 const auto& globalPos = scvf.ipGlobal();
319 const auto& localPos = element.geometry().local(globalPos);
326 [[deprecated(
"This method will be removed after release (3.11). Use effectiveViscosity(..., ipData) instead!")]]
328 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
329 const SubControlVolume<freeFlowMomentumIndex>& scv,
330 const bool considerPreviousTimeStep =
false)
const
332 return this->
effectiveViscosity(element, fvGeometry, ipData(fvGeometry, scv), considerPreviousTimeStep);
338 template <
class IpData>
340 const FVElementGeometry<freeFlowMomentumIndex>& fvGeometry,
342 const bool considerPreviousTimeStep =
false)
const
344 assert(!(considerPreviousTimeStep && !this->isTransient_));
350 massFvGeometry.bind(element);
351 const auto context = makeMomentumCouplingContext_(massFvGeometry, sol);
352 const auto& gridVarsCache = subDomainGridVars_(Dune::index_constant<freeFlowMassIndex>{}, considerPreviousTimeStep);
356 const auto eIdx = fvGeometry.elementIndex();
357 const auto& scv = massFvGeometry.scv(eIdx);
359 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
361 return variables.viscosity();
367 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
368 const auto& localBasis = massFvGeometry.feLocalBasis();
369 std::vector<ShapeValue> shapeValues;
370 localBasis.evaluateFunction(ipData.
local(), shapeValues);
373 for (
const auto& scv :
scvs(massFvGeometry))
375 const auto& variables = context.vars(gridVarsCache, massFvGeometry, scv);
376 mu += variables.viscosity()*shapeValues[scv.localDofIndex()][0];
382 DUNE_THROW(Dune::NotImplemented,
383 "Viscosity interpolation for discretization scheme " << MassDiscretizationMethod{}
391 const SubControlVolumeFace<freeFlowMassIndex>& scvf)
const
395 fvGeometry.bindElement(element);
397 const auto& localBasis = fvGeometry.feLocalBasis();
399 std::vector<ShapeValue> shapeValues;
400 const auto ipLocal = element.geometry().local(scvf.ipGlobal());
401 localBasis.evaluateFunction(ipLocal, shapeValues);
405 for (
const auto& localDof :
localDofs(fvGeometry))
414 VelocityVector
elementVelocity(
const FVElementGeometry<freeFlowMassIndex>& fvGeometry)
const
417 momentumFvGeometry.bindElement(fvGeometry.element());
419 const auto& localBasis = momentumFvGeometry.feLocalBasis();
423 std::vector<ShapeValue> shapeValues;
424 localBasis.evaluateFunction(referenceElement(fvGeometry.element()).position(0,0), shapeValues);
426 for (
const auto& localDof :
localDofs(momentumFvGeometry))
435 template <
class IpData>
436 VelocityVector
velocity(
const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
438 const bool considerPreviousTimeStep =
false)
const
440 assert(!(considerPreviousTimeStep && !this->isTransient_));
442 const auto& element = fvGeometry.element();
445 momentumFvGeometry.bindElement(fvGeometry.element());
458 template<std::
size_t j>
460 const Element<freeFlowMomentumIndex>& elementI,
461 const SubControlVolume<freeFlowMomentumIndex>& scvI,
462 Dune::index_constant<j> domainJ)
const
463 {
return emptyStencil_; }
480 const Element<freeFlowMassIndex>& elementI,
481 Dune::index_constant<freeFlowMomentumIndex> domainJ)
const
484 return massAndEnergyToMomentumStencils_[eIdx];
496 const Element<freeFlowMomentumIndex>& elementI,
497 Dune::index_constant<freeFlowMassIndex> domainJ)
const
500 return momentumToMassAndEnergyStencils_[eIdx];
529 template<std::
size_t i, std::
size_t j,
class LocalAssemblerI>
531 const LocalAssemblerI& localAssemblerI,
532 Dune::index_constant<j> domainJ,
533 std::size_t dofIdxGlobalJ,
534 const PrimaryVariables<j>& priVarsJ,
537 this->
curSol(domainJ)[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
539 if constexpr (GridVariablesCache<freeFlowMassIndex>::cachingEnabled)
546 const auto& deflectedElement =
problem.gridGeometry().element(dofIdxGlobalJ);
549 fvGeometry.bind(deflectedElement);
550 const auto& scv = fvGeometry.scv(dofIdxGlobalJ);
552 if constexpr (Concept::FVGridVariables<GridVariables<freeFlowMassIndex>>)
553 subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{},
true, scv).update(std::move(elemSol),
problem, deflectedElement, scv);
555 subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{},
true, ipData(fvGeometry, scv)).update(std::move(elemSol),
problem, deflectedElement, scv);
564 const auto deflectedElementIdx =
problem.gridGeometry().elementMapper().index(localAssemblerI.element());
565 const auto& deflectedElement =
problem.gridGeometry().element(deflectedElementIdx);
568 fvGeometry.bind(deflectedElement);
571 for (
const auto& scv :
scvs(fvGeometry))
573 if (scv.dofIndex() == dofIdxGlobalJ)
574 if constexpr (Concept::FVGridVariables<GridVariables<freeFlowMassIndex>>)
575 this->subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{},
true, scv).update(std::move(elemSol),
problem, deflectedElement, scv);
577 this->subDomainVariables_(Dune::index_constant<freeFlowMassIndex>{},
true, scv).update(std::move(elemSol),
problem, fvGeometry, ipData(fvGeometry, scv));
582 DUNE_THROW(Dune::NotImplemented,
583 "Context update for discretization scheme " << MassDiscretizationMethod{}
614 template<std::
size_t i,
class AssembleElementFunc>
617 if (elementSets_.empty())
618 DUNE_THROW(Dune::InvalidStateException,
"Call computeColorsForAssembly before assembling in parallel!");
625 for (
const auto& elements : elementSets_)
629 const auto element = grid.entity(elements[eIdx]);
630 assembleElement(element);
636 template<std::
size_t i>
637 const auto& subDomainGridVars_(Dune::index_constant<i> domainIdx,
bool current)
const
639 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
640 return current ? gridVars_(domainIdx).curGridVolVars()
641 : gridVars_(domainIdx).prevGridVolVars();
643 return current ? gridVars_(domainIdx).curGridVars()
644 : gridVars_(domainIdx).prevGridVars();
647 template<std::
size_t i>
648 auto& subDomainGridVars_(Dune::index_constant<i> domainIdx,
bool current)
650 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
651 return current ? gridVars_(domainIdx).curGridVolVars()
652 : gridVars_(domainIdx).prevGridVolVars();
654 return current ? gridVars_(domainIdx).curGridVars()
655 : gridVars_(domainIdx).prevGridVars();
658 template<std::
size_t i,
class Scv>
659 auto& subDomainVariables_(Dune::index_constant<i> domainIdx,
bool current,
const Scv& scv)
661 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
662 return subDomainGridVars_(domainIdx, current).volVars(scv);
664 return subDomainGridVars_(domainIdx, current).variables(scv);
667 template<std::
size_t i,
class Scv>
668 const auto& subDomainVariables_(Dune::index_constant<i> domainIdx,
bool current,
const Scv& scv)
const
670 if constexpr (Concept::FVGridVariables<GridVariables<i>>)
671 return subDomainGridVars_(domainIdx, current).volVars(scv);
673 return subDomainGridVars_(domainIdx, current).variables(scv);
676 template<
class SolutionVector>
677 auto makeMomentumCouplingContext_(
const FVElementGeometry<freeFlowMassIndex>& fvGeometry,
678 const SolutionVector& sol)
const
680 if constexpr (GridVariablesCache<freeFlowMassIndex>::cachingEnabled)
681 return MomentumCouplingContextGlobalCaching{};
683 return MomentumCouplingContextNoCaching{
elementSolution(fvGeometry.element(), sol, fvGeometry.gridGeometry())};
690 template<std::
size_t i>
691 const GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
const
693 if (std::get<i>(gridVariables_))
694 return *std::get<i>(gridVariables_);
696 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
703 template<std::
size_t i>
704 GridVariables<i>& gridVars_(Dune::index_constant<i> domainIdx)
706 if (std::get<i>(gridVariables_))
707 return *std::get<i>(gridVariables_);
709 DUNE_THROW(Dune::InvalidStateException,
"The gridVariables pointer was not set. Use setGridVariables() before calling this function");
712 void computeCouplingStencils_()
716 auto momentumFvGeometry =
localView(momentumGridGeometry);
717 auto massFvGeometry =
localView(massGridGeometry);
719 massAndEnergyToMomentumStencils_.clear();
720 massAndEnergyToMomentumStencils_.resize(massGridGeometry.gridView().size(0));
722 momentumToMassAndEnergyStencils_.clear();
723 momentumToMassAndEnergyStencils_.resize(momentumGridGeometry.gridView().size(0));
725 assert(massAndEnergyToMomentumStencils_.size() == momentumToMassAndEnergyStencils_.size());
726 bool hasCubeElements =
false;
728 for (
const auto& element : elements(momentumGridGeometry.gridView()))
731 hasCubeElements =
true;
733 momentumFvGeometry.bindElement(element);
734 massFvGeometry.bindElement(element);
735 const auto eIdx = momentumFvGeometry.elementIndex();
737 for (
const auto& localDof :
localDofs(momentumFvGeometry))
738 massAndEnergyToMomentumStencils_[eIdx].push_back(localDof.dofIndex());
741 for (
const auto& scv :
scvs(massFvGeometry))
742 momentumToMassAndEnergyStencils_[eIdx].push_back(scv.dofIndex());
748 if(hasCubeElements && !GridGeometry<freeFlowMomentumIndex>::enableHybridCVFE)
750 std::cerr <<
"Warning: Coupled Navier-Stokes problem on cube elements uses non-hybrid pq1bubble. "
751 <<
"The hybrid variant is recommended because it implements two bubble functions for stability reasons."
757 CouplingStencilType emptyStencil_;
758 std::vector<CouplingStencilType> momentumToMassAndEnergyStencils_;
759 std::vector<CouplingStencilType> massAndEnergyToMomentumStencils_;
762 GridVariablesTuple gridVariables_;
764 const SolutionVector* prevSol_;
767 std::deque<std::vector<ElementSeed<freeFlowMomentumIndex>>> elementSets_;
773template<class Traits, class DiscretizationMethod = typename Detail::MomentumDiscretizationMethod<Traits>::type>
777template<
class Traits,
class D>
779{
using type = std::true_type; };
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:31
const LocalPosition & local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:46
The interface of the coupling manager for free flow systems.
Definition: couplingmanager_cvfe.hh:54
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given interpolation point.
Definition: couplingmanager_cvfe.hh:218
VelocityVector faceVelocity(const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
Returns the velocity at a given sub control volume face.
Definition: couplingmanager_cvfe.hh:390
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_cvfe.hh:144
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_cvfe.hh:169
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given sub control volume.
Definition: couplingmanager_cvfe.hh:327
static constexpr auto freeFlowMomentumIndex
Definition: couplingmanager_cvfe.hh:57
static constexpr auto pressureIdx
Definition: couplingmanager_cvfe.hh:136
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_cvfe.hh:479
const CouplingStencilType & couplingStencil(Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, 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_cvfe.hh:495
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_cvfe.hh:235
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the density at a given position.
Definition: couplingmanager_cvfe.hh:261
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given sub control volume face.
Definition: couplingmanager_cvfe.hh:313
VelocityVector elementVelocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry) const
Returns the velocity at the element center.
Definition: couplingmanager_cvfe.hh:414
Scalar density(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the density at a given sub control volume.
Definition: couplingmanager_cvfe.hh:249
Scalar effectiveViscosity(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the effective viscosity at a given position.
Definition: couplingmanager_cvfe.hh:339
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given sub control volume.
Definition: couplingmanager_cvfe.hh:206
void assembleMultithreaded(Dune::index_constant< i > domainId, AssembleElementFunc &&assembleElement) const
Execute assembly kernel in parallel.
Definition: couplingmanager_cvfe.hh:615
Scalar pressure(const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
Returns the pressure at a given sub control volume face.
Definition: couplingmanager_cvfe.hh:192
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_cvfe.hh:157
void computeColorsForAssembly()
Compute colors for multithreaded assembly.
Definition: couplingmanager_cvfe.hh:593
static constexpr auto freeFlowMassIndex
Definition: couplingmanager_cvfe.hh:58
VelocityVector velocity(const FVElementGeometry< freeFlowMassIndex > &fvGeometry, const IpData &ipData, const bool considerPreviousTimeStep=false) const
Returns the velocity at an interpolation point.
Definition: couplingmanager_cvfe.hh:436
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_cvfe.hh:459
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:332
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:297
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:319
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:348
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:229
typename Traits::template TupleOfSharedPtr< SubSolutionVector > SolutionVectorStorage
the type in which the solution vector is stored in the manager
Definition: multidomain/couplingmanager.hh:60
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Classes representing interpolation point data for control-volume finite element schemes.
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
PrimaryVariables evalSolutionAtLocalPos(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::LocalCoordinate &localPos, bool ignoreState=false)
Interpolates a given cvfe element solution at a given local position. Uses the finite element cache o...
Definition: evalsolution.hh:173
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_cvfe.hh:530
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
Class representing dofs on elements for control-volume finite element schemes.
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.
constexpr FCDiamond fcdiamond
Definition: method.hh:162
constexpr CCTpfa cctpfa
Definition: method.hh:154
constexpr Box box
Definition: method.hh:156
constexpr PQ1Bubble pq1bubble
Definition: method.hh:158
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:241
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50
Parallel for loop (multithreading)
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
std::true_type type
Definition: couplingmanager_cvfe.hh:779
Definition: couplingmanager_cvfe.hh:774