3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_ONEEQ_STAGGERED_FLUXVARIABLES_HH
26
27#include <numeric>
34
35namespace Dumux {
36
43// forward declaration
44template<class TypeTag, class BaseFluxVariables, class DiscretizationMethod>
45class OneEqFluxVariablesImpl;
46
47template<class TypeTag, class BaseFluxVariables>
48class OneEqFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethods::Staggered>
49: public BaseFluxVariables
50{
51 using ParentType = BaseFluxVariables;
52
54
55 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
56 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
57 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
58
59 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
60 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
61
62 using GridFaceVariables = typename GridVariables::GridFaceVariables;
63 using ElementFaceVariables = typename GridFaceVariables::LocalView;
64 using FaceVariables = typename GridFaceVariables::FaceVariables;
65
69 using FVElementGeometry = typename GridGeometry::LocalView;
70 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
71 using Extrusion = Extrusion_t<GridGeometry>;
72 using GridView = typename GridGeometry::GridView;
74 using Element = typename GridView::template Codim<0>::Entity;
76 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
77
78 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
79
80public:
81
85 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
86 const Element &element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const ElementFaceVariables& elemFaceVars,
90 const SubControlVolumeFace &scvf,
91 const FluxVariablesCache& fluxVarsCache)
92 {
93 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
94 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
95
96 // calculate advective flux
97 auto upwindTermK = [](const auto& volVars)
98 {
99 return volVars.viscosityTilde() * volVars.density();
100 };
101
102 flux[viscosityTildeEqIdx]
103 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
104
105 // calculate diffusive flux
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()];
110
111 // effective diffusion coefficients
112 Scalar insideCoeff = (insideVolVars.viscosity() + insideVolVars.viscosityTilde() * insideVolVars.density())
113 / insideVolVars.sigma();
114 Scalar outsideCoeff = (outsideVolVars.viscosity() + outsideVolVars.viscosityTilde() * outsideVolVars.density())
115 / outsideVolVars.sigma();
116
117 // scale by extrusion factor
118 insideCoeff *= insideVolVars.extrusionFactor();
119 outsideCoeff *= outsideVolVars.extrusionFactor();
120
121
122 Scalar coeff = 0.0;
123 Scalar distance = 0.0;
124 if (scvf.boundary())
125 {
126 coeff = insideCoeff;
127 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
128 }
129 else
130 {
131 // average and distance
132 coeff = arithmeticMean(insideCoeff, outsideCoeff,
133 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
134 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
135 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
136 }
137
138 const auto bcTypes = problem.boundaryTypes(element, scvf);
139 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::viscosityTildeEqIdx)
140 || bcTypes.isSymmetry())))
141 {
142 flux[viscosityTildeEqIdx]
143 += coeff / distance
144 * (insideVolVars.viscosityTilde() - outsideVolVars.viscosityTilde())
145 * Extrusion::area(scvf);
146 }
147 return flux;
148 }
149};
150
151} // end namespace
152
153#endif
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
Definition: adapt.hh:29
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.