version 3.11-dev
Loading...
Searching...
No Matches
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-FileCopyrightText: 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
94template<class Problem,
95 class FVElementGeometry,
96 class ElementVolumeVariables,
97 class IpData>
99{
100 using Element = typename FVElementGeometry::Element;
101 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
102
103 static constexpr int dim = FVElementGeometry::GridGeometry::GridView::dimension;
104 static constexpr int dimWorld = FVElementGeometry::GridGeometry::GridView::dimensionworld;
105
106 using Tensor = Dune::FieldMatrix<typename GlobalPosition::value_type, dim, dimWorld>;
107
108public:
109
112 const Problem& problem,
113 const FVElementGeometry& fvGeometry,
114 const ElementVolumeVariables& elemVolVars,
115 const IpData& ipData
116 )
117 : problem_(problem)
118 , fvGeometry_(fvGeometry)
119 , elemVolVars_(elemVolVars)
120 , ipData_(ipData)
121 , velocity_(0.0)
122 , gradVelocity_(0.0)
123 {
124 calculateVelocity_();
125 calculateGradVelocity();
126 }
127
128 const Problem& problem() const
129 { return problem_; }
130
131 const Element& element() const
132 { return fvGeometry_.element(); }
133
134 const FVElementGeometry& fvGeometry() const
135 { return fvGeometry_; }
136
137 const ElementVolumeVariables& elemVolVars() const
138 { return elemVolVars_; }
139
140 const GlobalPosition& velocity() const
141 { return velocity_; }
142
143 const Tensor& gradVelocity() const
144 { return gradVelocity_; }
145
146private:
147 void calculateVelocity_()
148 {
149 const auto& shapeValues = ipData_.shapeValues();
150 for (const auto& localDof : localDofs(fvGeometry_))
151 velocity_.axpy(shapeValues[localDof.index()][0], elemVolVars_[localDof.index()].velocity());
152 }
153
154 void calculateGradVelocity()
155 {
156 for (const auto& localDof : localDofs(fvGeometry_))
157 {
158 const auto& volVars = elemVolVars_[localDof.index()];
159 for (int dir = 0; dir < dim; ++dir)
160 gradVelocity_[dir].axpy(volVars.velocity(dir), ipData_.gradN(localDof.index()));
161 }
162 }
163
164 const Problem& problem_;
165 const FVElementGeometry& fvGeometry_;
166 const ElementVolumeVariables& elemVolVars_;
167 const IpData& ipData_;
168 GlobalPosition velocity_;
169 Tensor gradVelocity_;
170};
171
176template<class GridGeometry, class NumEqVector>
178{
179 using GridView = typename GridGeometry::GridView;
180 using Element = typename GridView::template Codim<0>::Entity;
181 using Scalar = typename NumEqVector::value_type;
182
183 using Extrusion = Extrusion_t<GridGeometry>;
184
185 static constexpr int dim = GridView::dimension;
186 static constexpr int dimWorld = GridView::dimensionworld;
187
188 using Tensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
189 static_assert(NumEqVector::dimension == dimWorld, "Wrong dimension of velocity vector");
190
191public:
195 template<class Context>
196 NumEqVector advectiveMomentumFlux(const Context& context) const
197 {
198 if (!context.problem().enableInertiaTerms())
199 return NumEqVector(0.0);
200
201 const auto& fvGeometry = context.fvGeometry();
202 const auto& elemVolVars = context.elemVolVars();
203 const auto& scvf = context.scvFace();
204 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
205 const auto& shapeValues = fluxVarCache.shapeValues();
206
207 // interpolate velocity at scvf
208 NumEqVector v(0.0);
209 for (const auto& localDof : localDofs(fvGeometry))
210 v.axpy(shapeValues[localDof.index()][0], elemVolVars[localDof.index()].velocity());
211
212 // get density from the problem
213 const Scalar density = context.problem().density(context.element(), context.fvGeometry(), fluxVarCache.ipData());
214
215 const auto vn = v*scvf.unitOuterNormal();
216 const auto& insideVolVars = elemVolVars[fvGeometry.scv(scvf.insideScvIdx())];
217 const auto& outsideVolVars = elemVolVars[fvGeometry.scv(scvf.outsideScvIdx())];
218 const auto upwindVelocity = vn > 0 ? insideVolVars.velocity() : outsideVolVars.velocity();
219 const auto downwindVelocity = vn > 0 ? outsideVolVars.velocity() : insideVolVars.velocity();
220 static const auto upwindWeight = getParamFromGroup<Scalar>(context.problem().paramGroup(), "Flux.UpwindWeight");
221 const auto advectiveTermIntegrand = density*vn * (upwindWeight * upwindVelocity + (1.0-upwindWeight)*downwindVelocity);
222
223 return advectiveTermIntegrand * Extrusion::area(fvGeometry, scvf) * insideVolVars.extrusionFactor();
224 }
225
229 template<class Context>
230 NumEqVector diffusiveMomentumFlux(const Context& context) const
231 {
232 const auto& element = context.element();
233 const auto& fvGeometry = context.fvGeometry();
234 const auto& elemVolVars = context.elemVolVars();
235 const auto& scvf = context.scvFace();
236 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
237
238 // interpolate velocity gradient at scvf
239 Tensor gradV(0.0);
240 for (const auto& localDof : localDofs(fvGeometry))
241 {
242 const auto& volVars = elemVolVars[localDof];
243 for (int dir = 0; dir < dim; ++dir)
244 gradV[dir].axpy(volVars.velocity(dir), fluxVarCache.gradN(localDof.index()));
245 }
246
247 // get viscosity from the problem
248 const auto mu = context.problem().effectiveViscosity(element, fvGeometry, fluxVarCache.ipData());
249
250 static const bool enableUnsymmetrizedVelocityGradient
251 = getParamFromGroup<bool>(context.problem().paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
252
253 // compute -mu*gradV*n*dA
254 NumEqVector diffusiveFlux = enableUnsymmetrizedVelocityGradient ?
255 mv(gradV, scvf.unitOuterNormal())
256 : mv(gradV + getTransposed(gradV),scvf.unitOuterNormal());
257
258 diffusiveFlux *= -mu;
259
260 static const bool enableDilatationTerm = getParamFromGroup<bool>(context.problem().paramGroup(), "FreeFlow.EnableDilatationTerm", false);
261 if (enableDilatationTerm)
262 diffusiveFlux += 2.0/3.0 * mu * trace(gradV) * scvf.unitOuterNormal();
263
264 diffusiveFlux *= Extrusion::area(fvGeometry, scvf) * elemVolVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
265 return diffusiveFlux;
266 }
267
268 template<class Context>
269 NumEqVector pressureContribution(const Context& context) const
270 {
271 const auto& element = context.element();
272 const auto& fvGeometry = context.fvGeometry();
273 const auto& elemVolVars = context.elemVolVars();
274 const auto& scvf = context.scvFace();
275 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
276
277 // The pressure force needs to take the extruded scvf area into account
278 const auto pressure = context.problem().pressure(element, fvGeometry, fluxVarCache.ipData());
279
280 // The pressure contribution calculated above might have a much larger numerical value compared to the viscous or inertial forces.
281 // This may lead to numerical inaccuracies due to loss of significance (cancellation) for the final residual value.
282 // In the end, we are only interested in a pressure gradient between the two relevant faces so we can
283 // subtract a constant reference value from the actual pressure contribution.
284 const auto referencePressure = context.problem().referencePressure();
285
286 NumEqVector pn(scvf.unitOuterNormal());
287 pn *= (pressure-referencePressure)*Extrusion::area(fvGeometry, scvf)*elemVolVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
288
289 return pn;
290 }
291};
292
297template<class GridGeometry, class NumEqVector>
299{
300 using GridView = typename GridGeometry::GridView;
301 using Element = typename GridView::template Codim<0>::Entity;
302 using Scalar = typename NumEqVector::value_type;
303
304 using Extrusion = Extrusion_t<GridGeometry>;
305
306 static constexpr int dim = GridView::dimension;
307 static constexpr int dimWorld = GridView::dimensionworld;
308
309 using Tensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
310 static_assert(NumEqVector::dimension == dimWorld, "Wrong dimension of velocity vector");
311
312public:
316 template<class Problem, class FVElementGeometry, class ElementVariables,
317 class SubControlVolumeFace, class VelocityVector>
319 const FVElementGeometry& fvGeometry,
320 const ElementVariables& elemVars,
321 const SubControlVolumeFace& scvf,
322 const VelocityVector& integratedVelocity) const
323 {
324 if (!problem.enableInertiaTerms())
325 return NumEqVector(0.0);
326
327 // get density from the problem
328 const Scalar density = problem.density(fvGeometry.element(), fvGeometry, ipData(fvGeometry, scvf));
329
330 const auto vn_integral = integratedVelocity*scvf.unitOuterNormal();
331 const auto& insideVolVars = elemVars[fvGeometry.scv(scvf.insideScvIdx())];
332 const auto& outsideVolVars = elemVars[fvGeometry.scv(scvf.outsideScvIdx())];
333 const auto upwindVelocity = vn_integral > 0 ? insideVolVars.velocity() : outsideVolVars.velocity();
334 const auto downwindVelocity = vn_integral > 0 ? outsideVolVars.velocity() : insideVolVars.velocity();
335 static const auto upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
336 const auto advectiveFlux = density*vn_integral * (upwindWeight * upwindVelocity + (1.0-upwindWeight)*downwindVelocity);
337
338 return advectiveFlux;
339 }
340
344 template<class Context, class IpData>
345 NumEqVector diffusiveMomentumFluxIntegrand(const Context& context, const IpData& ipData) const
346 {
347 const auto& element = context.element();
348 const auto& fvGeometry = context.fvGeometry();
349
350 // get viscosity from the problem
351 const auto mu = context.problem().effectiveViscosity(element, fvGeometry, ipData);
352
353 static const bool enableUnsymmetrizedVelocityGradient
354 = getParamFromGroup<bool>(context.problem().paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
355
356 // compute -mu*gradV*n*dA
357 const auto& gradV = context.gradVelocity();
358 NumEqVector diffusiveFluxIntegrand = enableUnsymmetrizedVelocityGradient ?
359 mv(gradV, ipData.unitOuterNormal())
360 : mv(gradV + getTransposed(gradV), ipData.unitOuterNormal());
361
362 diffusiveFluxIntegrand *= -mu;
363
364 static const bool enableDilatationTerm = getParamFromGroup<bool>(context.problem().paramGroup(), "FreeFlow.EnableDilatationTerm", false);
365 if (enableDilatationTerm)
366 diffusiveFluxIntegrand += 2.0/3.0 * mu * trace(gradV) * ipData.unitOuterNormal();
367
368 return diffusiveFluxIntegrand;
369 }
370
371 template<class Context, class IpData>
372 NumEqVector pressureFluxIntegrand(const Context& context, const IpData& ipData) const
373 {
374 const auto& element = context.element();
375 const auto& fvGeometry = context.fvGeometry();
376
377 // The pressure force integrand
378 const auto pressure = context.problem().pressure(element, fvGeometry, ipData);
379
380 // The pressure contribution calculated above might have a much larger numerical value compared to the viscous or inertial forces.
381 // This may lead to numerical inaccuracies due to loss of significance (cancellation) for the final residual value.
382 // Subtracting a constant reference value from the actual pressure contribution helps improving cancellation issues.
383 const auto referencePressure = context.problem().referencePressure();
384
385 NumEqVector pn(ipData.unitOuterNormal());
386 pn *= (pressure-referencePressure);
387
388 return pn;
389 }
390};
391
392} // end namespace Dumux
393
394#endif
The flux variables class for the Navier-Stokes model using control-volume finite element schemes.
Definition flux.hh:178
NumEqVector advectiveMomentumFlux(const Context &context) const
Returns the advective momentum flux.
Definition flux.hh:196
NumEqVector diffusiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition flux.hh:230
NumEqVector pressureContribution(const Context &context) const
Definition flux.hh:269
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
The flux function class for the Navier-Stokes model using control-volume finite element schemes.
Definition flux.hh:299
NumEqVector pressureFluxIntegrand(const Context &context, const IpData &ipData) const
Definition flux.hh:372
NumEqVector diffusiveMomentumFluxIntegrand(const Context &context, const IpData &ipData) const
Returns the diffusive momentum flux due to viscous forces.
Definition flux.hh:345
NumEqVector advectiveMomentumFluxIntegral(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const SubControlVolumeFace &scvf, const VelocityVector &integratedVelocity) const
Returns the advective momentum flux contribution for a given integrated velocity at the face.
Definition flux.hh:318
const FVElementGeometry & fvGeometry() const
Definition flux.hh:134
const Element & element() const
Definition flux.hh:131
const Problem & problem() const
Definition flux.hh:128
const GlobalPosition & velocity() const
Definition flux.hh:140
const Tensor & gradVelocity() const
Definition flux.hh:143
NavierStokesMomentumFluxFunctionContext(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const IpData &ipData)
Initialize the flux variables storing some temporary pointers.
Definition flux.hh:111
const ElementVolumeVariables & elemVolVars() const
Definition flux.hh:137
Some exceptions thrown in DuMux.
Helper classes to compute the integration elements.
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
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:829
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition math.hh:800
Dune::FieldMatrix< Scalar, n, m > getTransposed(const Dune::FieldMatrix< Scalar, m, n > &M)
Transpose a FieldMatrix.
Definition math.hh:712
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
Define some often used mathematical functions.
Definition adapt.hh:17
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition localdof.hh:50
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:236
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.