All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
kOmegaSSTDES.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "kOmegaSSTDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicMomentumTransportModel>
38 tmp<volScalarField::Internal>
40 {
41  return sqrt(this->k_())/(this->betaStar_*this->omega_());
42 }
43 
44 
45 template<class BasicMomentumTransportModel>
47 (
48  const volScalarField::Internal& F1,
49  const volScalarField::Internal& F2
50 ) const
51 {
52  switch (FSST_)
53  {
54  case 0:
55  return max(Lt()/(CDES_*this->delta()()), scalar(1));
56  case 1:
57  return max(Lt()*(1 - F1)/(CDES_*this->delta()()), scalar(1));
58  case 2:
59  return max(Lt()*(1 - F2)/(CDES_*this->delta()()), scalar(1));
60  default:
62  << "Incorrect FSST = " << FSST_ << ", should be 0, 1 or 2"
63  << exit(FatalError);
64  return F1;
65  }
66 }
67 
68 
69 template<class BasicMomentumTransportModel>
72 (
73  const volScalarField::Internal& F1,
74  const volScalarField::Internal& F2
75 ) const
76 {
77  return this->betaStar_*this->omega_()*FDES(F1, F2);
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
83 template<class BasicMomentumTransportModel>
85 (
86  const alphaField& alpha,
87  const rhoField& rho,
88  const volVectorField& U,
89  const surfaceScalarField& alphaRhoPhi,
90  const surfaceScalarField& phi,
91  const viscosity& viscosity,
92  const word& type
93 )
94 :
95  kOmegaSST
96  <
97  LESeddyViscosity<BasicMomentumTransportModel>,
98  BasicMomentumTransportModel
99  >
100  (
101  type,
102  alpha,
103  rho,
104  U,
105  alphaRhoPhi,
106  phi,
107  viscosity
108  ),
109 
110  CDES_
111  (
113  (
114  "CDES",
115  this->coeffDict_,
116  0.61
117  )
118  ),
119  FSST_(this->coeffDict_.lookupOrDefault("FSST", 2))
120 {
121  if (type == typeName)
122  {
123  this->printCoeffs(type);
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
130 template<class BasicMomentumTransportModel>
132 {
133  if
134  (
135  kOmegaSST
136  <
137  LESeddyViscosity<BasicMomentumTransportModel>,
138  BasicMomentumTransportModel
139  >::read()
140  )
141  {
142  CDES_.readIfPresent(this->coeffDict());
143  this->coeffDict().readIfPresent("FSST", FSST_);
144 
145  return true;
146  }
147  else
148  {
149  return false;
150  }
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace LESModels
157 } // End namespace Foam
158 
159 // ************************************************************************* //
scalar delta
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
Definition: kOmegaSSTDES.C:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
U
Definition: pEqn.H:72
error FatalError
#define F1(B, C, D)
Definition: SHA1.C:173
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmegaSSTDES.H:120
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
phi
Definition: correctPhi.H:3
kOmegaSSTDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTDES.C:85
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmegaSSTDES.H:121
virtual tmp< volScalarField::Internal > FDES(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
The DES dissipation-rate multiplier with options zonal filtering.
Definition: kOmegaSSTDES.C:47
tmp< volScalarField::Internal > Lt() const
Return the turbulent length-scale.
Definition: kOmegaSSTDES.C:39
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:131