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-2024 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 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
64 template<class TurbulenceThermophysicalTransportModel>
65 bool
67 {
68  if
69  (
71  ::read()
72  )
73  {
74  Sct_.read(this->coeffDict());
75 
76  return true;
77  }
78  else
79  {
80  return false;
81  }
82 }
83 
84 
85 template<class TurbulenceThermophysicalTransportModel>
88 {
90  (
92  (
94  (
95  "q",
96  this->momentumTransport().alphaRhoPhi().group()
97  ),
98  -fvc::interpolate(this->alpha()*this->kappaEff())
99  *fvc::snGrad(this->thermo().T())
100  )
101  );
102 
103  const PtrList<volScalarField>& Y = this->thermo().Y();
104  const volScalarField& p = this->thermo().p();
105  const volScalarField& T = this->thermo().T();
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(this->thermo().hsi(i, p, T));
122 
123  hGradY += fvc::interpolate(hi)*fvc::snGrad(Y[i]);
124  }
125 
126  tmpq.ref() -=
128  (
129  this->alpha()
130  *(
131  this->thermo().kappa()/this->thermo().Cp()
132  + (this->Prt_/Sct_)*this->alphat()
133  )
134  )*hGradY;
135  }
136 
137  return tmpq;
138 }
139 
140 
141 template<class TurbulenceThermophysicalTransportModel>
144 (
145  const label patchi
146 ) const
147 {
148  tmp<scalarField> tmpq
149  (
150  - (
151  this->alpha().boundaryField()[patchi]
152  *this->kappaEff(patchi)
153  *this->thermo().T().boundaryField()[patchi].snGrad()
154  )
155  );
156 
157  const PtrList<volScalarField>& Y = this->thermo().Y();
158  const fvPatchScalarField& p = this->thermo().p().boundaryField()[patchi];
159  const fvPatchScalarField& T = this->thermo().T().boundaryField()[patchi];
160 
161  if (Y.size())
162  {
163  scalarField hGradY(tmpq->size(), scalar(0));
164 
165  forAll(Y, i)
166  {
167  const scalarField hi(this->thermo().hsi(i, p, T));
168 
169  hGradY += hi*Y[i].boundaryField()[patchi].snGrad();
170  }
171 
172  tmpq.ref() -=
173  (
174  this->alpha().boundaryField()[patchi]
175  *(
176  this->thermo().kappa().boundaryField()[patchi]
177  /this->thermo().Cp().boundaryField()[patchi]
178  + (this->Prt_.value()/Sct_.value())*this->alphat(patchi)
179  )
180  )*hGradY;
181  }
182 
183  return tmpq;
184 }
185 
186 
187 template<class TurbulenceThermophysicalTransportModel>
190 (
192 ) const
193 {
194  tmp<fvScalarMatrix> tmpDivq
195  (
196  fvm::Su
197  (
198  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T()),
199  he
200  )
201  );
202 
203  const PtrList<volScalarField>& Y = this->thermo().Y();
204 
205  tmpDivq.ref() -=
206  fvm::laplacianCorrection(this->alpha()*this->alphaEff(), he);
207 
208  surfaceScalarField hGradY
209  (
211  (
212  "hGradY",
213  he.mesh(),
214  dimensionedScalar(he.dimensions()/dimLength, 0)
215  )
216  );
217 
218  forAll(Y, i)
219  {
220  const volScalarField hi
221  (
222  this->thermo().hsi(i, this->thermo().p(), this->thermo().T())
223  );
224 
225  hGradY += fvc::interpolate(hi)*fvc::snGrad(Y[i]);
226  }
227 
228  tmpDivq.ref() -=
229  fvc::div
230  (
232  (
233  this->alpha()
234  *(
235  this->thermo().kappa()/this->thermo().Cp()
236  + (this->Prt_/Sct_)*this->alphat()
237  )
238  )*hGradY*he.mesh().magSf()
239  );
240 
241  return tmpDivq;
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 } // End namespace turbulenceThermophysicalTransportModels
248 } // End namespace Foam
249 
250 // ************************************************************************* //
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:433
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
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:197
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....
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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.
label patchi
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
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimEnergy
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimLength
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
thermo he()
volScalarField & p
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31