kOmegaSSTSato.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 "kOmegaSSTSato.H"
28 #include "twoPhaseSystem.H"
29 #include "dragModel.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  gasTurbulencePtr_(NULL),
66 
67  Cmub_
68  (
70  (
71  "Cmub",
72  this->coeffDict_,
73  0.6
74  )
75  )
76 {
77  if (type == typeName)
78  {
79  this->printCoeffs(type);
80  }
81 }
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
86 template<class BasicTurbulenceModel>
88 {
90  {
91  Cmub_.readIfPresent(this->coeffDict());
92 
93  return true;
94  }
95  else
96  {
97  return false;
98  }
99 }
100 
101 template<class BasicTurbulenceModel>
103 <
104  typename BasicTurbulenceModel::transportModel
105 >&
107 {
108  if (!gasTurbulencePtr_)
109  {
110  const volVectorField& U = this->U_;
111 
112  const transportModel& liquid = this->transport();
113  const twoPhaseSystem& fluid =
114  refCast<const twoPhaseSystem>(liquid.fluid());
115  const transportModel& gas = fluid.otherPhase(liquid);
116 
117  gasTurbulencePtr_ =
118  &U.db()
120  (
122  (
124  gas.name()
125  )
126  );
127  }
128 
129  return *gasTurbulencePtr_;
130 }
131 
132 
133 template<class BasicTurbulenceModel>
135 {
137  this->gasTurbulence();
138 
140  (
141  pow(this->betaStar_, 0.25)*this->y_*sqrt(this->k_)/this->nu()
142  );
143 
144  this->nut_ =
145  this->a1_*this->k_
146  /max
147  (
148  this->a1_*this->omega_,
149  this->F23()*sqrt(2.0)*mag(symm(fvc::grad(this->U_)))
150  )
151  + sqr(1 - exp(-yPlus/16.0))
152  *Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
153  *(mag(this->U_ - gasTurbulence.U()));
154 
155  this->nut_.correctBoundaryConditions();
156 }
157 
158 
159 template<class BasicTurbulenceModel>
161 {
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace RASModels
169 } // End namespace Foam
170 
171 // ************************************************************************* //
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::RASModels::kOmegaSSTSato::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTSato.H:162
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::twoPhaseSystem
Class which solves the volume fraction equations for two phases.
Definition: twoPhaseSystem.H:51
Foam::RASModels::kOmegaSST
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Definition: kOmegaSST.H:126
Foam::PhaseCompressibleTurbulenceModel
Templated abstract base class for multiphase compressible turbulence models.
Definition: PhaseCompressibleTurbulenceModel.H:51
Foam::RASModels::kOmegaSSTSato::kOmegaSSTSato
kOmegaSSTSato(const kOmegaSSTSato &)
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::RASModels::kOmegaSSTSato::read
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTSato.C:87
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::alpha
const alphaField & alpha() const
Access function to phase fraction.
Definition: TurbulenceModel.H:148
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::RASModels::kOmegaSSTSato::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTSato.C:160
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
kOmegaSSTSato.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::RASModels::kOmegaSSTSato::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTSato.H:163
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::transport
const transportModel & transport() const
Access function to incompressible transport model.
Definition: TurbulenceModel.H:154
Foam::RASModels::kOmegaSSTSato::correctNut
virtual void correctNut()
Definition: kOmegaSSTSato.C:134
Foam::kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:385
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::RASModels::kOmegaSSTSato::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTSato.H:161
Foam::RASModels::kOmegaSSTSato::gasTurbulence
const PhaseCompressibleTurbulenceModel< typename BasicTurbulenceModel::transportModel > & gasTurbulence() const
Return the turbulence model for the gas phase.
Definition: kOmegaSSTSato.C:106
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::yPlus
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
Definition: yPlus.H:110
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:142