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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2014-2015 OpenFOAM Foundation
9 Copyright (C) 2018-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "buoyantKEpsilon.H"
30#include "gravityMeshObject.H"
31#include "fvcGrad.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace RASModels
39{
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43template<class BasicTurbulenceModel>
45(
46 const alphaField& alpha,
47 const rhoField& rho,
48 const volVectorField& U,
49 const surfaceScalarField& alphaRhoPhi,
51 const transportModel& transport,
52 const word& propertiesName,
53 const word& type
54)
55:
56 kEpsilon<BasicTurbulenceModel>
57 (
58 alpha,
59 rho,
60 U,
61 alphaRhoPhi,
62 phi,
63 transport,
64 propertiesName,
65 type
66 ),
67
68 Cg_
69 (
70 dimensioned<scalar>::getOrAddToDict
71 (
72 "Cg",
73 this->coeffDict_,
74 1.0
75 )
76 )
77{
78 if (type == typeName)
79 {
80 this->printCoeffs(type);
81 }
82}
83
84
85// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86
87template<class BasicTurbulenceModel>
89{
91 {
92 Cg_.readIfPresent(this->coeffDict());
93
94 return true;
95 }
96
97 return false;
98}
99
100
101template<class BasicTurbulenceModel>
104{
106 meshObjects::gravity::New(this->mesh_.time());
107
108 return
109 (Cg_*this->Cmu_)*this->alpha_*this->k_*(g & fvc::grad(this->rho_))
110 /(this->epsilon_ + this->epsilonMin_);
111}
112
113
114template<class BasicTurbulenceModel>
117{
119 meshObjects::gravity::New(this->mesh_.time());
120
121 if (mag(g.value()) > SMALL)
122 {
123 return -fvm::SuSp(Gcoef(), this->k_);
124 }
125
127}
128
129
130template<class BasicTurbulenceModel>
133{
135 meshObjects::gravity::New(this->mesh_.time());
136
137 if (mag(g.value()) > SMALL)
138 {
139 vector gHat(g.value()/mag(g.value()));
140
141 volScalarField v(gHat & this->U_);
143 (
144 mag(this->U_ - gHat*v)
145 + dimensionedScalar("SMALL", dimVelocity, SMALL)
146 );
147
148 return -fvm::SuSp(this->C1_*tanh(mag(v)/u)*Gcoef(), this->epsilon_);
149 }
150
152}
153
154
155// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156
157} // End namespace RASModels
158} // End namespace Foam
159
160// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
surfaceScalarField & phi
Additional buoyancy generation/dissipation term applied to the k and epsilon equations of the standar...
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > Gcoef() const
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:92
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:69
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:54
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Calculate the gradient of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimVelocity
dimensionedScalar tanh(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
volScalarField & alpha