buoyantKEpsilon.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 | Copyright (C) 2014-2015 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 "buoyantKEpsilon.H"
28 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasicTurbulenceModel>
42 (
43  const alphaField& alpha,
44  const rhoField& rho,
45  const volVectorField& U,
46  const surfaceScalarField& alphaRhoPhi,
47  const surfaceScalarField& phi,
48  const transportModel& transport,
49  const word& propertiesName,
50  const word& type
51 )
52 :
54  (
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  transport,
61  propertiesName,
62  type
63  ),
64 
65  Cg_
66  (
68  (
69  "Cg",
70  this->coeffDict_,
71  1.0
72  )
73  )
74 {
75  if (type == typeName)
76  {
78  this->printCoeffs(type);
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class BasicTurbulenceModel>
87 {
89  {
90  Cg_.readIfPresent(this->coeffDict());
91 
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 
101 template<class BasicTurbulenceModel>
104 {
106  this->mesh_.objectRegistry::template
107  lookupObject<uniformDimensionedVectorField>("g");
108 
109  return
110  (Cg_*this->Cmu_)*this->alpha_*this->k_*(g & fvc::grad(this->rho_))
111  /(this->epsilon_ + this->epsilonMin_);
112 }
113 
114 
115 template<class BasicTurbulenceModel>
118 {
120  this->mesh_.objectRegistry::template
121  lookupObject<uniformDimensionedVectorField>("g");
122 
123  if (mag(g.value()) > SMALL)
124  {
125  return -fvm::SuSp(Gcoef(), this->k_);
126  }
127  else
128  {
130  }
131 }
132 
133 
134 template<class BasicTurbulenceModel>
137 {
139  this->mesh_.objectRegistry::template
140  lookupObject<uniformDimensionedVectorField>("g");
141 
142  if (mag(g.value()) > SMALL)
143  {
144  vector gHat(g.value()/mag(g.value()));
145 
146  volScalarField v(gHat & this->U_);
148  (
149  mag(this->U_ - gHat*v)
150  + dimensionedScalar("SMALL", dimVelocity, SMALL)
151  );
152 
153  return -fvm::SuSp(this->C1_*tanh(mag(v)/u)*Gcoef(), this->epsilon_);
154  }
155  else
156  {
158  }
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace RASModels
165 } // End namespace Foam
166 
167 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:49
A class for handling words, derived from string.
Definition: word.H:59
BasicTurbulenceModel::rhoField rhoField
Namespace for OpenFOAM.
tmp< volScalarField > Gcoef() const
Calculate the gradient of the given field.
virtual tmp< fvScalarMatrix > epsilonSource() const
Generic dimensioned Type class.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicTurbulenceModel::alphaField alphaField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Macros for easy insertion into run-time selection tables.
virtual void correctNut()
Definition: kEpsilon.C:39
dimensionedScalar tanh(const dimensionedScalar &ds)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
surfaceScalarField & phi
virtual bool read()
Re-read model coefficients if they have changed.
Dimensioned<Type> registered with the database as a registered IOobject which has the functionality o...
const dimensionedVector & g
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Additional buoyancy generation/dissipation term applied to the k and epsilon equations of the standar...
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:64
virtual tmp< fvScalarMatrix > kSource() const
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
U
Definition: pEqn.H:82
const Type & value() const
Return const reference to value.