24#ifndef DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
37template<
class TypeTag,
class BaseFluxVariables, DiscretizationMethod discMethod>
38class OneEqFluxVariablesImpl;
45template<
class TypeTag,
class BaseFluxVariables>
47:
public BaseFluxVariables
49 using ParentType = BaseFluxVariables;
53 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
54 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
55 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
57 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
58 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
60 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
61 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
62 using FaceVariables =
typename GridFaceVariables::FaceVariables;
68 using Element =
typename GridView::template Codim<0>::Entity;
71 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
74 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
82 const Element &element,
83 const FVElementGeometry& fvGeometry,
84 const ElementVolumeVariables& elemVolVars,
85 const ElementFaceVariables& elemFaceVars,
86 const SubControlVolumeFace &scvf,
87 const FluxVariablesCache& fluxVarsCache)
89 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
90 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
93 auto upwindTermK = [](
const auto& volVars)
95 return volVars.viscosityTilde();
98 flux[viscosityTildeEqIdx]
99 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
102 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
103 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
104 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
105 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
108 Scalar insideCoeff = (insideVolVars.kinematicViscosity() + insideVolVars.viscosityTilde())
109 / insideVolVars.sigma();
110 Scalar outsideCoeff = (outsideVolVars.kinematicViscosity() + outsideVolVars.viscosityTilde())
111 / outsideVolVars.sigma();
114 insideCoeff *= insideVolVars.extrusionFactor();
115 outsideCoeff *= outsideVolVars.extrusionFactor();
119 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
120 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
121 Scalar distance = 0.0;
126 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
131 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
134 const auto bcTypes = problem.boundaryTypes(element, scvf);
135 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::viscosityTildeEqIdx)
136 || bcTypes.isSymmetry())))
138 flux[viscosityTildeEqIdx]
140 * (insideVolVars.viscosityTilde() - outsideVolVars.viscosityTilde())
Base class for the flux variables living on a sub control volume face.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:49
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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:81
Declares all properties used in Dumux.