12#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FLUXVARIABLES_HH
13#define DUMUX_NAVIERSTOKES_MOMENTUM_CVFE_FLUXVARIABLES_HH
15#include <dune/common/fmatrix.hh>
34template<
class Problem,
35 class FVElementGeometry,
36 class ElementVolumeVariables,
37 class ElementFluxVariablesCache>
40 using Element =
typename FVElementGeometry::Element;
41 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
50 const SubControlVolumeFace& scvf
63 {
return fvGeometry_.element(); }
65 const SubControlVolumeFace&
scvFace()
const
69 {
return fvGeometry_; }
72 {
return elemVolVars_; }
75 {
return elemFluxVarsCache_; }
78 const Problem& problem_;
79 const FVElementGeometry& fvGeometry_;
80 const ElementVolumeVariables& elemVolVars_;
81 const ElementFluxVariablesCache& elemFluxVarsCache_;
82 const SubControlVolumeFace& scvf_;
94template<
class Problem,
95 class FVElementGeometry,
96 class ElementVolumeVariables,
100 using Element =
typename FVElementGeometry::Element;
101 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
103 static constexpr int dim = FVElementGeometry::GridGeometry::GridView::dimension;
104 static constexpr int dimWorld = FVElementGeometry::GridGeometry::GridView::dimensionworld;
106 using Tensor = Dune::FieldMatrix<typename GlobalPosition::value_type, dim, dimWorld>;
124 calculateVelocity_();
125 calculateGradVelocity();
132 {
return fvGeometry_.element(); }
135 {
return fvGeometry_; }
138 {
return elemVolVars_; }
141 {
return velocity_; }
144 {
return gradVelocity_; }
147 void calculateVelocity_()
149 const auto& shapeValues = ipData_.shapeValues();
150 for (
const auto& localDof :
localDofs(fvGeometry_))
151 velocity_.axpy(shapeValues[localDof.index()][0], elemVolVars_[localDof.index()].velocity());
154 void calculateGradVelocity()
156 for (
const auto& localDof :
localDofs(fvGeometry_))
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()));
164 const Problem& problem_;
165 const FVElementGeometry& fvGeometry_;
166 const ElementVolumeVariables& elemVolVars_;
167 const IpData& ipData_;
168 GlobalPosition velocity_;
169 Tensor gradVelocity_;
176template<
class Gr
idGeometry,
class NumEqVector>
179 using GridView =
typename GridGeometry::GridView;
180 using Element =
typename GridView::template Codim<0>::Entity;
181 using Scalar =
typename NumEqVector::value_type;
185 static constexpr int dim = GridView::dimension;
186 static constexpr int dimWorld = GridView::dimensionworld;
188 using Tensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
189 static_assert(NumEqVector::dimension == dimWorld,
"Wrong dimension of velocity vector");
195 template<
class Context>
198 if (!context.problem().enableInertiaTerms())
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();
209 for (
const auto& localDof :
localDofs(fvGeometry))
210 v.axpy(shapeValues[localDof.index()][0], elemVolVars[localDof.index()].velocity());
213 const Scalar density = context.problem().density(context.element(), context.fvGeometry(), fluxVarCache.ipData());
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();
221 const auto advectiveTermIntegrand = density*vn * (upwindWeight * upwindVelocity + (1.0-upwindWeight)*downwindVelocity);
223 return advectiveTermIntegrand * Extrusion::area(fvGeometry, scvf) * insideVolVars.extrusionFactor();
229 template<
class Context>
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];
240 for (
const auto& localDof :
localDofs(fvGeometry))
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()));
248 const auto mu = context.problem().effectiveViscosity(element, fvGeometry, fluxVarCache.ipData());
250 static const bool enableUnsymmetrizedVelocityGradient
251 =
getParamFromGroup<bool>(context.problem().paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
254 NumEqVector diffusiveFlux = enableUnsymmetrizedVelocityGradient ?
255 mv(gradV, scvf.unitOuterNormal())
258 diffusiveFlux *= -mu;
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();
264 diffusiveFlux *= Extrusion::area(fvGeometry, scvf) * elemVolVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
265 return diffusiveFlux;
268 template<
class Context>
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];
278 const auto pressure = context.problem().pressure(element, fvGeometry, fluxVarCache.ipData());
284 const auto referencePressure = context.problem().referencePressure();
287 pn *= (pressure-referencePressure)*Extrusion::area(fvGeometry, scvf)*elemVolVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
297template<
class Gr
idGeometry,
class NumEqVector>
300 using GridView =
typename GridGeometry::GridView;
301 using Element =
typename GridView::template Codim<0>::Entity;
302 using Scalar =
typename NumEqVector::value_type;
306 static constexpr int dim = GridView::dimension;
307 static constexpr int dimWorld = GridView::dimensionworld;
309 using Tensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
310 static_assert(NumEqVector::dimension == dimWorld,
"Wrong dimension of velocity vector");
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
324 if (!problem.enableInertiaTerms())
328 const Scalar density = problem.density(fvGeometry.element(), fvGeometry, ipData(fvGeometry, scvf));
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();
336 const auto advectiveFlux = density*vn_integral * (upwindWeight * upwindVelocity + (1.0-upwindWeight)*downwindVelocity);
338 return advectiveFlux;
344 template<
class Context,
class IpData>
347 const auto& element = context.element();
348 const auto& fvGeometry = context.fvGeometry();
351 const auto mu = context.problem().effectiveViscosity(element, fvGeometry, ipData);
353 static const bool enableUnsymmetrizedVelocityGradient
354 =
getParamFromGroup<bool>(context.problem().paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
357 const auto& gradV = context.gradVelocity();
358 NumEqVector diffusiveFluxIntegrand = enableUnsymmetrizedVelocityGradient ?
359 mv(gradV, ipData.unitOuterNormal())
362 diffusiveFluxIntegrand *= -mu;
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();
368 return diffusiveFluxIntegrand;
371 template<
class Context,
class IpData>
374 const auto& element = context.element();
375 const auto& fvGeometry = context.fvGeometry();
378 const auto pressure = context.problem().pressure(element, fvGeometry, ipData);
383 const auto referencePressure = context.problem().referencePressure();
386 pn *= (pressure-referencePressure);
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.
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.