version 3.8
freeflow/rans/twoeq/lowrekepsilon/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_LOWREKEPSILON_LOCAL_RESIDUAL_HH
13#define DUMUX_STAGGERED_LOWREKEPSILON_LOCAL_RESIDUAL_HH
14
15#include <dune/common/hybridutilities.hh>
19
20namespace Dumux {
21
27// forward declaration
28template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod>
29class LowReKEpsilonResidualImpl;
30
31template<class TypeTag, class BaseLocalResidual>
32class LowReKEpsilonResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
33: public BaseLocalResidual
34{
35 using ParentType = BaseLocalResidual;
36
38
39 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
40 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
41 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
42
43 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
44 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
46
47 using GridFaceVariables = typename GridVariables::GridFaceVariables;
48 using ElementFaceVariables = typename GridFaceVariables::LocalView;
49
53 using Element = typename GridView::template Codim<0>::Entity;
54 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
55 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
56 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
59
60 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
61 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
62
63public:
64 using ParentType::ParentType;
65
67 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
68 const SubControlVolume& scv,
69 const VolumeVariables& volVars) const
70 {
71 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
72
73 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy() * volVars.density();
74 storage[dissipationEqIdx] = volVars.dissipationTilde() * volVars.density();
75
76 return storage;
77 }
78
79 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
80 const Element &element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const ElementFaceVariables& elemFaceVars,
84 const SubControlVolume &scv) const
85 {
86 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
87 elemVolVars, elemFaceVars, scv);
88
89 const auto& volVars = elemVolVars[scv];
90
91 // production
92 source[turbulentKineticEnergyEqIdx] += 2.0 * volVars.dynamicEddyViscosity()
93 * volVars.stressTensorScalarProduct();
94 source[dissipationEqIdx] += volVars.cOneEpsilon() * volVars.fOne()
95 * volVars.dissipationTilde() / volVars.turbulentKineticEnergy()
96 * 2.0 * volVars.dynamicEddyViscosity()
97 * volVars.stressTensorScalarProduct();
98
99 // destruction
100 source[turbulentKineticEnergyEqIdx] -= volVars.dissipationTilde() * volVars.density();
101 source[dissipationEqIdx] -= volVars.cTwoEpsilon() * volVars.fTwo() * volVars.density()
102 * volVars.dissipationTilde() * volVars.dissipationTilde()
103 / volVars.turbulentKineticEnergy();
104
105 // dampening functions
106 source[turbulentKineticEnergyEqIdx] -= volVars.dValue() * volVars.density();
107 source[dissipationEqIdx] += volVars.eValue() * volVars.density();
108
109 return source;
110 }
111};
112} // end namespace Dumux
113
114#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/lowrekepsilon/staggered/localresidual.hh:79
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/lowrekepsilon/staggered/localresidual.hh:67
Element-wise calculation of the residual for low-Reynolds k-epsilon models using the staggered discre...
Definition: freeflow/rans/twoeq/lowrekepsilon/localresidual.hh:24
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