nonUnityLewisEddyDiffusivity.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-2023 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 "fvcDiv.H"
28 #include "fvcLaplacian.H"
29 #include "fvcSnGrad.H"
30 #include "fvmSup.H"
31 #include "surfaceInterpolate.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace turbulenceThermophysicalTransportModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class TurbulenceThermophysicalTransportModel>
45 (
46  const momentumTransportModel& momentumTransport,
47  const thermoModel& thermo
48 )
49 :
50  unityLewisEddyDiffusivity<TurbulenceThermophysicalTransportModel>
51  (
52  typeName,
53  momentumTransport,
54  thermo,
55  false
56  ),
57 
58  Sct_("Sct", dimless, this->coeffDict_)
59 {
60  this->printCoeffs(typeName);
61 }
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
66 template<class TurbulenceThermophysicalTransportModel>
67 bool
69 {
70  if
71  (
73  ::read()
74  )
75  {
76  Sct_.read(this->coeffDict());
77 
78  return true;
79  }
80  else
81  {
82  return false;
83  }
84 }
85 
86 
87 template<class TurbulenceThermophysicalTransportModel>
90 {
92  (
94  (
96  (
97  "q",
98  this->momentumTransport().alphaRhoPhi().group()
99  ),
100  -fvc::interpolate(this->alpha()*this->kappaEff())
101  *fvc::snGrad(this->thermo().T())
102  )
103  );
104 
105  const PtrList<volScalarField>& Y = this->thermo().Y();
106 
107  if (Y.size())
108  {
109  surfaceScalarField hGradY
110  (
112  (
113  "hGradY",
114  Y[0].mesh(),
116  )
117  );
118 
119  forAll(Y, i)
120  {
121  const volScalarField hi
122  (
123  this->thermo().hsi(i, this->thermo().p(), this->thermo().T())
124  );
125 
126  hGradY += fvc::interpolate(hi)*fvc::snGrad(Y[i]);
127  }
128 
129  tmpq.ref() -=
131  (
132  this->alpha()
133  *(
134  this->thermo().kappa()/this->thermo().Cp()
135  + (this->Prt_/Sct_)*this->alphat()
136  )
137  )*hGradY;
138  }
139 
140  return tmpq;
141 }
142 
143 
144 template<class TurbulenceThermophysicalTransportModel>
147 (
149 ) const
150 {
151  tmp<fvScalarMatrix> tmpDivq
152  (
153  fvm::Su
154  (
155  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
156  he
157  )
158  );
159 
160  const PtrList<volScalarField>& Y = this->thermo().Y();
161 
162  tmpDivq.ref() -=
163  fvm::laplacianCorrection(this->alpha()*this->alphaEff(), he);
164 
165  surfaceScalarField hGradY
166  (
168  (
169  "hGradY",
170  he.mesh(),
171  dimensionedScalar(he.dimensions()/dimLength, 0)
172  )
173  );
174 
175  forAll(Y, i)
176  {
177  const volScalarField hi
178  (
179  this->thermo().hsi(i, this->thermo().p(), this->thermo().T())
180  );
181 
182  hGradY += fvc::interpolate(hi)*fvc::snGrad(Y[i]);
183  }
184 
185  tmpDivq.ref() -=
186  fvc::div
187  (
189  (
190  this->alpha()
191  *(
192  this->thermo().kappa()/this->thermo().Cp()
193  + (this->Prt_/Sct_)*this->alphat()
194  )
195  )*hGradY*he.mesh().magSf()
196  );
197 
198  return tmpDivq;
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace turbulenceThermophysicalTransportModels
205 } // End namespace Foam
206 
207 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
nonUnityLewisEddyDiffusivity(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow....
Calculate the divergence of the given field.
Calculate the laplacian of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
Definition: fvmLaplacian.C:340
Namespace for OpenFOAM.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const dimensionSet dimEnergy
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
thermo he()
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31