Loading...
Searching...
No Matches
kOmegaSST.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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019-2023 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 "kOmegaSST.H"
31#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
40namespace RASVariables
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47
48
49// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50
52{
54 if (solverControl_.average())
55 {
56 GMean_.reset
57 (
59 (
61 (
62 "GMean",
63 mesh_.time().timeName(),
64 mesh_,
68 ),
69 mesh_,
71 )
72 );
73 }
74}
75
76
77// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78
80(
81 const fvMesh& mesh,
82 const solverControl& SolverControl
83)
84:
85 RASModelVariables(mesh, SolverControl)
86{
87 TMVar1BaseName_ = "k";
88 TMVar2BaseName_ = "omega";
89
93
94 distPtr_.ref(wallDist::New(mesh_).y().constCast());
95
98}
99
100
101// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102
104{
106 (
108 (
110 TMVar2().internalField().group()
111 )
112 );
113 // Recompute G and modify values next to the walls
114 // Ideally, grad(U) should be cached to avoid the overhead
115 const volVectorField& U = turbModel.U();
118 (
119 IOobject::scopedName(this->type(), "GbyNu"),
120 (tgradU() && devTwoSymm(tgradU()))
121 );
122
123 // NB: leave tmp registered (for correctBoundaryConditions)
125 (
126 turbModel.GName(),
128 nutRefInst()*GbyNu0
129 );
130
131 // Use correctBoundaryConditions instead of updateCoeffs to avoid
132 // messing with updateCoeffs in the next iteration of omegaEqn
135 return tG;
136}
137
138
139// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140
142{
144 {
146 << "Using GMean" << endl;
150 << "Using instantaneous G" << endl;
151 return computeG();
152}
153
154
156{
159 {
160 const label iAverageIter = solverControl_.averageIter();
161 scalar avIter(iAverageIter);
162 scalar oneOverItP1 = 1./(avIter + 1);
163 scalar mult = avIter*oneOverItP1;
164 GMean_() = GMean_()*mult + computeG()*oneOverItP1;
165 }
166}
167
168
170(
172)
173{
174 // The presence of G is required to update the boundary value of omega
175 const volVectorField& U = turbulence.U();
176 const volScalarField S2(2*magSqr(symm(fvc::grad(U))));
177
178 volScalarField G(turbulence.GName(), nutRef() * S2);
180}
181
182
183// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184
185} // End namespace RASVariables
186} // End namespace incompressible
187} // End namespace Foam
188
189// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
void correctBoundaryConditions()
Correct boundary field.
@ REGISTER
Request registration (bool: true).
@ 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.
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
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
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
correct bounday conditions of turbulent fields
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 & nutRef() const
const volScalarField & TMVar2Inst() const
tmp< volScalarField::Internal > computeG()
Definition kOmegaSST.C:96
virtual tmp< volScalarField::Internal > G()
Return the turbulence production term.
Definition kOmegaSST.C:134
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
Correct boundary conditions of turbulent fields.
Definition kOmegaSST.C:163
kOmegaSST(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
Definition kOmegaSST.C:73
virtual void computeMeanFields()
Compute mean fields on the fly.
Definition kOmegaSST.C:148
autoPtr< volScalarField::Internal > GMean_
Average of the production term.
Definition kOmegaSST.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 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.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
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.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)