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-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 "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_
64  (
65  dimensioned<scalar>::lookupOrAddToDict
66  (
67  "Cg",
68  this->coeffDict_,
69  1.0
70  )
71  )
72 {
73  if (type == typeName)
74  {
75  this->printCoeffs(type);
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class BasicMomentumTransportModel>
84 {
86  {
87  Cg_.readIfPresent(this->coeffDict());
88 
89  return true;
90  }
91  else
92  {
93  return false;
94  }
95 }
96 
97 
98 template<class BasicMomentumTransportModel>
101 {
103  this->mesh_.objectRegistry::template
104  lookupObject<uniformDimensionedVectorField>("g");
105 
106  return
107  (Cg_*this->Cmu_)*this->alpha_*this->k_*(g & fvc::grad(this->rho_))
108  /(this->epsilon_ + this->epsilonMin_);
109 }
110 
111 
112 template<class BasicMomentumTransportModel>
115 {
117  this->mesh_.objectRegistry::template
118  lookupObject<uniformDimensionedVectorField>("g");
119 
120  if (mag(g.value()) > small)
121  {
122  return -fvm::SuSp(Gcoef(), this->k_);
123  }
124  else
125  {
127  }
128 }
129 
130 
131 template<class BasicMomentumTransportModel>
134 {
136  this->mesh_.objectRegistry::template
137  lookupObject<uniformDimensionedVectorField>("g");
138 
139  if (mag(g.value()) > small)
140  {
141  vector gHat(g.value()/mag(g.value()));
142 
143  volScalarField v(gHat & this->U_);
145  (
146  mag(this->U_ - gHat*v)
148  );
149 
150  return -fvm::SuSp(this->C1_*tanh(mag(v)/u)*Gcoef(), this->epsilon_);
151  }
152  else
153  {
155  }
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 } // End namespace RASModels
162 } // End namespace Foam
163 
164 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
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
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
Definition: kEpsilon.C:65
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:50
Generic dimensioned Type class.
const Type & value() const
Return const reference to 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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimVelocity
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488