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 SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
102 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
104 static constexpr int dim = FVElementGeometry::GridGeometry::GridView::dimension;
105 static constexpr int dimWorld = FVElementGeometry::GridGeometry::GridView::dimensionworld;
107 using Tensor = Dune::FieldMatrix<typename GlobalPosition::value_type, dim, dimWorld>;
128 {
return fvGeometry_.element(); }
131 {
return fvGeometry_; }
134 {
return elemVolVars_; }
141 GlobalPosition v(0.0);
142 const auto& shapeValues = ipData_.shapeValues();
143 for (
const auto& localDof :
localDofs(fvGeometry_))
144 v.axpy(shapeValues[localDof.index()][0], elemVolVars_[localDof.index()].velocity());
152 for (
const auto& localDof :
localDofs(fvGeometry_))
154 const auto& volVars = elemVolVars_[localDof.index()];
155 for (
int dir = 0; dir < dim; ++dir)
156 gradV[dir].axpy(volVars.velocity(dir), ipData_.gradN(localDof.index()));
162 const Problem& problem_;
163 const FVElementGeometry& fvGeometry_;
164 const ElementVolumeVariables& elemVolVars_;
165 const IpData& ipData_;
172template<
class Gr
idGeometry,
class NumEqVector>
175 using GridView =
typename GridGeometry::GridView;
176 using Element =
typename GridView::template Codim<0>::Entity;
177 using Scalar =
typename NumEqVector::value_type;
181 static constexpr int dim = GridView::dimension;
182 static constexpr int dimWorld = GridView::dimensionworld;
184 using Tensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
185 static_assert(NumEqVector::dimension == dimWorld,
"Wrong dimension of velocity vector");
191 template<
class Context>
194 if (!context.problem().enableInertiaTerms())
197 const auto& fvGeometry = context.fvGeometry();
198 const auto& elemVolVars = context.elemVolVars();
199 const auto& scvf = context.scvFace();
200 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
201 const auto& shapeValues = fluxVarCache.shapeValues();
205 for (
const auto& localDof :
localDofs(fvGeometry))
206 v.axpy(shapeValues[localDof.index()][0], elemVolVars[localDof.index()].velocity());
209 const Scalar
density = context.problem().density(context.element(), context.fvGeometry(), fluxVarCache.ipData());
211 const auto vn = v*scvf.unitOuterNormal();
212 const auto& insideVolVars = elemVolVars[fvGeometry.scv(scvf.insideScvIdx())];
213 const auto& outsideVolVars = elemVolVars[fvGeometry.scv(scvf.outsideScvIdx())];
214 const auto upwindVelocity = vn > 0 ? insideVolVars.velocity() : outsideVolVars.velocity();
215 const auto downwindVelocity = vn > 0 ? outsideVolVars.velocity() : insideVolVars.velocity();
216 static const auto upwindWeight = getParamFromGroup<Scalar>(context.problem().paramGroup(),
"Flux.UpwindWeight");
217 const auto advectiveTermIntegrand =
density*vn * (upwindWeight * upwindVelocity + (1.0-upwindWeight)*downwindVelocity);
219 return advectiveTermIntegrand * Extrusion::area(fvGeometry, scvf) * insideVolVars.extrusionFactor();
225 template<
class Context>
228 const auto& element = context.element();
229 const auto& fvGeometry = context.fvGeometry();
230 const auto& elemVolVars = context.elemVolVars();
231 const auto& scvf = context.scvFace();
232 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
236 for (
const auto& localDof :
localDofs(fvGeometry))
238 const auto& volVars = elemVolVars[localDof];
239 for (
int dir = 0; dir < dim; ++dir)
240 gradV[dir].axpy(volVars.velocity(dir), fluxVarCache.gradN(localDof.index()));
244 const auto mu = context.problem().effectiveViscosity(element, fvGeometry, fluxVarCache.ipData());
246 static const bool enableUnsymmetrizedVelocityGradient
247 = getParamFromGroup<bool>(context.problem().paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
250 NumEqVector diffusiveFlux = enableUnsymmetrizedVelocityGradient ?
251 mv(gradV, scvf.unitOuterNormal())
254 diffusiveFlux *= -mu;
256 static const bool enableDilatationTerm = getParamFromGroup<bool>(context.problem().paramGroup(),
"FreeFlow.EnableDilatationTerm",
false);
257 if (enableDilatationTerm)
258 diffusiveFlux += 2.0/3.0 * mu *
trace(gradV) * scvf.unitOuterNormal();
260 diffusiveFlux *= Extrusion::area(fvGeometry, scvf) * elemVolVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
261 return diffusiveFlux;
264 template<
class Context>
267 const auto& element = context.element();
268 const auto& fvGeometry = context.fvGeometry();
269 const auto& elemVolVars = context.elemVolVars();
270 const auto& scvf = context.scvFace();
271 const auto& fluxVarCache = context.elemFluxVarsCache()[scvf];
274 const auto pressure = context.problem().pressure(element, fvGeometry, fluxVarCache.ipData());
280 const auto referencePressure = context.problem().referencePressure();
283 pn *= (
pressure-referencePressure)*Extrusion::area(fvGeometry, scvf)*elemVolVars[fvGeometry.scv(scvf.insideScvIdx())].extrusionFactor();
The flux variables class for the Navier-Stokes model using control-volume finite element schemes.
Definition: flux.hh:174
NumEqVector advectiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:192
NumEqVector diffusiveMomentumFlux(const Context &context) const
Returns the diffusive momentum flux due to viscous forces.
Definition: flux.hh:226
NumEqVector pressureContribution(const Context &context) const
Definition: flux.hh:265
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
Context for interpolating data on interpolation points.
Definition: flux.hh:99
const FVElementGeometry & fvGeometry() const
Definition: flux.hh:130
Tensor gradVelocity() const
Definition: flux.hh:149
const Element & element() const
Definition: flux.hh:127
const Problem & problem() const
Definition: flux.hh:124
NavierStokesMomentumFluxFunctionContext(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const IpData &ipData)
Initialize the flux variables storing some temporary pointers.
Definition: flux.hh:112
const IpData & ipData() const
Definition: flux.hh:136
GlobalPosition velocity() const
Definition: flux.hh:139
const ElementVolumeVariables & elemVolVars() const
Definition: flux.hh:133
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: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
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.