Loading...
Searching...
No Matches
kEpsilon.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2020 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "kEpsilon.H"
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
39namespace RASVariables
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51{
53 if (solverControl_.average())
54 {
55 GMean_.reset
56 (
58 (
60 (
61 "GMean",
62 mesh_.time().timeName(),
63 mesh_,
66 ),
67 mesh_,
69 )
70 );
71 }
72}
73
74
75
76// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77
79(
80 const fvMesh& mesh,
81 const solverControl& SolverControl
82)
83:
84 RASModelVariables(mesh, SolverControl)
85{
86 TMVar1BaseName_ = "k";
87 TMVar2BaseName_ = "epsilon";
88
92
95}
96
97
98// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99
101{
103 (
105 (
107 TMVar2().internalField().group()
108 )
109 );
110 // Recompute G and modify values next to the walls
111 // Ideally, grad(U) should be cached to avoid the overhead
112 const volVectorField& U = turbModel.U();
115 (
116 IOobject::scopedName(this->type(), "GbyNu"),
117 (tgradU() && devTwoSymm(tgradU()))
118 );
119
120 // NB: leave tmp registered (for correctBoundaryConditions)
121 auto tG =
123 (
124 turbModel.GName(),
125 nutRefInst()*GbyNu0
126 );
127 // Use correctBoundaryConditions instead of updateCoeffs to avoid
128 // messing with updateCoeffs in the next iteration of omegaEqn
131 return tG;
132}
133
134
135// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136
138{
140 {
142 << "Using GMean" << endl;
146 << "Using instantaneous G" << endl;
147 return computeG();
148}
149
150
152{
155 {
156 const label iAverageIter = solverControl_.averageIter();
157 scalar avIter(iAverageIter);
158 scalar oneOverItP1 = 1./(avIter + 1);
159 scalar mult = avIter*oneOverItP1;
160 GMean_() = GMean_()*mult + computeG()*oneOverItP1;
161 }
162}
163
164
165
166// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167
168} // End namespace RASVariables
169} // End namespace incompressible
170} // End namespace Foam
171
172// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
DimensionedField< scalar, volMesh > Internal
void correctBoundaryConditions()
Correct boundary field.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
const volScalarField & TMVar2() const
RASModelVariables(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
const volScalarField & nutRefInst() const
virtual void computeMeanFields()
Compute mean fields on the fly.
const volScalarField & TMVar2Inst() const
tmp< volScalarField::Internal > computeG()
Definition kEpsilon.C:93
virtual tmp< volScalarField::Internal > G()
Return the turbulence production term.
Definition kEpsilon.C:130
kEpsilon(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
Definition kEpsilon.C:72
virtual void computeMeanFields()
Compute mean fields on the fly.
Definition kEpsilon.C:144
autoPtr< volScalarField::Internal > GMean_
Average of the production term.
Definition kEpsilon.H:66
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Base class for solver control classes.
bool doAverageIter() const
Whether or not to add fields of the current iteration to the average fields.
label & averageIter()
Return average iteration index reference.
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
static const word propertiesName
Default name of the turbulence properties dictionary.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
const dimensionSet dimArea(sqr(dimLength))
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.