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-2020 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 "phaseSystem.H"
29 #include "dragModel.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace LESModels
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 transportModel& transport,
49  const word& type
50 )
51 :
53  (
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  transport,
60  type
61  ),
62 
63  gasTurbulencePtr_(nullptr),
64 
65  alphaInversion_
66  (
68  (
69  "alphaInversion",
70  this->coeffDict_,
71  0.3
72  )
73  ),
74 
75  Cp_
76  (
78  (
79  "Cp",
80  this->coeffDict_,
81  this->Ck_.value()
82  )
83  ),
84 
85  Cmub_
86  (
88  (
89  "Cmub",
90  this->coeffDict_,
91  0.6
92  )
93  )
94 {
95  if (type == typeName)
96  {
97  this->printCoeffs(type);
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class BasicMomentumTransportModel>
106 {
108  {
109  alphaInversion_.readIfPresent(this->coeffDict());
110  Cp_.readIfPresent(this->coeffDict());
111  Cmub_.readIfPresent(this->coeffDict());
112 
113  return true;
114  }
115  else
116  {
117  return false;
118  }
119 }
120 
121 
122 template<class BasicMomentumTransportModel>
124 <
125  typename BasicMomentumTransportModel::transportModel
126 >&
128 {
129  if (!gasTurbulencePtr_)
130  {
131  const volVectorField& U = this->U_;
132 
133  const transportModel& liquid = this->transport();
134  const phaseSystem& fluid = liquid.fluid();
135  const transportModel& gas = fluid.otherPhase(liquid);
136 
137  gasTurbulencePtr_ =
138  &U.db().lookupObject
139  <
141  >
142  (
144  (
145  momentumTransportModel::typeName,
146  gas.name()
147  )
148  );
149  }
150 
151  return *gasTurbulencePtr_;
152 }
153 
154 
155 template<class BasicMomentumTransportModel>
157 {
159  gasTurbulence = this->gasTurbulence();
160 
161  this->nut_ =
162  this->Ck_*sqrt(this->k_)*this->delta()
163  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
164  *(mag(this->U_ - gasTurbulence.U()));
165 
166  this->nut_.correctBoundaryConditions();
167  fv::options::New(this->mesh_).correct(this->nut_);
168 }
169 
170 
171 template<class BasicMomentumTransportModel>
173 {
175  gasTurbulence = this->gasTurbulence();
176 
177  const transportModel& liquid = this->transport();
178  const phaseSystem& fluid = liquid.fluid();
179  const transportModel& gas = fluid.otherPhase(liquid);
180 
181  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
182 
183  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
184 
185  tmp<volScalarField> bubbleG
186  (
187  Cp_*sqr(magUr)*drag.K()/liquid.rho()
188  );
189 
190  return bubbleG;
191 }
192 
193 
194 template<class BasicMomentumTransportModel>
197 {
198  const volVectorField& U = this->U_;
199  const alphaField& alpha = this->alpha_;
200  const rhoField& rho = this->rho_;
201 
202  const momentumTransportModel& gasTurbulence = this->gasTurbulence();
203 
204  return
205  (
206  max(alphaInversion_ - alpha, scalar(0))
207  *rho
208  *min
209  (
210  this->Ce_*sqrt(gasTurbulence.k())/this->delta(),
211  1.0/U.time().deltaT()
212  )
213  );
214 }
215 
216 
217 template<class BasicMomentumTransportModel>
219 {
220  const alphaField& alpha = this->alpha_;
221  const rhoField& rho = this->rho_;
222 
224  gasTurbulence = this->gasTurbulence();
225 
226  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
227 
228  return
229  alpha*rho*bubbleG()
230  + phaseTransferCoeff*gasTurbulence.k()
231  - fvm::Sp(phaseTransferCoeff, this->k_);
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace LESModels
238 } // End namespace Foam
239 
240 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
Templated abstract base class for multiphase compressible turbulence models.
scalar delta
const transportModel & transport() const
Access function to incompressible transport model.
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)
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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:192
One equation eddy-viscosity model.
Definition: kEqn.H:71
phi
Definition: pEqn.H:104
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:105
A class for handling words, derived from string.
Definition: word.H:59
NicenoKEqn(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: NicenoKEqn.C:42
static word groupName(Name name, const word &group)
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:172
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:196
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:112
phaseSystem & fluid
Definition: createFields.H:11
virtual tmp< fvScalarMatrix > kSource() const
Definition: NicenoKEqn.C:218
BasicMomentumTransportModel::rhoField rhoField
Definition: kEqn.H:98
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
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:328
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:322
BasicMomentumTransportModel::transportModel transportModel
Definition: kEqn.H:99
virtual void correctNut()
Definition: NicenoKEqn.C:156
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
BasicMomentumTransportModel::alphaField alphaField
Definition: kEqn.H:97
Namespace for OpenFOAM.