version 3.8
flux.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_NAVIERSTOKES_MOMENTUM_CVFE_FLUXVARIABLES_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FLUXVARIABLES_HH
14
15#include <dune/common/fmatrix.hh>
16
17#include <dumux/common/math.hh>
20
22
23namespace Dumux {
24
34template<class Problem,
35 class FVElementGeometry,
36 class ElementVolumeVariables,
37 class ElementFluxVariablesCache>
39{
40 using Element = typename FVElementGeometry::Element;
41 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
42public:
43
46 const Problem& problem,
47 const FVElementGeometry& fvGeometry,
48 const ElementVolumeVariables& elemVolVars,
49 const ElementFluxVariablesCache& elemFluxVarsCache,
50 const SubControlVolumeFace& scvf
51 )
52 : problem_(problem)
53 , fvGeometry_(fvGeometry)
54 , elemVolVars_(elemVolVars)
55 , elemFluxVarsCache_(elemFluxVarsCache)
56 , scvf_(scvf)
57 {}
58
59 const Problem& problem() const
60 { return problem_; }
61
62 const Element& element() const
63 { return fvGeometry_.element(); }
64
65 const SubControlVolumeFace& scvFace() const
66 { return scvf_; }
67
68 const FVElementGeometry& fvGeometry() const
69 { return fvGeometry_; }
70
71 const ElementVolumeVariables& elemVolVars() const
72 { return elemVolVars_; }
73
74 const ElementFluxVariablesCache& elemFluxVarsCache() const
75 { return elemFluxVarsCache_; }
76
77private:
78 const Problem& problem_;
79 const FVElementGeometry& fvGeometry_;
80 const ElementVolumeVariables& elemVolVars_;
81 const ElementFluxVariablesCache& elemFluxVarsCache_;
82 const SubControlVolumeFace& scvf_;
83};
84
89template<class GridGeometry, class NumEqVector>
91{
92 using GridView = typename GridGeometry::GridView;
93 using Element = typename GridView::template Codim<0>::Entity;
94 using Scalar = typename NumEqVector::value_type;
95
96 using Extrusion = Extrusion_t<GridGeometry>;
97
98 static constexpr int dim = GridView::dimension;
99 static constexpr int dimWorld = GridView::dimensionworld;
100
101 using Tensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
102 static_assert(NumEqVector::dimension == dimWorld, "Wrong dimension of velocity vector");
103
104public:
108 template<class Context>
109 NumEqVector advectiveMomentumFlux(const Context& context) const
110 {
111 if (!context.problem().enableInertiaTerms())
112 return NumEqVector(0.0);
113
114 const auto& fvGeometry = context.fvGeometry();
115 const auto& elemVolVars = context.elemVolVars();
116 const auto& scvf = context.scvFace();
117 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
118 const auto& shapeValues = fluxVarCache.shapeValues();
119
120 // interpolate velocity at scvf
121 NumEqVector v(0.0);
122 for (const auto& scv : scvs(fvGeometry))
123 v.axpy(shapeValues[scv.indexInElement()][0], elemVolVars[scv].velocity());
124
125 // get density from the problem
126 const Scalar density = context.problem().density(context.element(), context.fvGeometry(), scvf);
127
128 const auto vn = v*scvf.unitOuterNormal();
129 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
130 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
131 const auto upwindVelocity = vn > 0 ? insideVolVars.velocity() : outsideVolVars.velocity();
132 const auto downwindVelocity = vn > 0 ? outsideVolVars.velocity() : insideVolVars.velocity();
133 static const auto upwindWeight = getParamFromGroup<Scalar>(context.problem().paramGroup(), "Flux.UpwindWeight");
134 const auto advectiveTermIntegrand = density*vn * (upwindWeight * upwindVelocity + (1.0-upwindWeight)*downwindVelocity);
135
136 return advectiveTermIntegrand * Extrusion::area(fvGeometry, scvf) * insideVolVars.extrusionFactor();
137 }
138
142 template<class Context>
143 NumEqVector diffusiveMomentumFlux(const Context& context) const
144 {
145 const auto& element = context.element();
146 const auto& fvGeometry = context.fvGeometry();
147 const auto& elemVolVars = context.elemVolVars();
148 const auto& scvf = context.scvFace();
149 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
150
151 // interpolate velocity gradient at scvf
152 Tensor gradV(0.0);
153 for (const auto& scv : scvs(fvGeometry))
154 {
155 const auto& volVars = elemVolVars[scv];
156 for (int dir = 0; dir < dim; ++dir)
157 gradV[dir].axpy(volVars.velocity(dir), fluxVarCache.gradN(scv.indexInElement()));
158 }
159
160 // get viscosity from the problem
161 const auto mu = context.problem().effectiveViscosity(element, fvGeometry, scvf);
162
163 static const bool enableUnsymmetrizedVelocityGradient
164 = getParamFromGroup<bool>(context.problem().paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
165
166 // compute -mu*gradV*n*dA
167 NumEqVector diffusiveFlux = enableUnsymmetrizedVelocityGradient ?
168 mv(gradV, scvf.unitOuterNormal())
169 : mv(gradV + getTransposed(gradV),scvf.unitOuterNormal());
170
171 diffusiveFlux *= -mu;
172
173 static const bool enableDilatationTerm = getParamFromGroup<bool>(context.problem().paramGroup(), "FreeFlow.EnableDilatationTerm", false);
174 if (enableDilatationTerm)
175 diffusiveFlux += 2.0/3.0 * mu * trace(gradV) * scvf.unitOuterNormal();
176
177 diffusiveFlux *= Extrusion::area(fvGeometry, scvf) * elemVolVars[scvf.insideScvIdx()].extrusionFactor();
178 return diffusiveFlux;
179 }
180
181 template<class Context>
182 NumEqVector pressureContribution(const Context& context) const
183 {
184 const auto& element = context.element();
185 const auto& fvGeometry = context.fvGeometry();
186 const auto& elemVolVars = context.elemVolVars();
187 const auto& scvf = context.scvFace();
188
189 // The pressure force needs to take the extruded scvf area into account
190 const auto pressure = context.problem().pressure(element, fvGeometry, scvf);
191
192 // The pressure contribution calculated above might have a much larger numerical value compared to the viscous or inertial forces.
193 // This may lead to numerical inaccuracies due to loss of significance (cancellation) for the final residual value.
194 // In the end, we are only interested in a pressure gradient between the two relevant faces so we can
195 // subtract a constant reference value from the actual pressure contribution.
196 const auto referencePressure = context.problem().referencePressure();
197
198 NumEqVector pn(scvf.unitOuterNormal());
199 pn *= (pressure-referencePressure)*Extrusion::area(fvGeometry, scvf)*elemVolVars[scvf.insideScvIdx()].extrusionFactor();
200
201 return pn;
202 }
203};
204
205} // end namespace Dumux
206
207#endif
The flux variables class for the Navier-Stokes model using control-volume finite element schemes.
Definition: flux.hh:91
NumEqVector advectiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:109
NumEqVector diffusiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:143
NumEqVector pressureContribution(const Context &context) const
Definition: flux.hh:182
Context for computing fluxes.
Definition: flux.hh:39
const Element & element() const
Definition: flux.hh:62
const FVElementGeometry & fvGeometry() const
Definition: flux.hh:68
NavierStokesMomentumFluxContext(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf)
Initialize the flux variables storing some temporary pointers.
Definition: flux.hh:45
const SubControlVolumeFace & scvFace() const
Definition: flux.hh:65
const ElementVolumeVariables & elemVolVars() const
Definition: flux.hh:71
const ElementFluxVariablesCache & elemFluxVarsCache() const
Definition: flux.hh:74
const Problem & problem() const
Definition: flux.hh:59
Some exceptions thrown in DuMux
Helper classes to compute the integration elements.
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
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:771
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition: math.hh:683
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
Define some often used mathematical functions.
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.