unityLewisFourier.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-2021 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::laminarThermophysicalTransportModels::unityLewisFourier
26 
27 Description
28  unityLewisFourier's energy gradient heat flux model for laminar flow.
29  Specie fluxes are computed assuming a unity turbulent Lewis number.
30 
31 SourceFiles
32  unityLewisFourier.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef unityLewisFourier_H
37 #define unityLewisFourier_H
38 
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace laminarThermophysicalTransportModels
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class unityLewisFourier Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 template<class BasicThermophysicalTransportModel>
54 :
56  <
57  BasicThermophysicalTransportModel
58  >
59 {
60 
61 public:
62 
63  typedef typename BasicThermophysicalTransportModel::alphaField
64  alphaField;
65 
68 
69  typedef typename BasicThermophysicalTransportModel::thermoModel
71 
72 
73  //- Runtime type information
74  TypeName("unityLewisFourier");
75 
76 
77  // Constructors
78 
79  //- Construct from a momentum transport model and a thermo model
81  (
82  const momentumTransportModel& momentumTransport,
83  const thermoModel& thermo
84  );
85 
86  //- Construct from a type name, a momentum transport model and a thermo
87  // model
89  (
90  const word& type,
91  const momentumTransportModel& momentumTransport,
92  const thermoModel& thermo
93  );
94 
95 
96  //- Destructor
97  virtual ~unityLewisFourier()
98  {}
99 
100 
101  // Member Functions
102 
103  //- Const access to the coefficients dictionary
104  virtual const dictionary& coeffDict() const;
105 
106  //- Read thermophysicalTransport dictionary
107  virtual bool read();
108 
109  //- Effective mass diffusion coefficient
110  // for a given specie mass-fraction [kg/m/s]
111  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const
112  {
113  return volScalarField::New
114  (
115  "DEff",
116  this->thermo().alpha()
117  );
118  }
119 
120  //- Effective mass diffusion coefficient
121  // for a given specie mass-fraction for patch [kg/m/s]
123  (
124  const volScalarField& Yi,
125  const label patchi
126  ) const
127  {
128  return this->thermo().alpha(patchi);
129  }
130 
131  //- Return the heat flux [W/m^2]
132  virtual tmp<surfaceScalarField> q() const;
133 
134  //- Return the source term for the energy equation
135  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
136 
137  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
138  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
139 
140  //- Return the source term for the given specie mass-fraction equation
141  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
142 
143  //- Correct the unityLewisFourier viscosity
144  virtual void correct();
145 };
146 
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 } // End namespace laminarThermophysicalTransportModels
151 } // End namespace Foam
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 #ifdef NoRepository
156  #include "unityLewisFourier.C"
157 #endif
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 #endif
162 
163 // ************************************************************************* //
Templated abstract base class for laminar thermophysical transport models.
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
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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< surfaceScalarField > q() const
Return the heat flux [W/m^2].
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
CompressibleMomentumTransportModel< dynamicTransportModel > momentumTransportModel
BasicThermophysicalTransportModel::thermoModel thermoModel
A class for handling words, derived from string.
Definition: word.H:59
unityLewisFourier(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
TypeName("unityLewisFourier")
Runtime type information.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
thermo he()
virtual bool read()
Read thermophysicalTransport dictionary.
label patchi
BasicThermophysicalTransportModel::alphaField alphaField
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
unityLewisFourier&#39;s energy gradient heat flux model for laminar flow. Specie fluxes are computed assu...
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Correct the unityLewisFourier viscosity.
Namespace for OpenFOAM.
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.