unityLewisEddyDiffusivity.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) 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 
27 #include "fvmLaplacian.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace turbulenceThermophysicalTransportModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class TurbulenceThermophysicalTransportModel>
41 {
42  alphat_ =
43  this->momentumTransport().rho()
44  *this->momentumTransport().nut()/Prt_;
45  alphat_.correctBoundaryConditions();
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 template<class TurbulenceThermophysicalTransportModel>
54 (
55  const momentumTransportModel& momentumTransport,
56  const thermoModel& thermo
57 )
58 :
60  (
61  typeName,
62  momentumTransport,
63  thermo,
64  false
65  )
66 {}
67 
68 
69 template<class TurbulenceThermophysicalTransportModel>
72 (
73  const word& type,
74  const momentumTransportModel& momentumTransport,
75  const thermoModel& thermo,
76  const bool allowDefaultPrt
77 )
78 :
79  TurbulenceThermophysicalTransportModel
80  (
81  type,
82  momentumTransport,
83  thermo
84  ),
85 
86  Prt_
87  (
88  allowDefaultPrt
90  (
91  "Prt",
92  this->coeffDict_,
93  1
94  )
96  (
97  "Prt",
98  dimless,
99  this->coeffDict_
100  )
101  ),
102 
103  alphat_
104  (
105  IOobject
106  (
108  (
109  "alphat",
110  this->momentumTransport().alphaRhoPhi().group()
111  ),
112  momentumTransport.time().timeName(),
113  momentumTransport.mesh(),
116  ),
117  momentumTransport.mesh()
118  )
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
124 template<class TurbulenceThermophysicalTransportModel>
126 {
128  {
129  Prt_.readIfPresent(this->coeffDict());
130 
131  return true;
132  }
133  else
134  {
135  return false;
136  }
137 }
138 
139 
140 template<class TurbulenceThermophysicalTransportModel>
143 {
145  (
147  (
148  "q",
149  this->momentumTransport().alphaRhoPhi().group()
150  ),
151  -fvc::interpolate(this->alphaEff()*this->alpha())
152  *fvc::snGrad(this->thermo().he())
153  );
154 }
155 
156 
157 template<class TurbulenceThermophysicalTransportModel>
160 (
161  volScalarField& he
162 ) const
163 {
164  return -fvm::laplacian(this->alpha()*this->alphaEff(), he);
165 }
166 
167 
168 template<class TurbulenceThermophysicalTransportModel>
171 (
172  const volScalarField& Yi
173 ) const
174 {
176  (
178  (
179  "j(" + Yi.name() + ')',
180  this->momentumTransport().alphaRhoPhi().group()
181  ),
182  -fvc::interpolate(this->DEff(Yi)*this->alpha())*fvc::snGrad(Yi)
183  );
184 }
185 
186 
187 template<class TurbulenceThermophysicalTransportModel>
190 (
191  volScalarField& Yi
192 ) const
193 {
194  return -fvm::laplacian(this->alpha()*this->DEff(Yi), Yi);
195 }
196 
197 
198 template<class TurbulenceThermophysicalTransportModel>
201 {
203  correctAlphat();
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 } // End namespace turbulenceThermophysicalTransportModels
210 } // End namespace Foam
211 
212 // ************************************************************************* //
const char *const group
Group name for atomic constants.
const word & name() const
Return name.
Definition: IOobject.H:303
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
Calculate the matrix for the laplacian of the field.
Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow...
const dimensionSet dimless
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
static word groupName(Name name, const word &group)
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
thermo he()
virtual void correct()
Correct the unityLewisEddyDiffusivity viscosity.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
unityLewisEddyDiffusivity(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.