continuousGasKEpsilon.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 "continuousGasKEpsilon.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 #include "dragModel.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicMomentumTransportModel>
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const viscosity& viscosity,
51  const word& type
52 )
53 :
55  (
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  viscosity,
62  type
63  ),
64 
65  liquidTurbulencePtr_(nullptr),
66 
67  nutEff_
68  (
69  IOobject
70  (
71  IOobject::groupName("nutEff", alphaRhoPhi.group()),
72  this->runTime_.timeName(),
73  this->mesh_,
76  ),
77  this->nut_
78  ),
79 
80  alphaInversion_
81  (
83  (
84  "alphaInversion",
85  this->coeffDict_,
86  0.7
87  )
88  )
89 {
90  if (type == typeName)
91  {
92  this->printCoeffs(type);
93  }
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class BasicMomentumTransportModel>
101 {
103  {
104  alphaInversion_.readIfPresent(this->coeffDict());
105 
106  return true;
107  }
108  else
109  {
110  return false;
111  }
112 }
113 
114 
115 template<class BasicMomentumTransportModel>
117 {
119 
120  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
121 
122  const phaseModel& gas = refCast<const phaseModel>(this->properties());
123  const phaseSystem& fluid = gas.fluid();
124  const phaseModel& liquid = fluid.otherPhase(gas);
125 
127  fluid.lookupInterfacialModel
129  (dispersedPhaseInterface(gas, liquid));
130 
131  volScalarField thetal(liquidTurbulence.k()/liquidTurbulence.epsilon());
132  volScalarField rhodv(gas.rho() + virtualMass.Cvm()*liquid.rho());
133  volScalarField thetag
134  (
135  (rhodv/(18*liquid.rho()*liquid.thermo().nu()))*sqr(gas.d())
136  );
137  volScalarField expThetar
138  (
139  min
140  (
141  exp(min(thetal/thetag, scalar(50))),
142  scalar(1)
143  )
144  );
145  volScalarField omega((1 - expThetar)/(1 + expThetar));
146 
147  nutEff_ = omega*liquidTurbulence.nut();
148  fvConstraints::New(this->mesh_).constrain(nutEff_);
149 }
150 
151 
152 template<class BasicMomentumTransportModel>
155 {
156  if (!liquidTurbulencePtr_)
157  {
158  const volVectorField& U = this->U_;
159 
160  const phaseModel& gas = refCast<const phaseModel>(this->properties());
161  const phaseSystem& fluid = gas.fluid();
162  const phaseModel& liquid = fluid.otherPhase(gas);
163 
164  liquidTurbulencePtr_ =
166  (
168  (
169  momentumTransportModel::typeName,
170  liquid.name()
171  )
172  );
173  }
174 
175  return *liquidTurbulencePtr_;
176 }
177 
178 
179 template<class BasicMomentumTransportModel>
182 {
183  volScalarField blend
184  (
185  max
186  (
187  min
188  (
189  (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
190  scalar(1)
191  ),
192  scalar(0)
193  )
194  );
195 
196  return volScalarField::New
197  (
198  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
199  blend*this->nut_
200  + (1.0 - blend)*rhoEff()*nutEff_/this->rho_
201  + this->nu()
202  );
203 }
204 
205 
206 template<class BasicMomentumTransportModel>
209 {
210  const phaseModel& gas = refCast<const phaseModel>(this->properties());
211  const phaseSystem& fluid = gas.fluid();
212  const phaseModel& liquid = fluid.otherPhase(gas);
213 
215  fluid.lookupInterfacialModel
217  (dispersedPhaseInterface(gas, liquid));
218 
219  return volScalarField::New
220  (
221  IOobject::groupName("rhoEff", this->alphaRhoPhi_.group()),
222  gas.rho() + (virtualMass.Cvm() + 3.0/20.0)*liquid.rho()
223  );
224 }
225 
226 
227 template<class BasicMomentumTransportModel>
230 {
231  const volVectorField& U = this->U_;
232  const alphaField& alpha = this->alpha_;
233  const rhoField& rho = this->rho_;
234 
235  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
236 
237  return
238  (
239  max(alphaInversion_ - alpha, scalar(0))
240  *rho
241  *min
242  (
243  liquidTurbulence.epsilon()/liquidTurbulence.k(),
244  1.0/U.time().deltaT()
245  )
246  );
247 }
248 
249 
250 template<class BasicMomentumTransportModel>
253 {
254  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
255  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
256 
257  return
258  phaseTransferCoeff*liquidTurbulence.k()
259  - fvm::Sp(phaseTransferCoeff, this->k_);
260 }
261 
262 
263 template<class BasicMomentumTransportModel>
266 {
267  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
268  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
269 
270  return
271  phaseTransferCoeff*liquidTurbulence.epsilon()
272  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
273 }
274 
275 
276 template<class BasicMomentumTransportModel>
279 {
280  tmp<volScalarField> tk(this->k());
281 
283  (
284  IOobject::groupName("R", this->alphaRhoPhi_.group()),
285  ((2.0/3.0)*I)*tk() - (nutEff_)*dev(twoSymm(fvc::grad(this->U_))),
286  tk().boundaryField().types()
287  );
288 }
289 
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace RASModels
294 } // End namespace Foam
295 
296 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
U
Definition: pEqn.H:72
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > rho() const =0
Return the density field.
const word & name() const
Definition: phaseModel.H:110
virtual bool read()
Re-read model coefficients if they have changed.
label k
Boltzmann constant.
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.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
static const Identity< scalar > I
Definition: Identity.H:93
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:83
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:68
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
A class for handling words, derived from string.
Definition: word.H:59
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
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))
BasicMomentumTransportModel::alphaField alphaField
Definition: kEpsilon.H:115
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
continuousGasKEpsilon(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.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
const momentumTransportModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
tmp< volScalarField > phaseTransferCoeff() const
const Time & time() const
Return time.
Definition: IOobject.C:318
const phaseSystem & fluid() const
Return the system to which this phase belongs.
BasicMomentumTransportModel::rhoField rhoField
Definition: kEpsilon.H:116
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:73
virtual void correctNut()
Definition: kEpsilon.C:41
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
Namespace for OpenFOAM.