3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porenetwork/1p/fluxvariablescache.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PNM_1P_FLUXVARIABLESCACHE_HH
25#define DUMUX_PNM_1P_FLUXVARIABLESCACHE_HH
26
28
29namespace Dumux::PoreNetwork {
30
36template<class AdvectionType>
38{
39 using Scalar = typename AdvectionType::Scalar;
40public:
41
42 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables>
43 void update(const Problem& problem,
44 const Element& element,
45 const FVElementGeometry& fvGeometry,
46 const ElementVolumeVariables& elemVolVars,
47 const typename FVElementGeometry::SubControlVolumeFace& scvf)
48 {
49 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
50 throatCrossSectionShape_ = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
51 throatShapeFactor_ = fvGeometry.gridGeometry().throatShapeFactor(eIdx);
52 throatCrossSectionalArea_ = problem.spatialParams().throatCrossSectionalArea(element, elemVolVars);
53 throatLength_ = problem.spatialParams().throatLength(element, elemVolVars);
54 throatInscribedRadius_ = problem.spatialParams().throatInscribedRadius(element, elemVolVars);
55 poreToPoreDistance_ = element.geometry().volume();
56
57 cache_.fill(problem, element, fvGeometry, scvf, elemVolVars, *this, /*phaseIdx*/0);
58 transmissibility_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, *this, /*phaseIdx*/0);
59 }
60
65 { return throatCrossSectionShape_; }
66
70 Scalar throatShapeFactor() const
71 { return throatShapeFactor_; }
72
76 Scalar transmissibility(const int phaseIdx = 0) const
77 { return transmissibility_; }
78
82 Scalar throatCrossSectionalArea(const int phaseIdx = 0) const
83 { return throatCrossSectionalArea_; }
84
88 Scalar throatLength() const
89 { return throatLength_; }
90
94 Scalar throatInscribedRadius() const
95 { return throatInscribedRadius_; }
96
100 Scalar poreToPoreDistance() const
101 { return poreToPoreDistance_; }
102
106 const auto& singlePhaseFlowVariables() const
107 { return cache_; }
108
109private:
110 Throat::Shape throatCrossSectionShape_;
111 Scalar throatShapeFactor_;
112 Scalar transmissibility_;
113 Scalar throatCrossSectionalArea_;
114 Scalar throatLength_;
115 Scalar throatInscribedRadius_;
116 Scalar poreToPoreDistance_;
117
118 typename AdvectionType::Transmissibility::SinglePhaseCache cache_;
119};
120
121} // end namespace Dumux::PoreNetwork
122
123#endif
This file contains functions related to calculate pore-throat properties.
Definition: discretization/porenetwork/fvelementgeometry.hh:33
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:38
Flux variables cache for the single-phase-flow PNM Store data required for flux calculation.
Definition: porenetwork/1p/fluxvariablescache.hh:38
const auto & singlePhaseFlowVariables() const
Returns the throats's cached flow variables for single-phase flow.
Definition: porenetwork/1p/fluxvariablescache.hh:106
Scalar transmissibility(const int phaseIdx=0) const
Returns the throats's transmissibility.
Definition: porenetwork/1p/fluxvariablescache.hh:76
Scalar poreToPoreDistance() const
Returns the throats's pore-to-pore-center distance.
Definition: porenetwork/1p/fluxvariablescache.hh:100
Scalar throatShapeFactor() const
Returns the throats's shape factor.
Definition: porenetwork/1p/fluxvariablescache.hh:70
Scalar throatInscribedRadius() const
Returns the throats's inscribed radius.
Definition: porenetwork/1p/fluxvariablescache.hh:94
Throat::Shape throatCrossSectionShape() const
Returns the throats's cross-sectional shape.
Definition: porenetwork/1p/fluxvariablescache.hh:64
Scalar throatLength() const
Returns the throats's length.
Definition: porenetwork/1p/fluxvariablescache.hh:88
Scalar throatCrossSectionalArea(const int phaseIdx=0) const
Returns the throats's cross-sectional area.
Definition: porenetwork/1p/fluxvariablescache.hh:82
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf)
Definition: porenetwork/1p/fluxvariablescache.hh:43