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-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 "fvcSnGrad.H"
28 #include "fvmLaplacian.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace turbulenceThermophysicalTransportModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class TurbulenceThermophysicalTransportModel>
42 {
43  alphat_ =
44  this->momentumTransport().rho()
45  *this->momentumTransport().nut()/Prt_;
46  alphat_.correctBoundaryConditions();
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class TurbulenceThermophysicalTransportModel>
55 (
56  const momentumTransportModel& momentumTransport,
57  const thermoModel& thermo,
58  const bool allowDefaultPrt
59 )
60 :
62  (
63  typeName,
64  momentumTransport,
65  thermo,
66  allowDefaultPrt
67  )
68 {}
69 
70 
71 template<class TurbulenceThermophysicalTransportModel>
74 (
75  const word& type,
76  const momentumTransportModel& momentumTransport,
77  const thermoModel& thermo,
78  const bool allowDefaultPrt
79 )
80 :
81  TurbulenceThermophysicalTransportModel
82  (
83  type,
84  momentumTransport,
85  thermo
86  ),
87 
88  Prt_
89  (
90  allowDefaultPrt
91  ? dimensioned<scalar>("Prt", dimless, this->coeffDict(), 1)
92  : dimensioned<scalar>("Prt", dimless, this->coeffDict())
93  ),
94 
95  alphat_
96  (
97  IOobject
98  (
99  IOobject::groupName
100  (
101  "alphat",
102  this->momentumTransport().alphaRhoPhi().group()
103  ),
104  momentumTransport.time().name(),
105  momentumTransport.mesh(),
106  IOobject::MUST_READ,
107  IOobject::AUTO_WRITE
108  ),
109  momentumTransport.mesh()
110  )
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
116 template<class TurbulenceThermophysicalTransportModel>
118 {
120  {
121  Prt_.readIfPresent(this->coeffDict());
122 
123  return true;
124  }
125  else
126  {
127  return false;
128  }
129 }
130 
131 
132 template<class TurbulenceThermophysicalTransportModel>
135 {
137  (
139  (
140  "q",
141  this->momentumTransport().alphaRhoPhi().group()
142  ),
143  -fvc::interpolate(this->alpha()*this->alphaEff())
144  *fvc::snGrad(this->thermo().he())
145  );
146 }
147 
148 
149 template<class TurbulenceThermophysicalTransportModel>
152 (
153  const label patchi
154 ) const
155 {
156  return
157  - (
158  this->alpha().boundaryField()[patchi]
159  *this->alphaEff(patchi)
160  *this->thermo().he().boundaryField()[patchi].snGrad()
161  );
162 }
163 
164 
165 template<class TurbulenceThermophysicalTransportModel>
168 (
170 ) const
171 {
172  return -fvm::laplacian(this->alpha()*this->alphaEff(), he);
173 }
174 
175 
176 template<class TurbulenceThermophysicalTransportModel>
179 (
180  const volScalarField& Yi
181 ) const
182 {
184  (
186  (
187  "j(" + Yi.name() + ')',
188  this->momentumTransport().alphaRhoPhi().group()
189  ),
190  -fvc::interpolate(this->DEff(Yi)*this->alpha())*fvc::snGrad(Yi)
191  );
192 }
193 
194 
195 template<class TurbulenceThermophysicalTransportModel>
198 (
199  const volScalarField& Yi,
200  const label patchi
201 ) const
202 {
203  return
204  - (
205  this->alpha().boundaryField()[patchi]
206  *this->DEff(Yi, patchi)
207  *Yi.boundaryField()[patchi].snGrad()
208  );
209 }
210 
211 
212 template<class TurbulenceThermophysicalTransportModel>
215 (
216  volScalarField& Yi
217 ) const
218 {
219  return -fvm::laplacian(this->alpha()*this->DEff(Yi), Yi);
220 }
221 
222 
223 template<class TurbulenceThermophysicalTransportModel>
225 predict()
226 {
227  TurbulenceThermophysicalTransportModel::predict();
228  correctAlphat();
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace turbulenceThermophysicalTransportModels
235 } // End namespace Foam
236 
237 // ************************************************************************* //
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
Generic dimensioned Type class.
A class for managing temporary objects.
Definition: tmp.H:55
Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow....
unityLewisEddyDiffusivity(const momentumTransportModel &momentumTransport, const thermoModel &thermo, const bool allowDefaultPrt=false)
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.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
virtual void predict()
Correct the unityLewisEddyDiffusivity viscosity.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const char *const group
Group name for atomic constants.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
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 dimless
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()
fluidMulticomponentThermo & thermo
Definition: createFields.H:31