24#ifndef DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
44template<
class TypeTag,
class BaseFluxVariables,
class DiscretizationMethod>
45class OneEqFluxVariablesImpl;
47template<
class TypeTag,
class BaseFluxVariables>
49:
public BaseFluxVariables
51 using ParentType = BaseFluxVariables;
55 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
56 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
57 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
59 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
60 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
62 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
63 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
64 using FaceVariables =
typename GridFaceVariables::FaceVariables;
69 using FVElementGeometry =
typename GridGeometry::LocalView;
70 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
72 using GridView =
typename GridGeometry::GridView;
74 using Element =
typename GridView::template Codim<0>::Entity;
78 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
86 const Element &element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const ElementFaceVariables& elemFaceVars,
90 const SubControlVolumeFace &scvf,
91 const FluxVariablesCache& fluxVarsCache)
93 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
94 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
97 auto upwindTermK = [](
const auto& volVars)
99 return volVars.viscosityTilde() * volVars.density();
102 flux[viscosityTildeEqIdx]
103 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
106 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
107 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
108 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
109 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
112 Scalar insideCoeff = (insideVolVars.viscosity() + insideVolVars.viscosityTilde() * insideVolVars.density())
113 / insideVolVars.sigma();
114 Scalar outsideCoeff = (outsideVolVars.viscosity() + outsideVolVars.viscosityTilde() * outsideVolVars.density())
115 / outsideVolVars.sigma();
118 insideCoeff *= insideVolVars.extrusionFactor();
119 outsideCoeff *= outsideVolVars.extrusionFactor();
127 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
133 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
134 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
135 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
138 const auto bcTypes = problem.boundaryTypes(element, scvf);
139 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::viscosityTildeEqIdx)
140 || bcTypes.isSymmetry())))
142 flux[viscosityTildeEqIdx]
144 * (insideVolVars.viscosityTilde() - outsideVolVars.viscosityTilde())
145 * Extrusion::area(scvf);
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:292
constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) arithmetic mean of two scalar values.
Definition: math.hh:50
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
The flux variables class for the one-equation model by Spalart-Allmaras using the staggered grid disc...
Definition: freeflow/rans/oneeq/fluxvariables.hh:35
CellCenterPrimaryVariables computeMassFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Computes the flux for the cell center residual.
Definition: freeflow/rans/oneeq/staggered/fluxvariables.hh:85
Declares all properties used in Dumux.