All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2022 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 #include "dragModel.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicMomentumTransportModel>
43 (
44  const alphaField& alpha,
45  const rhoField& rho,
46  const volVectorField& U,
47  const surfaceScalarField& alphaRhoPhi,
48  const surfaceScalarField& phi,
49  const viscosity& viscosity,
50  const word& type
51 )
52 :
54  (
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  viscosity,
61  type
62  ),
63 
64  gasTurbulencePtr_(nullptr),
65 
66  alphaInversion_
67  (
69  (
70  "alphaInversion",
71  this->coeffDict_,
72  0.3
73  )
74  ),
75 
76  Cp_
77  (
79  (
80  "Cp",
81  this->coeffDict_,
82  this->Ck_.value()
83  )
84  ),
85 
86  Cmub_
87  (
89  (
90  "Cmub",
91  this->coeffDict_,
92  0.6
93  )
94  )
95 {
96  if (type == typeName)
97  {
98  this->printCoeffs(type);
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class BasicMomentumTransportModel>
107 {
109  {
110  alphaInversion_.readIfPresent(this->coeffDict());
111  Cp_.readIfPresent(this->coeffDict());
112  Cmub_.readIfPresent(this->coeffDict());
113 
114  return true;
115  }
116  else
117  {
118  return false;
119  }
120 }
121 
122 
123 template<class BasicMomentumTransportModel>
126 {
127  if (!gasTurbulencePtr_)
128  {
129  const volVectorField& U = this->U_;
130 
131  const phaseModel& liquid =
132  refCast<const phaseModel>(this->properties());
133  const phaseSystem& fluid = liquid.fluid();
134  const phaseModel& gas = fluid.otherPhase(liquid);
135 
136  gasTurbulencePtr_ =
137  &U.db().lookupObject
138  <
140  >
141  (
143  (
144  momentumTransportModel::typeName,
145  gas.name()
146  )
147  );
148  }
149 
150  return *gasTurbulencePtr_;
151 }
152 
153 
154 template<class BasicMomentumTransportModel>
156 {
157  const phaseCompressibleMomentumTransportModel& gasTurbulence =
158  this->gasTurbulence();
159 
160  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
161  const phaseSystem& fluid = liquid.fluid();
162  const phaseModel& gas = fluid.otherPhase(liquid);
163 
164  this->nut_ =
165  this->Ck_*sqrt(this->k_)*this->delta()
166  + Cmub_*gas.d()*gasTurbulence.alpha()
167  *(mag(this->U_ - gasTurbulence.U()));
168 
169  this->nut_.correctBoundaryConditions();
170  fvConstraints::New(this->mesh_).constrain(this->nut_);
171 }
172 
173 
174 template<class BasicMomentumTransportModel>
176 {
177  const phaseCompressibleMomentumTransportModel& gasTurbulence =
178  this->gasTurbulence();
179 
180  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
181  const phaseSystem& fluid = liquid.fluid();
182  const phaseModel& gas = fluid.otherPhase(liquid);
183 
185  fluid.lookupInterfacialModel<dragModels::dispersedDragModel>
186  (dispersedPhaseInterface(gas, liquid));
187 
188  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
189 
190  tmp<volScalarField> bubbleG
191  (
192  Cp_*sqr(magUr)*drag.K()/liquid.rho()
193  );
194 
195  return bubbleG;
196 }
197 
198 
199 template<class BasicMomentumTransportModel>
202 {
203  const volVectorField& U = this->U_;
204  const alphaField& alpha = this->alpha_;
205  const rhoField& rho = this->rho_;
206 
207  const momentumTransportModel& gasTurbulence = this->gasTurbulence();
208 
209  return
210  (
211  max(alphaInversion_ - alpha, scalar(0))
212  *rho
213  *min
214  (
215  this->Ce_*sqrt(gasTurbulence.k())/this->delta(),
216  1.0/U.time().deltaT()
217  )
218  );
219 }
220 
221 
222 template<class BasicMomentumTransportModel>
224 {
225  const alphaField& alpha = this->alpha_;
226  const rhoField& rho = this->rho_;
227 
228  const phaseCompressibleMomentumTransportModel& gasTurbulence =
229  this->gasTurbulence();
230 
231  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
232 
233  return
234  alpha*rho*bubbleG()
235  + phaseTransferCoeff*gasTurbulence.k()
236  - fvm::Sp(phaseTransferCoeff, this->k_);
237 }
238 
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 } // End namespace LESModels
243 } // End namespace Foam
244 
245 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
scalar delta
Templated abstract base class for multiphase compressible turbulence models.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
U
Definition: pEqn.H:72
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const alphaField & alpha() const
Access function to phase fraction.
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< volScalarField > rho() const =0
Return the density field.
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
const word & name() const
Definition: phaseModel.H:110
Generic dimensioned Type class.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
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.
One equation eddy-viscosity model.
Definition: kEqn.H:71
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:68
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:106
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:175
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:201
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
virtual tmp< fvScalarMatrix > kSource() const
Definition: NicenoKEqn.C:223
Class to represent a interface between phases where one phase is considered dispersed within the othe...
NicenoKEqn(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.
Definition: NicenoKEqn.C:43
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
BasicMomentumTransportModel::rhoField rhoField
Definition: kEqn.H:98
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.
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:318
const phaseSystem & fluid() const
Return the system to which this phase belongs.
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 > &)
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
virtual void correctNut()
Definition: NicenoKEqn.C:155
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
BasicMomentumTransportModel::alphaField alphaField
Definition: kEqn.H:97
Namespace for OpenFOAM.