NicenoKEqn.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) 2013-2018 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 "NicenoKEqn.H"
27 #include "fvOptions.H"
28 #include "twoPhaseSystem.H"
29 #include "dragModel.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace LESModels
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_(nullptr),
66 
67  alphaInversion_
68  (
70  (
71  "alphaInversion",
72  this->coeffDict_,
73  0.3
74  )
75  ),
76 
77  Cp_
78  (
80  (
81  "Cp",
82  this->coeffDict_,
83  this->Ck_.value()
84  )
85  ),
86 
87  Cmub_
88  (
90  (
91  "Cmub",
92  this->coeffDict_,
93  0.6
94  )
95  )
96 {
97  if (type == typeName)
98  {
99  this->printCoeffs(type);
100  }
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class BasicTurbulenceModel>
108 {
110  {
111  alphaInversion_.readIfPresent(this->coeffDict());
112  Cp_.readIfPresent(this->coeffDict());
113  Cmub_.readIfPresent(this->coeffDict());
114 
115  return true;
116  }
117  else
118  {
119  return false;
120  }
121 }
122 
123 
124 template<class BasicTurbulenceModel>
126 <
127  typename BasicTurbulenceModel::transportModel
128 >&
130 {
131  if (!gasTurbulencePtr_)
132  {
133  const volVectorField& U = this->U_;
134 
135  const transportModel& liquid = this->transport();
136  const twoPhaseSystem& fluid =
137  refCast<const twoPhaseSystem>(liquid.fluid());
138  const transportModel& gas = fluid.otherPhase(liquid);
139 
140  gasTurbulencePtr_ =
141  &U.db()
143  (
145  (
147  gas.name()
148  )
149  );
150  }
151 
152  return *gasTurbulencePtr_;
153 }
154 
155 
156 template<class BasicTurbulenceModel>
158 {
160  this->gasTurbulence();
161 
162  this->nut_ =
163  this->Ck_*sqrt(this->k_)*this->delta()
164  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
165  *(mag(this->U_ - gasTurbulence.U()));
166 
167  this->nut_.correctBoundaryConditions();
168  fv::options::New(this->mesh_).correct(this->nut_);
169 
170  BasicTurbulenceModel::correctNut();
171 }
172 
173 
174 template<class BasicTurbulenceModel>
176 {
178  this->gasTurbulence();
179 
180  const transportModel& liquid = this->transport();
181  const twoPhaseSystem& fluid =
182  refCast<const twoPhaseSystem>(liquid.fluid());
183  const transportModel& gas = fluid.otherPhase(liquid);
184 
185  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
186 
187  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
188 
189  tmp<volScalarField> bubbleG
190  (
191  Cp_*sqr(magUr)*drag.K()/liquid.rho()
192  );
193 
194  return bubbleG;
195 }
196 
197 
198 template<class BasicTurbulenceModel>
201 {
202  const volVectorField& U = this->U_;
203  const alphaField& alpha = this->alpha_;
204  const rhoField& rho = this->rho_;
205 
206  const turbulenceModel& gasTurbulence = this->gasTurbulence();
207 
208  return
209  (
210  max(alphaInversion_ - alpha, scalar(0))
211  *rho
212  *min
213  (
214  this->Ce_*sqrt(gasTurbulence.k())/this->delta(),
215  1.0/U.time().deltaT()
216  )
217  );
218 }
219 
220 
221 template<class BasicTurbulenceModel>
223 {
224  const alphaField& alpha = this->alpha_;
225  const rhoField& rho = this->rho_;
226 
228  this->gasTurbulence();
229 
230  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
231 
232  return
233  alpha*rho*bubbleG()
234  + phaseTransferCoeff*gasTurbulence.k()
235  - fvm::Sp(phaseTransferCoeff, this->k_);
236 }
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 } // End namespace LESModels
242 } // End namespace Foam
243 
244 // ************************************************************************* //
scalar delta
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
BasicTurbulenceModel::alphaField alphaField
Definition: NicenoKEqn.H:115
surfaceScalarField & phi
multiphaseSystem & fluid
Definition: createFields.H:11
const volVectorField & U() const
Access function to velocity field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicTurbulenceModel::rhoField rhoField
Definition: NicenoKEqn.H:116
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Abstract base class for turbulence models (RAS, LES and laminar).
Info<< "Reading strained laminar flame speed field Su\"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Definition: createFields.H:182
One equation eddy-viscosity model.
Definition: kEqn.H:71
Class which solves the volume fraction equations for two phases.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const transportModel & transport() const
Access function to incompressible transport model.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:107
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:175
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:200
virtual tmp< fvScalarMatrix > kSource() const
Definition: NicenoKEqn.C:222
NicenoKEqn(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: NicenoKEqn.C:42
Templated abstract base class for multiphase compressible turbulence models.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
U
Definition: pEqn.H:72
One-equation SGS model for the continuous phase in a two-phase system including bubble-generated turb...
Definition: NicenoKEqn.H:72
const Time & time() const
Return time.
Definition: IOobject.C:360
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:354
BasicTurbulenceModel::transportModel transportModel
Definition: NicenoKEqn.H:117
virtual tmp< volScalarField > K(const volScalarField &Ur) const =0
The dragfunction K used in the momentum eq.
virtual void correctNut()
Definition: NicenoKEqn.C:157
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
Namespace for OpenFOAM.