3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/oneeq/staggered/localresidual.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_STAGGERED_ONEEQ_LOCAL_RESIDUAL_HH
25#define DUMUX_STAGGERED_ONEEQ_LOCAL_RESIDUAL_HH
26
27#include <dune/common/hybridutilities.hh>
31
32namespace Dumux {
33
40// forward declaration
41template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod>
42class OneEqResidualImpl;
43
44template<class TypeTag, class BaseLocalResidual>
45class OneEqResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
46: public BaseLocalResidual
47{
48 using ParentType = BaseLocalResidual;
49
51
52 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
53 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
54 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
55
56 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
57 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
59
60 using GridFaceVariables = typename GridVariables::GridFaceVariables;
61 using ElementFaceVariables = typename GridFaceVariables::LocalView;
62
66 using Element = typename GridView::template Codim<0>::Entity;
67 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
68 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
69 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
72
73 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
74
75public:
76 using ParentType::ParentType;
77
79 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars) const
82 {
83 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
84 storage[viscosityTildeEqIdx] = volVars.viscosityTilde() * volVars.density();
85 return storage;
86 }
87
88 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
89 const Element &element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const ElementFaceVariables& elemFaceVars,
93 const SubControlVolume &scv) const
94 {
95 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
96 elemVolVars, elemFaceVars, scv);
97
98 const auto& volVars = elemVolVars[scv];
99
100 source[viscosityTildeEqIdx] += volVars.cb1() * (1.0 - volVars.ft2())
101 * volVars.stressTensorScalarProductTilde()
102 * volVars.viscosityTilde() * volVars.density();
103
104 source[viscosityTildeEqIdx] -= (volVars.cw1() * volVars.fW()
105 - volVars.cb1() * volVars.ft2() / problem.karmanConstant() / problem.karmanConstant())
106 * volVars.viscosityTilde() * volVars.viscosityTilde()
107 / volVars.wallDistance() / volVars.wallDistance() * volVars.density();;
108
109 for (unsigned int axisIdx = 0; axisIdx < ModelTraits::dim(); ++axisIdx)
110 {
111 source[viscosityTildeEqIdx] += volVars.cb2() / volVars.sigma()
112 * volVars.storedViscosityTildeGradient()[axisIdx]
113 * volVars.storedViscosityTildeGradient()[axisIdx]
114 * volVars.density();
115 }
116
117 return source;
118 }
119};
120} // end namespace Dumux
121
122#endif
The available discretization methods in Dumux.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Element-wise calculation of the residual for one-equation turbulence models using the staggered discr...
Definition: freeflow/rans/oneeq/localresidual.hh:36
CellCenterPrimaryVariables computeStorageForCellCenter(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate fluxes entering or leaving the cell center control volume.
Definition: freeflow/rans/oneeq/staggered/localresidual.hh:79
CellCenterPrimaryVariables computeSourceForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolume &scv) const
Definition: freeflow/rans/oneeq/staggered/localresidual.hh:88
Declares all properties used in Dumux.