GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/flux/box/forchheimerslaw.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 29 31 93.5%
Functions: 1 2 50.0%
Branches: 21 36 58.3%

Line Branch Exec Source
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 // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5 // SPDX-License-Identifier: GPL-3.0-or-later
6 //
7 /*!
8 * \file
9 * \ingroup BoxFlux
10 * \brief Forchheimers's law for the box method
11 */
12 #ifndef DUMUX_DISCRETIZATION_BOX_FORCHHEIMERS_LAW_HH
13 #define DUMUX_DISCRETIZATION_BOX_FORCHHEIMERS_LAW_HH
14
15 #include <dune/common/fvector.hh>
16 #include <dune/common/fmatrix.hh>
17
18 #include <dumux/common/math.hh>
19 #include <dumux/common/parameters.hh>
20 #include <dumux/common/properties.hh>
21 #include <dumux/common/typetraits/typetraits.hh>
22
23 #include <dumux/discretization/method.hh>
24 #include <dumux/discretization/extrusion.hh>
25 #include <dumux/flux/cvfe/darcyslaw.hh>
26 #include <dumux/flux/facetensoraverage.hh>
27
28 namespace Dumux {
29
30 // forward declarations
31 template<class TypeTag, class ForchheimerVelocity, class DiscretizationMethod>
32 class ForchheimersLawImplementation;
33
34 /*!
35 * \ingroup BoxFlux
36 * \brief Forchheimer's law for box scheme
37 * \note Forchheimer's law is specialized for network and surface grids (i.e. if grid dim < dimWorld)
38 * \tparam Scalar the scalar type for scalar physical quantities
39 * \tparam GridGeometry the grid geometry
40 * \tparam ForchheimerVelocity class for the calculation of the Forchheimer velocity
41 * \tparam isNetwork whether we are computing on a network grid embedded in a higher world dimension
42 */
43 template<class Scalar, class GridGeometry, class ForchheimerVelocity>
44 class BoxForchheimersLaw;
45
46 /*!
47 * \ingroup BoxFlux
48 * \brief Forchheimer's law for box scheme
49 */
50 template <class TypeTag, class ForchheimerVelocity>
51 class ForchheimersLawImplementation<TypeTag, ForchheimerVelocity, DiscretizationMethods::Box>
52 : public BoxForchheimersLaw<GetPropType<TypeTag, Properties::Scalar>,
53 GetPropType<TypeTag, Properties::GridGeometry>,
54 ForchheimerVelocity>
55 {};
56
57 /*!
58 * \ingroup BoxFlux
59 * \brief Specialization of the BoxForchheimersLaw
60 */
61 template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
62 class BoxForchheimersLaw
63 {
64 using ThisType = BoxForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity>;
65 using FVElementGeometry = typename GridGeometry::LocalView;
66 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
67 using Extrusion = Extrusion_t<GridGeometry>;
68 using GridView = typename GridGeometry::GridView;
69 using Element = typename GridView::template Codim<0>::Entity;
70
71 using DimWorldVector = typename ForchheimerVelocity::DimWorldVector;
72 using DimWorldMatrix = typename ForchheimerVelocity::DimWorldMatrix;
73
74 using DarcysLaw = CVFEDarcysLaw<ScalarType, GridGeometry>;
75
76 public:
77 //! state the scalar type of the law
78 using Scalar = ScalarType;
79
80 using DiscretizationMethod = DiscretizationMethods::Box;
81 //! state the discretization method this implementation belongs to
82 static constexpr DiscretizationMethod discMethod{};
83
84 /*! \brief Compute the advective flux of a phase across
85 * the given sub-control volume face using the Forchheimer equation.
86 *
87 * The flux is given in N*m, and can be converted
88 * into a volume flux (m^3/s) or mass flux (kg/s) by applying an upwind scheme
89 * for the mobility or the product of density and mobility, respectively.
90 */
91 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
92 1360 static Scalar flux(const Problem& problem,
93 const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const ElementVolumeVariables& elemVolVars,
96 const SubControlVolumeFace& scvf,
97 int phaseIdx,
98 const ElementFluxVarsCache& elemFluxVarsCache)
99 {
100
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1360 times.
1360 const auto& fluxVarCache = elemFluxVarsCache[scvf];
101
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 1360 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1360 times.
2720 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
102
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 1360 times.
✓ Branch 2 taken 1360 times.
✗ Branch 3 not taken.
1360 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
103
1/2
✓ Branch 0 taken 1360 times.
✗ Branch 1 not taken.
1360 const auto& insideVolVars = elemVolVars[insideScv];
104
1/2
✓ Branch 0 taken 1360 times.
✗ Branch 1 not taken.
1360 const auto& outsideVolVars = elemVolVars[outsideScv];
105
106 1360 auto insideK = insideVolVars.permeability();
107 1360 auto outsideK = outsideVolVars.permeability();
108
109 // scale with correct extrusion factor
110 1360 insideK *= insideVolVars.extrusionFactor();
111 1360 outsideK *= outsideVolVars.extrusionFactor();
112
113
2/4
✓ Branch 0 taken 1360 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1360 times.
✗ Branch 3 not taken.
2720 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
114
5/8
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1359 times.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
1360 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
115
116 1360 const auto& shapeValues = fluxVarCache.shapeValues();
117
118 // evaluate gradP - rho*g at integration point
119 1360 DimWorldVector gradP(0.0);
120 1360 Scalar rho(0.0);
121
4/4
✓ Branch 0 taken 5440 times.
✓ Branch 1 taken 1360 times.
✓ Branch 2 taken 5440 times.
✓ Branch 3 taken 1360 times.
13600 for (auto&& scv : scvs(fvGeometry))
122 {
123
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 5440 times.
5440 const auto& volVars = elemVolVars[scv];
124
125
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 5440 times.
5440 if (enableGravity)
126 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
127
128 // the global shape function gradient
129 21760 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
130 }
131
132
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 1360 times.
1360 if (enableGravity)
133 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
134
135 1360 DimWorldVector darcyVelocity = mv(K, gradP);
136 1360 darcyVelocity *= -1;
137
138 5440 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
139 4080 darcyVelocity *= ForchheimerVelocity::UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, darcyVelocity*scvf.unitOuterNormal(), phaseIdx);
140
141 1360 const auto velocity = ForchheimerVelocity::velocity(fvGeometry,
142 elemVolVars,
143 scvf,
144 phaseIdx,
145 calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf),
146 darcyVelocity);
147
148 1360 Scalar flux = velocity * scvf.unitOuterNormal();
149 1360 flux *= Extrusion::area(fvGeometry, scvf);
150 1360 return flux;
151 }
152
153 //! The flux variables cache has to be bound to an element prior to flux calculations
154 //! During the binding, the transmissibility will be computed and stored using the method below.
155 template<class Problem, class ElementVolumeVariables>
156 static Scalar calculateTransmissibility(const Problem& problem,
157 const Element& element,
158 const FVElementGeometry& fvGeometry,
159 const ElementVolumeVariables& elemVolVars,
160 const SubControlVolumeFace& scvf)
161 {
162 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
163 }
164
165 /*! \brief Returns the harmonic mean of \f$\sqrt{K_0}\f$ and \f$\sqrt{K_1}\f$.
166 *
167 * This is a specialization for scalar-valued permeabilities which returns a tensor with identical diagonal entries.
168 */
169 template<class Problem, class ElementVolumeVariables>
170 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf)
173 {
174 1360 return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
175 }
176 };
177
178 } // end namespace Dumux
179
180 #endif
181