buoyantKEpsilon.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) 2014-2024 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 BasicMomentumTransportModel>
42 (
43  const alphaField& alpha,
44  const rhoField& rho,
45  const volVectorField& U,
46  const surfaceScalarField& alphaRhoPhi,
47  const surfaceScalarField& phi,
48  const viscosity& viscosity,
49  const word& type
50 )
51 :
52  kEpsilon<BasicMomentumTransportModel>
53  (
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  viscosity,
60  type
61  ),
62 
63  Cg_("Cg", this->coeffDict(), 1.0)
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 template<class BasicMomentumTransportModel>
71 {
73  {
74  Cg_.readIfPresent(this->coeffDict());
75 
76  return true;
77  }
78  else
79  {
80  return false;
81  }
82 }
83 
84 
85 template<class BasicMomentumTransportModel>
88 {
90  this->mesh_.objectRegistry::template
91  lookupObject<uniformDimensionedVectorField>("g");
92 
93  return
94  (Cg_*this->Cmu_)*this->alpha_*this->k_*(g & fvc::grad(this->rho_))
95  /this->epsilon_;
96 }
97 
98 
99 template<class BasicMomentumTransportModel>
102 {
104  this->mesh_.objectRegistry::template
105  lookupObject<uniformDimensionedVectorField>("g");
106 
107  if (mag(g.value()) > small)
108  {
109  return -fvm::SuSp(Gcoef(), this->k_);
110  }
111  else
112  {
114  }
115 }
116 
117 
118 template<class BasicMomentumTransportModel>
121 {
123  this->mesh_.objectRegistry::template
124  lookupObject<uniformDimensionedVectorField>("g");
125 
126  if (mag(g.value()) > small)
127  {
128  vector gHat(g.value()/mag(g.value()));
129 
130  volScalarField v(gHat & this->U_);
132  (
133  mag(this->U_ - gHat*v)
135  );
136 
137  return -fvm::SuSp(this->C1_*tanh(mag(v)/u)*Gcoef(), this->epsilon_);
138  }
139  else
140  {
142  }
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 } // End namespace RASModels
149 } // End namespace Foam
150 
151 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
tmp< volScalarField > Gcoef() const
buoyantKEpsilon(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.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
virtual bool read()
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: kEpsilon.C:74
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: kEpsilon.C:59
Type & value()
Return a reference to the value.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Calculate the gradient of the given field.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
dimensionedScalar tanh(const dimensionedScalar &ds)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488