version 3.8
freeflow/rans/twoeq/kepsilon/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// 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_STAGGERED_KEPSILON_LOCAL_RESIDUAL_HH
13#define DUMUX_STAGGERED_KEPSILON_LOCAL_RESIDUAL_HH
14
15#include <dune/common/hybridutilities.hh>
19
20namespace Dumux {
21
27// forward declaration
28template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod>
30
31template<class TypeTag, class BaseLocalResidual>
32class KEpsilonResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
33: public BaseLocalResidual
34{
35 using ParentType = BaseLocalResidual;
36 friend class StaggeredLocalResidual<TypeTag>;
37
39
40 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
41 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
42 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
43
44 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
45 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
47
48 using GridFaceVariables = typename GridVariables::GridFaceVariables;
49 using ElementFaceVariables = typename GridFaceVariables::LocalView;
50
54 using Element = typename GridView::template Codim<0>::Entity;
55 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
56 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
57 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
60
61 using CellCenterResidual = CellCenterPrimaryVariables;
62
63 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
64 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
65
66public:
67 using ParentType::ParentType;
68
69 // account for the offset of the cell center privars within the PrimaryVariables container
70 static constexpr auto cellCenterOffset = ModelTraits::numEq() - CellCenterPrimaryVariables::dimension;
71 static_assert(cellCenterOffset == ModelTraits::dim(), "cellCenterOffset must equal dim for staggered NavierStokes");
72
74 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
75 const SubControlVolume& scv,
76 const VolumeVariables& volVars) const
77 {
78 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
79
80 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy() * volVars.density();
81 storage[dissipationEqIdx] = volVars.dissipation() * volVars.density();
82
83 return storage;
84 }
85
86 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
87 const Element &element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const ElementFaceVariables& elemFaceVars,
91 const SubControlVolume &scv) const
92 {
93 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
94 elemVolVars, elemFaceVars, scv);
95
96 const auto& volVars = elemVolVars[scv];
97 Scalar turbulentKineticEnergy = volVars.turbulentKineticEnergy();
98 Scalar dissipation = volVars.dissipation();
99
100 // production
101 // turbulence production is equal to dissipation -> exclude both terms (according to local equilibrium hypothesis, see FLUENT)
102 if (!volVars.isMatchingPoint())
103 {
104 source[turbulentKineticEnergyEqIdx] += 2.0 * volVars.dynamicEddyViscosity()
105 * volVars.stressTensorScalarProduct();
106 }
107 source[dissipationEqIdx] += volVars.cOneEpsilon()
108 * dissipation / turbulentKineticEnergy
109 * 2.0 * volVars.dynamicEddyViscosity()
110 * volVars.stressTensorScalarProduct();
111
112 // destruction
113 // turbulence production is equal to dissipation -> exclude both terms (according to local equilibrium hypothesis, see FLUENT)
114 if (!volVars.isMatchingPoint())
115 {
116 source[turbulentKineticEnergyEqIdx] -= dissipation * volVars.density();
117 }
118 source[dissipationEqIdx] -= volVars.cTwoEpsilon()
119 * dissipation * dissipation
120 / turbulentKineticEnergy * volVars.density();
121
122 return source;
123 }
124};
125} // end namespace Dumux
126
127#endif
CellCenterPrimaryVariables computeSourceForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolume &scv) const
Definition: freeflow/rans/twoeq/kepsilon/staggered/localresidual.hh:86
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/twoeq/kepsilon/staggered/localresidual.hh:74
Element-wise calculation of the residual for k-epsilon models using the staggered discretization.
Definition: freeflow/rans/twoeq/kepsilon/localresidual.hh:24
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:29
Defines all properties used in Dumux.
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