Fickian.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) 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::turbulenceThermophysicalTransportModels::Fickian
26 
27 Description
28  Base class for multi-component Fickian based temperature gradient heat
29  flux models with optional Soret thermal diffusion of species.
30 
31  The mixture diffusion coefficients are specified as Function2<scalar>s of
32  pressure and temperature but independent of composition.
33 
34  The heat flux source is implemented as an implicit energy correction to the
35  temperature gradient based flux source. At convergence the energy
36  correction is 0.
37 
38 SourceFiles
39  Fickian.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "Function2.H"
44 
45 #ifndef Fickian_H
46 #define Fickian_H
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class Fickian Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 template<class BasicThermophysicalTransportModel>
58 class Fickian
59 :
60  public BasicThermophysicalTransportModel
61 {
62  // Private data
63 
64  // Model coefficients
65 
66  //- Choice between mass diffusion coefficient functions w.r.t
67  // mixture and mixing binary mass diffusion coefficient functions
68  bool mixtureDiffusionCoefficients_;
69 
70  //- Array of specie binary mass diffusion coefficient functions
71  // [m^2/s]
73 
74  //- List of specie mass diffusion coefficient functions
75  // w.r.t the mixture [m^2/s]
76  PtrList<Function2<scalar>> DmFuncs_;
77 
78  //- List of specie Soret thermal diffusion coefficient
79  // functions [kg/m/s]
80  PtrList<Function2<scalar>> DTFuncs_;
81 
82  //- Mass diffusion coefficients in the mixture
84 
85 
86 public:
87 
88  typedef typename BasicThermophysicalTransportModel::alphaField
89  alphaField;
90 
91  typedef typename
94 
95  typedef typename BasicThermophysicalTransportModel::thermoModel
97 
98 
99  // Constructors
100 
101  //- Construct from a momentum transport model and a thermo model
102  Fickian
103  (
104  const word& type,
105  const momentumTransportModel& momentumTransport,
106  const thermoModel& thermo
107  );
108 
109 
110  //- Destructor
111  virtual ~Fickian()
112  {}
113 
114 
115  // Member Functions
116 
117  //- Read thermophysicalTransport dictionary
118  virtual bool read();
119 
120  //- Effective mass diffusion coefficient
121  // for a given specie mass-fraction [kg/m/s]
122  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const;
123 
124  //- Effective mass diffusion coefficient
125  // for a given specie mass-fraction for patch [kg/m/s]
126  virtual tmp<scalarField> DEff
127  (
128  const volScalarField& Yi,
129  const label patchi
130  ) const;
131 
132  //- Return the heat flux [W/m^2]
133  virtual tmp<surfaceScalarField> q() const;
134 
135  //- Return the source term for the energy equation
136  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
137 
138  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
139  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
140 
141  //- Return the source term for the given specie mass-fraction equation
142  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
143 
144  //- Update the diffusion coefficients
145  virtual void correct();
146 };
147 
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 } // End namespace Foam
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 #ifdef NoRepository
156  #include "Fickian.C"
157 #endif
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 #endif
162 
163 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fickian.C:295
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fickian.C:402
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
BasicThermophysicalTransportModel::thermoModel thermoModel
Definition: Fickian.H:95
virtual ~Fickian()
Destructor.
Definition: Fickian.H:110
Fickian(const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
Definition: Fickian.C:43
compressibleMomentumTransportModel momentumTransportModel
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fickian.C:218
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
Definition: Fickian.H:92
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fickian.C:187
A class for handling words, derived from string.
Definition: word.H:59
virtual void correct()
Update the diffusion coefficients.
Definition: Fickian.C:438
thermo he()
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
Definition: Fickian.C:369
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: Fickian.C:76
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
BasicThermophysicalTransportModel::alphaField alphaField
Definition: Fickian.H:88