version 3.8
freeflow/rans/oneeq/staggered/fluxvariables.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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
13#define DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
14
15#include <numeric>
22
23namespace Dumux {
24
31// forward declaration
32template<class TypeTag, class BaseFluxVariables, class DiscretizationMethod>
33class OneEqFluxVariablesImpl;
34
35template<class TypeTag, class BaseFluxVariables>
36class OneEqFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethods::Staggered>
37: public BaseFluxVariables
38{
39 using ParentType = BaseFluxVariables;
40
42
43 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
44 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
45 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
46
47 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
48 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
49
50 using GridFaceVariables = typename GridVariables::GridFaceVariables;
51 using ElementFaceVariables = typename GridFaceVariables::LocalView;
52 using FaceVariables = typename GridFaceVariables::FaceVariables;
53
57 using FVElementGeometry = typename GridGeometry::LocalView;
58 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
59 using Extrusion = Extrusion_t<GridGeometry>;
60 using GridView = typename GridGeometry::GridView;
62 using Element = typename GridView::template Codim<0>::Entity;
64 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
65
66 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
67
68public:
69
73 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
74 const Element &element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const ElementFaceVariables& elemFaceVars,
78 const SubControlVolumeFace &scvf,
79 const FluxVariablesCache& fluxVarsCache)
80 {
81 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
82 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
83
84 // calculate advective flux
85 auto upwindTermK = [](const auto& volVars)
86 {
87 return volVars.viscosityTilde() * volVars.density();
88 };
89
90 flux[viscosityTildeEqIdx]
91 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermK);
92
93 // calculate diffusive flux
94 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
95 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
96 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
97 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
98
99 // effective diffusion coefficients
100 Scalar insideCoeff = (insideVolVars.viscosity() + insideVolVars.viscosityTilde() * insideVolVars.density())
101 / insideVolVars.sigma();
102 Scalar outsideCoeff = (outsideVolVars.viscosity() + outsideVolVars.viscosityTilde() * outsideVolVars.density())
103 / outsideVolVars.sigma();
104
105 // scale by extrusion factor
106 insideCoeff *= insideVolVars.extrusionFactor();
107 outsideCoeff *= outsideVolVars.extrusionFactor();
108
109
110 Scalar coeff = 0.0;
111 Scalar distance = 0.0;
112 if (scvf.boundary())
113 {
114 coeff = insideCoeff;
115 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
116 }
117 else
118 {
119 // average and distance
120 coeff = arithmeticMean(insideCoeff, outsideCoeff,
121 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
122 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
123 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
124 }
125
126 const auto bcTypes = problem.boundaryTypes(element, scvf);
127 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::viscosityTildeEqIdx)
128 || bcTypes.isSymmetry())))
129 {
130 flux[viscosityTildeEqIdx]
131 += coeff / distance
132 * (insideVolVars.viscosityTilde() - outsideVolVars.viscosityTilde())
133 * Extrusion::area(fvGeometry, scvf);
134 }
135 return flux;
136 }
137};
138
139} // end namespace
140
141#endif
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:73
The flux variables class for the one-equation model by Spalart-Allmaras using the staggered grid disc...
Definition: freeflow/rans/oneeq/fluxvariables.hh:23
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Base class for the flux variables living on a sub control volume face.
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:38
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:282
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166