version 3.8
box/forchheimerslaw.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
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>
22
27
28namespace Dumux {
29
30// forward declarations
31template<class TypeTag, class ForchheimerVelocity, class DiscretizationMethod>
32class ForchheimersLawImplementation;
33
43template<class Scalar, class GridGeometry, class ForchheimerVelocity>
44class BoxForchheimersLaw;
45
50template <class TypeTag, class ForchheimerVelocity>
51class ForchheimersLawImplementation<TypeTag, ForchheimerVelocity, DiscretizationMethods::Box>
52: public BoxForchheimersLaw<GetPropType<TypeTag, Properties::Scalar>,
53 GetPropType<TypeTag, Properties::GridGeometry>,
54 ForchheimerVelocity>
55{};
56
61template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
63{
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
75
76public:
78 using Scalar = ScalarType;
79
82 static constexpr DiscretizationMethod discMethod{};
83
91 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
92 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 const auto& fluxVarCache = elemFluxVarsCache[scvf];
101 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
102 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
103 const auto& insideVolVars = elemVolVars[insideScv];
104 const auto& outsideVolVars = elemVolVars[outsideScv];
105
106 auto insideK = insideVolVars.permeability();
107 auto outsideK = outsideVolVars.permeability();
108
109 // scale with correct extrusion factor
110 insideK *= insideVolVars.extrusionFactor();
111 outsideK *= outsideVolVars.extrusionFactor();
112
113 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
114 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
115
116 const auto& shapeValues = fluxVarCache.shapeValues();
117
118 // evaluate gradP - rho*g at integration point
119 DimWorldVector gradP(0.0);
120 Scalar rho(0.0);
121 for (auto&& scv : scvs(fvGeometry))
122 {
123 const auto& volVars = elemVolVars[scv];
124
125 if (enableGravity)
126 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
127
128 // the global shape function gradient
129 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
130 }
131
132 if (enableGravity)
133 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
134
135 DimWorldVector darcyVelocity = mv(K, gradP);
136 darcyVelocity *= -1;
137
138 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
139 darcyVelocity *= ForchheimerVelocity::UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, darcyVelocity*scvf.unitOuterNormal(), phaseIdx);
140
141 const auto velocity = ForchheimerVelocity::velocity(fvGeometry,
142 elemVolVars,
143 scvf,
144 phaseIdx,
145 calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf),
146 darcyVelocity);
147
148 Scalar flux = velocity * scvf.unitOuterNormal();
149 flux *= Extrusion::area(fvGeometry, scvf);
150 return flux;
151 }
152
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
169 template<class Problem, class ElementVolumeVariables>
170 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf)
173 {
174 return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
175 }
176};
177
178} // end namespace Dumux
179
180#endif
Forchheimer's law for box scheme.
Definition: box/forchheimerslaw.hh:63
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: box/forchheimerslaw.hh:156
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: box/forchheimerslaw.hh:170
static constexpr DiscretizationMethod discMethod
state the discretization method this implementation belongs to
Definition: box/forchheimerslaw.hh:82
ScalarType Scalar
state the scalar type of the law
Definition: box/forchheimerslaw.hh:78
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux of a phase across the given sub-control volume face using the Forchheimer ...
Definition: box/forchheimerslaw.hh:92
Darcy's law for control-volume finite element schemes.
Definition: flux/cvfe/darcyslaw.hh:36
Forchheimer's law For a detailed description see dumux/flow/forchheimerslaw.hh.
Definition: forchheimervelocity.hh:39
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:126
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:50
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Definition: forchheimervelocity.hh:55
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:51
forward declare
Definition: forchheimerslaw_fwd.hh:27
Defines all properties used in Dumux.
Type traits.
Helper classes to compute the integration elements.
A free function to average a Tensor at an interface.
Specialization of Darcy's Law for control-volume finite element schemes.
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:800
Define some often used mathematical functions.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: method.hh:46