3.1-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>
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, class BaseFluxVariables, DiscretizationMethod discMethod>
38class OneEqFluxVariablesImpl;
39
45template<class TypeTag, class BaseFluxVariables>
46class OneEqFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethod::staggered>
47: public BaseFluxVariables
48{
49 using ParentType = BaseFluxVariables;
50
52
53 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
54 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
55 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
56
57 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
58 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
59
60 using GridFaceVariables = typename GridVariables::GridFaceVariables;
61 using ElementFaceVariables = typename GridFaceVariables::LocalView;
62 using FaceVariables = typename GridFaceVariables::FaceVariables;
63
68 using Element = typename GridView::template Codim<0>::Entity;
69 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
71 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
72 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
73
74 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
75
76public:
77
81 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
82 const Element &element,
83 const FVElementGeometry& fvGeometry,
84 const ElementVolumeVariables& elemVolVars,
85 const ElementFaceVariables& elemFaceVars,
86 const SubControlVolumeFace &scvf,
87 const FluxVariablesCache& fluxVarsCache)
88 {
89 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
90 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
91
92 // calculate advective flux
93 auto upwindTermK = [](const auto& volVars)
94 {
95 return volVars.viscosityTilde();
96 };
97
98 flux[viscosityTildeEqIdx]
99 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
100
101 // calculate diffusive flux
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()];
106
107 // effective diffusion coefficients
108 Scalar insideCoeff = (insideVolVars.kinematicViscosity() + insideVolVars.viscosityTilde())
109 / insideVolVars.sigma();
110 Scalar outsideCoeff = (outsideVolVars.kinematicViscosity() + outsideVolVars.viscosityTilde())
111 / outsideVolVars.sigma();
112
113 // scale by extrusion factor
114 insideCoeff *= insideVolVars.extrusionFactor();
115 outsideCoeff *= outsideVolVars.extrusionFactor();
116
117 // average and distance
118 Scalar coeff = arithmeticMean(insideCoeff, outsideCoeff,
119 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
120 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
121 Scalar distance = 0.0;
122
123 // adapt boundary handling
124 if (scvf.boundary())
125 {
126 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
127 coeff = insideCoeff;
128 }
129 else
130 {
131 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
132 }
133
134 const auto bcTypes = problem.boundaryTypes(element, scvf);
135 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::viscosityTildeEqIdx)
136 || bcTypes.isSymmetry())))
137 {
138 flux[viscosityTildeEqIdx]
139 += coeff / distance
140 * (insideVolVars.viscosityTilde() - outsideVolVars.viscosityTilde())
141 * scvf.area();
142 }
143 return flux;
144 }
145};
146
147} // end namespace
148
149#endif
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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.