All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2022 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  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
86  // Used for the implicit energy correction on the temperature laplacian
87  virtual tmp<volScalarField> alphaEff() const
88  {
89  return this->thermo().kappa()/this->thermo().Cpv() + alphat();
90  }
91 
92 
93 public:
94 
95  typedef typename TurbulenceThermophysicalTransportModel::alphaField
96  alphaField;
97 
98  typedef typename
101 
102  typedef typename TurbulenceThermophysicalTransportModel::thermoModel
103  thermoModel;
104 
105 
106  //- Runtime type information
107  TypeName("unityLewisEddyDiffusivity");
108 
109 
110  // Constructors
111 
112  //- Construct from a momentum transport model and a thermo model
114  (
115  const momentumTransportModel& momentumTransport,
116  const thermoModel& thermo,
117  const bool allowDefaultPrt = false
118  );
119 
120  //- Construct from a type name, a momentum transport model and a thermo
121  // model, and whether to default the turbulent Prandtl number to one
122  // if it is not specified
124  (
125  const word& type,
126  const momentumTransportModel& momentumTransport,
127  const thermoModel& thermo,
128  const bool allowDefaultPrt = false
129  );
130 
131 
132  //- Destructor
134  {}
135 
136 
137  // Member Functions
138 
139  //- Read thermophysicalTransport dictionary
140  virtual bool read();
141 
142  //- Turbulent thermal diffusivity for enthalpy [kg/m/s]
143  virtual tmp<volScalarField> alphat() const
144  {
145  return alphat_;
146  }
147 
148  //- Turbulent thermal diffusivity for enthalpy for a patch [kg/m/s]
149  virtual tmp<scalarField> alphat(const label patchi) const
150  {
151  return alphat()().boundaryField()[patchi];
152  }
153 
154  //- Effective thermal turbulent conductivity
155  // of mixture [W/m/K]
156  virtual tmp<volScalarField> kappaEff() const
157  {
158  return this->thermo().kappa() + this->thermo().Cp()*alphat();
159  }
160 
161  //- Effective thermal turbulent conductivity
162  // of mixture for patch [W/m/K]
163  virtual tmp<scalarField> kappaEff(const label patchi) const
164  {
165  return
166  this->thermo().kappa().boundaryField()[patchi]
167  + this->thermo().Cp().boundaryField()[patchi]*alphat(patchi);
168  }
169 
170  //- Effective mass diffusion coefficient
171  // for a given specie mass-fraction [kg/m/s]
172  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const
173  {
174  return volScalarField::New
175  (
176  "DEff",
177  this->thermo().kappa()/this->thermo().Cp() + alphat()
178  );
179  }
180 
181  //- Effective mass diffusion coefficient
182  // for a given specie mass-fraction for patch [kg/m/s]
183  virtual tmp<scalarField> DEff
184  (
185  const volScalarField& Yi,
186  const label patchi
187  ) const
188  {
189  return
190  this->thermo().kappa().boundaryField()[patchi]
191  /this->thermo().Cp().boundaryField()[patchi]
192  + alphat(patchi);
193  }
194 
195  //- Return the heat flux [W/m^2]
196  virtual tmp<surfaceScalarField> q() const;
197 
198  //- Return the source term for the energy equation
199  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
200 
201  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
202  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
203 
204  //- Return the source term for the given specie mass-fraction equation
205  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
206 
207  //- Correct the unityLewisEddyDiffusivity viscosity
208  virtual void predict();
209 };
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace turbulenceThermophysicalTransportModels
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #ifdef NoRepository
220  #include "unityLewisEddyDiffusivity.C"
221 #endif
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #endif
226 
227 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
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,.
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 > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
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