Fourier.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::Fourier
26 
27 Description
28  Fourier's temperature gradient heat flux model
29  for single specie laminar flow.
30 
31  The heat flux source is implemented as an implicit energy correction to the
32  temperature gradient based flux source. At convergence the energy
33  correction is 0.
34 
35 SourceFiles
36  Fourier.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef Fourier_H
41 #define Fourier_H
42 
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace laminarThermophysicalTransportModels
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class Fourier Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class BasicThermophysicalTransportModel>
57 class Fourier
58 :
60  <
61  BasicThermophysicalTransportModel
62  >
63 {
64 
65 public:
66 
67  typedef typename BasicThermophysicalTransportModel::alphaField
68  alphaField;
69 
72 
73  typedef typename BasicThermophysicalTransportModel::thermoModel
75 
76 
77  //- Runtime type information
78  TypeName("Fourier");
79 
80 
81  // Constructors
82 
83  //- Construct from components
84  Fourier
85  (
86  const momentumTransportModel& momentumTransport,
87  const thermoModel& thermo
88  );
89 
90 
91  //- Destructor
92  virtual ~Fourier()
93  {}
94 
95 
96  // Member Functions
97 
98  //- Const access to the coefficients dictionary
99  virtual const dictionary& coeffDict() const;
100 
101  //- Read thermophysicalTransport dictionary
102  virtual bool read();
103 
104  //- Effective mass diffusion coefficient
105  // for a given specie mass-fraction [kg/m/s]
106  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const;
107 
108  //- Effective mass diffusion coefficient
109  // for a given specie mass-fraction for patch [kg/m/s]
110  virtual tmp<scalarField> DEff
111  (
112  const volScalarField& Yi,
113  const label patchi
114  ) const;
115 
116  //- Return the heat flux [W/m^2]
117  virtual tmp<surfaceScalarField> q() const;
118 
119  //- Return the source term for the energy equation
120  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
121 
122  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
123  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
124 
125  //- Return the source term for the given specie mass-fraction equation
126  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
127 
128  //- Correct the Fourier viscosity
129  virtual void correct();
130 };
131 
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 } // End namespace laminarThermophysicalTransportModels
136 } // End namespace Foam
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 #ifdef NoRepository
141  #include "Fourier.C"
142 #endif
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 #endif
147 
148 // ************************************************************************* //
Templated abstract base class for laminar thermophysical transport models.
fluidReactionThermo & thermo
Definition: createFields.H:28
Fourier(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from components.
Definition: Fourier.C:41
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
Definition: Fourier.H:70
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: Fourier.C:59
BasicThermophysicalTransportModel::thermoModel thermoModel
Definition: Fourier.H:73
compressibleMomentumTransportModel momentumTransportModel
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fourier.C:156
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: Fourier.C:66
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
Definition: Fourier.C:140
virtual void correct()
Correct the Fourier viscosity.
Definition: Fourier.C:169
thermo he()
label patchi
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fourier.C:75
TypeName("Fourier")
Runtime type information.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fourier.C:108
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fourier.C:126
A class for managing temporary objects.
Definition: PtrList.H:53
Fourier&#39;s temperature gradient heat flux model for single specie laminar flow.
Definition: Fourier.H:56
BasicThermophysicalTransportModel::alphaField alphaField
Definition: Fourier.H:67
Namespace for OpenFOAM.