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 
86 public:
87 
88  typedef typename TurbulenceThermophysicalTransportModel::alphaField
89  alphaField;
90 
91  typedef typename
94 
95  typedef typename TurbulenceThermophysicalTransportModel::thermoModel
97 
98 
99  //- Runtime type information
100  TypeName("unityLewisEddyDiffusivity");
101 
102 
103  // Constructors
104 
105  //- Construct from a momentum transport model and a thermo model
107  (
108  const momentumTransportModel& momentumTransport,
109  const thermoModel& thermo
110  );
111 
112  //- Construct from a type name, a momentum transport model and a thermo
113  // model, and whether to default the turbulent Prandtl number to one
114  // if it is not specified
116  (
117  const word& type,
118  const momentumTransportModel& momentumTransport,
119  const thermoModel& thermo,
120  const bool allowDefaultPrt = false
121  );
122 
123 
124  //- Destructor
126  {}
127 
128 
129  // Member Functions
130 
131  //- Read thermophysicalTransport dictionary
132  virtual bool read();
133 
134  //- Turbulent thermal diffusivity for enthalpy [kg/m/s]
135  virtual tmp<volScalarField> alphat() const
136  {
137  return alphat_;
138  }
139 
140  //- Turbulent thermal diffusivity for enthalpy for a patch [kg/m/s]
141  virtual tmp<scalarField> alphat(const label patchi) const
142  {
143  return alphat()().boundaryField()[patchi];
144  }
145 
146  //- Effective thermal turbulent diffusivity for temperature
147  // of mixture [W/m/K]
148  virtual tmp<volScalarField> kappaEff() const
149  {
150  return this->thermo().kappaEff(alphat());
151  }
152 
153  //- Effective thermal turbulent diffusivity for temperature
154  // of mixture for patch [W/m/K]
155  virtual tmp<scalarField> kappaEff(const label patchi) const
156  {
157  return this->thermo().kappaEff
158  (
159  alphat(patchi),
160  patchi
161  );
162  }
163 
164  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
165  virtual tmp<volScalarField> alphaEff() const
166  {
167  return this->thermo().alphaEff(alphat());
168  }
169 
170  //- Effective thermal turbulent diffusivity of mixture
171  // for patch [kg/m/s]
172  virtual tmp<scalarField> alphaEff(const label patchi) const
173  {
174  return this->thermo().alphaEff
175  (
176  alphat(patchi),
177  patchi
178  );
179  }
180 
181  //- Effective mass diffusion coefficient
182  // for a given specie mass-fraction [kg/m/s]
183  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const
184  {
185  return volScalarField::New
186  (
187  "DEff",
188  this->thermo().kappa()/this->thermo().Cp() + alphat()
189  );
190  }
191 
192  //- Effective mass diffusion coefficient
193  // for a given specie mass-fraction for patch [kg/m/s]
195  (
196  const volScalarField& Yi,
197  const label patchi
198  ) const
199  {
200  return
201  this->thermo().kappa().boundaryField()[patchi]
202  /this->thermo().Cp().boundaryField()[patchi]
203  + alphat(patchi);
204  }
205 
206  //- Return the heat flux [W/m^2]
207  virtual tmp<surfaceScalarField> q() const;
208 
209  //- Return the source term for the energy equation
210  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
211 
212  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
213  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
214 
215  //- Return the source term for the given specie mass-fraction equation
216  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
217 
218  //- Correct the unityLewisEddyDiffusivity viscosity
219  virtual void correct();
220 };
221 
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 } // End namespace turbulenceThermophysicalTransportModels
226 } // End namespace Foam
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 #ifdef NoRepository
231  #include "unityLewisEddyDiffusivity.C"
232 #endif
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 #endif
237 
238 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
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 > kappaEff() const
Effective thermal turbulent diffusivity for temperature.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Eddy-diffusivity based energy gradient heat flux model for RAS or LES of turbulent flow...
virtual tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
compressibleMomentumTransportModel momentumTransportModel
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
thermo he()
virtual void correct()
Correct the unityLewisEddyDiffusivity viscosity.
label patchi
unityLewisEddyDiffusivity(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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField alphat_
Turbulent thermal diffusivity of enthalpy [kg/m/s].
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
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.
TypeName("unityLewisEddyDiffusivity")
Runtime type information.
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.