unityLewisEddyDiffusivity.H
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 Class
25  Foam::turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity
26 
27 Description
28  Eddy-diffusivity based energy gradient heat flux model for RAS or LES
29  of turbulent flow. Specie fluxes are computed assuming a unity turbulent
30  Lewis number.
31 
32 Usage
33  \verbatim
34  LES
35  {
36  model unityLewisEddyDiffusivity;
37  Prt 0.85;
38  }
39  \endverbatim
40 
41 SourceFiles
42  unityLewisEddyDiffusivity.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef unityLewisEddyDiffusivity_H
47 #define unityLewisEddyDiffusivity_H
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace turbulenceThermophysicalTransportModels
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class unityLewisEddyDiffusivity Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class TurbulenceThermophysicalTransportModel>
62 :
63  public TurbulenceThermophysicalTransportModel
64 {
65 
66 protected:
67 
68  // Protected data
69 
70  // Model coefficients
71 
72  //- Turbulent Prandtl number []
74 
75  // Fields
76 
77  //- Turbulent thermal diffusivity of enthalpy [kg/m/s]
79 
80 
81  // Protected Member Functions
82 
83  virtual void correctAlphat();
84 
85  // Used for the implicit energy correction on the temperature laplacian
87  {
88  return volScalarField::New
89  (
90  "alphaEff",
91  this->thermo().kappa()/this->thermo().Cpv() + alphat()
92  );
93  }
94 
95 
96 public:
97 
98  typedef typename TurbulenceThermophysicalTransportModel::alphaField
99  alphaField;
100 
101  typedef typename
104 
105  typedef typename TurbulenceThermophysicalTransportModel::thermoModel
106  thermoModel;
107 
108 
109  //- Runtime type information
110  TypeName("unityLewisEddyDiffusivity");
111 
112 
113  // Constructors
114 
115  //- Construct from a momentum transport model and a thermo model
117  (
118  const momentumTransportModel& momentumTransport,
119  const thermoModel& thermo,
120  const bool allowDefaultPrt = false
121  );
122 
123  //- Construct from a type name, a momentum transport model and a thermo
124  // model, and whether to default the turbulent Prandtl number to one
125  // if it is not specified
127  (
128  const word& type,
129  const momentumTransportModel& momentumTransport,
130  const thermoModel& thermo,
131  const bool allowDefaultPrt = false
132  );
133 
134 
135  //- Destructor
137  {}
138 
139 
140  // Member Functions
141 
142  //- Read thermophysicalTransport dictionary
143  virtual bool read();
144 
145  //- Turbulent thermal diffusivity for enthalpy [kg/m/s]
146  virtual tmp<volScalarField> alphat() const
147  {
148  return alphat_;
149  }
150 
151  //- Turbulent thermal diffusivity for enthalpy for a patch [kg/m/s]
152  virtual tmp<scalarField> alphat(const label patchi) const
153  {
154  return alphat()().boundaryField()[patchi];
155  }
156 
157  //- Effective thermal turbulent conductivity
158  // of mixture [W/m/K]
159  virtual tmp<volScalarField> kappaEff() const
160  {
161  return this->thermo().kappa() + this->thermo().Cp()*alphat();
162  }
163 
164  //- Effective thermal turbulent conductivity
165  // of mixture for patch [W/m/K]
166  virtual tmp<scalarField> kappaEff(const label patchi) const
167  {
168  return
169  this->thermo().kappa().boundaryField()[patchi]
170  + this->thermo().Cp().boundaryField()[patchi]*alphat(patchi);
171  }
172 
173  //- Effective thermal turbulent diffusivity
174  // of mixture for a patch [kg/m/s]
175  virtual tmp<scalarField> alphaEff(const label patchi) const
176  {
177  return
178  this->thermo().kappa().boundaryField()[patchi]
179  /this->thermo().Cpv().boundaryField()[patchi]
180  + alphat(patchi);
181  }
182 
183  //- Effective mass diffusion coefficient
184  // for a given specie mass-fraction [kg/m/s]
185  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const
186  {
187  return volScalarField::New
188  (
189  "DEff",
190  this->thermo().kappa()/this->thermo().Cp() + alphat()
191  );
192  }
193 
194  //- Effective mass diffusion coefficient
195  // for a given specie mass-fraction for patch [kg/m/s]
196  virtual tmp<scalarField> DEff
197  (
198  const volScalarField& Yi,
199  const label patchi
200  ) const
201  {
202  return
203  this->thermo().kappa().boundaryField()[patchi]
204  /this->thermo().Cp().boundaryField()[patchi]
205  + alphat(patchi);
206  }
207 
208  //- Return the heat flux [W/m^2]
209  virtual tmp<surfaceScalarField> q() const;
210 
211  //- Return the patch heat flux [W/m^2]
212  virtual tmp<scalarField> q(const label patchi) const;
213 
214  //- Return the source term for the energy equation
215  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
216 
217  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
218  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
219 
220  //- Return the specie flux
221  // for the given specie mass-fraction for patch [kg/m^2/s]
222  virtual tmp<scalarField> j
223  (
224  const volScalarField& Yi,
225  const label patchi
226  ) const;
227 
228  //- Return the source term for the given specie mass-fraction equation
229  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
230 
231  //- Correct the unityLewisEddyDiffusivity viscosity
232  virtual void predict();
233 };
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace turbulenceThermophysicalTransportModels
239 } // End namespace Foam
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #ifdef NoRepository
244  #include "unityLewisEddyDiffusivity.C"
245 #endif
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #endif
250 
251 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
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,.
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.
volScalarField alphat_
Turbulent thermal diffusivity of enthalpy [kg/m/s].
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual tmp< volScalarField > alphat() const
Turbulent thermal diffusivity for enthalpy [kg/m/s].
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].
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
TypeName("unityLewisEddyDiffusivity")
Runtime type information.
virtual tmp< volScalarField > kappaEff() const
Effective thermal turbulent conductivity.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
compressibleMomentumTransportModel momentumTransportModel
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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
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