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-2020 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 
36 namespace Foam
37 {
38 namespace incompressible
39 {
40 namespace RASVariables
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
54  if (solverControl_.average())
55  {
56  GMean_.reset
57  (
59  (
60  IOobject
61  (
62  "GMean",
63  mesh_.time().timeName(),
64  mesh_,
67  ),
68  mesh_,
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_ = "omega";
88 
92 
93  distPtr_.ref(const_cast<volScalarField&>(wallDist::New(mesh_).y()));
94 
97 }
98 
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();
113  tmp<volTensorField> tgradU = fvc::grad(U);
115  (
116  this->type() + ":GbyNu",
117  (tgradU() && dev(twoSymm(tgradU())))
118  );
119  auto tG =
121  (
122  turbModel.GName(),
123  nutRefInst()*GbyNu0
124  );
125  // Use correctBoundaryConditions instead of updateCoeffs to avoid
126  // messing with updateCoeffs in the next iteration of omegaEqn
129  return tG;
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 {
138  {
139  DebugInfo
140  << "Using GMean" << endl;
142  }
143  DebugInfo
144  << "Using instantaneous G" << endl;
145  return computeG();
146 }
147 
148 
150 {
153  {
154  const label iAverageIter = solverControl_.averageIter();
155  scalar avIter(iAverageIter);
156  scalar oneOverItP1 = 1./(avIter + 1);
157  scalar mult = avIter*oneOverItP1;
158  GMean_() = GMean_()*mult + computeG()*oneOverItP1;
159  }
160 }
161 
162 
164 (
166 )
167 {
168  // The presence of G is required to update the boundary value of omega
169  const volVectorField& U = turbulence.U();
170  const volScalarField S2(2*magSqr(symm(fvc::grad(U))));
171 
172  volScalarField G(turbulence.GName(), nutRef() * S2);
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace RASVariables
180 } // End namespace incompressible
181 } // End namespace Foam
182 
183 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
virtual tmp< volScalarField::Internal > G()
Return the turbulence production term.
Definition: kOmegaSST.C:128
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...
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 volVectorField & U() const
Access function to velocity field.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
virtual void computeMeanFields()
Compute mean fields on the fly.
Definition: kOmegaSST.C:142
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
kOmegaSST(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
Definition: kOmegaSST.C:72
const volScalarField & TMVar2() const
const volScalarField & nutRef() const
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
word GName() const
Helper function to return the name of the turbulence G field.
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
constexpr const char *const group
Group name for atomic constants.
Templated abstract base class for single-phase incompressible turbulence models.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const word propertiesName
Default name of the turbulence properties dictionary.
Base class for solver control classes.
Definition: solverControl.H:45
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
Correct boundary conditions of turbulent fields.
Definition: kOmegaSST.C:157
#define DebugInfo
Report an information message using Foam::Info.
autoPtr< volScalarField::Internal > GMean_
Average of the production term.
Definition: kOmegaSST.H:66
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
bool doAverageIter() const
Whether or not to add fields of the current iteration to the average fields.
const volScalarField & TMVar2Inst() const
dimensionedScalar pow3(const dimensionedScalar &ds)
const volScalarField & nutRefInst() const
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(RASModelVariables, kEpsilon, dictionary)
bool average() const
Whether averaging is enabled or not.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
label & averageIter()
Return average iteration index reference.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
tmp< volScalarField::Internal > computeG()
Definition: kOmegaSST.C:93
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
correct bounday conditions of turbulent fields
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
virtual void computeMeanFields()
Compute mean fields on the fly.