12#ifndef DUMUX_PNM_1P_FLUXVARIABLESCACHE_HH 
   13#define DUMUX_PNM_1P_FLUXVARIABLESCACHE_HH 
   24template<
class AdvectionType>
 
   27    using Scalar = 
typename AdvectionType::Scalar;
 
   32    template<
class Problem, 
class Element, 
class FVElementGeometry, 
class ElementVolumeVariables>
 
   34                const Element& element,
 
   35                const FVElementGeometry& fvGeometry,
 
   36                const ElementVolumeVariables& elemVolVars,
 
   37                const typename FVElementGeometry::SubControlVolumeFace& scvf)
 
   39        const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
 
   40        throatCrossSectionShape_ = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
 
   41        throatShapeFactor_ = fvGeometry.gridGeometry().throatShapeFactor(eIdx);
 
   42        throatCrossSectionalArea_ = problem.spatialParams().throatCrossSectionalArea(element, elemVolVars);
 
   43        throatLength_ = problem.spatialParams().throatLength(element, elemVolVars);
 
   44        throatInscribedRadius_ = problem.spatialParams().throatInscribedRadius(element, elemVolVars);
 
   45        poreToPoreDistance_ = element.geometry().volume();
 
   47        cache_.fill(problem, element, fvGeometry, scvf, elemVolVars, *
this, 0);
 
   48        transmissibility_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, *
this, 0);
 
 
   55    { 
return throatCrossSectionShape_; }
 
 
   61    { 
return throatShapeFactor_; }
 
 
   67    { 
return transmissibility_; }
 
 
   73    { 
return throatCrossSectionalArea_; }
 
 
   79    { 
return throatLength_; }
 
 
   85    { 
return throatInscribedRadius_; }
 
 
   91    { 
return poreToPoreDistance_; }
 
 
  101    Scalar throatShapeFactor_;
 
  102    Scalar transmissibility_;
 
  103    Scalar throatCrossSectionalArea_;
 
  104    Scalar throatLength_;
 
  105    Scalar throatInscribedRadius_;
 
  106    Scalar poreToPoreDistance_;
 
  108    typename AdvectionType::Transmissibility::SinglePhaseCache cache_;
 
 
Flux variables cache for the single-phase-flow PNM Store data required for flux calculation.
Definition porenetwork/1p/fluxvariablescache.hh:26
const auto & singlePhaseFlowVariables() const
Returns the throats's cached flow variables for single-phase flow.
Definition porenetwork/1p/fluxvariablescache.hh:96
Scalar transmissibility(const int phaseIdx=0) const
Returns the throats's transmissibility.
Definition porenetwork/1p/fluxvariablescache.hh:66
Scalar poreToPoreDistance() const
Returns the throats's pore-to-pore-center distance.
Definition porenetwork/1p/fluxvariablescache.hh:90
Scalar throatShapeFactor() const
Returns the throats's shape factor.
Definition porenetwork/1p/fluxvariablescache.hh:60
Scalar throatInscribedRadius() const
Returns the throats's inscribed radius.
Definition porenetwork/1p/fluxvariablescache.hh:84
Throat::Shape throatCrossSectionShape() const
Returns the throats's cross-sectional shape.
Definition porenetwork/1p/fluxvariablescache.hh:54
static bool constexpr isSolDependent
Definition porenetwork/1p/fluxvariablescache.hh:30
Scalar throatLength() const
Returns the throats's length.
Definition porenetwork/1p/fluxvariablescache.hh:78
Scalar throatCrossSectionalArea(const int phaseIdx=0) const
Returns the throats's cross-sectional area.
Definition porenetwork/1p/fluxvariablescache.hh:72
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf)
Definition porenetwork/1p/fluxvariablescache.hh:33
Shape
Collection of different pore-throat shapes.
Definition throatproperties.hh:26
Definition discretization/porenetwork/fvelementgeometry.hh:24
This file contains functions related to calculate pore-throat properties.