MaxwellStefan.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::MaxwellStefan
26 
27 Description
28  Base class for multi-component Maxwell Stefan generalised Fick's law
29  diffusion coefficients based temperature gradient heat flux model with
30  optional Soret thermal diffusion of species.
31 
32  The binary mass diffusion coefficients are specified as Function2<scalar>s
33  of pressure and temperature but independent of composition.
34 
35  The heat flux source is implemented as an implicit energy correction to the
36  temperature gradient based flux source. At convergence the energy
37  correction is 0.
38 
39  References:
40  \verbatim
41  Taylor, R., & Krishna, R. (1993).
42  Multicomponent mass transfer (Vol. 2).
43  John Wiley & Sons.
44 
45  Merk, H. J. (1959).
46  The macroscopic equations for simultaneous heat and mass transfer
47  in isotropic, continuous and closed systems.
48  Applied Scientific Research,
49  Section A, 8(1), 73-99.
50  \endverbatim
51 
52 SourceFiles
53  MaxwellStefan.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #include "Function2.H"
58 #include "LUscalarMatrix.H"
59 
60 #ifndef MaxwellStefan_H
61 #define MaxwellStefan_H
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class MaxwellStefan Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 template<class BasicThermophysicalTransportModel>
73 class MaxwellStefan
74 :
75  public BasicThermophysicalTransportModel
76 {
77  // Private data
78 
79  // Model coefficients
80 
81  //- Array of specie binary mass diffusion coefficient functions
82  // [m^2/s]
84 
85  //- List of specie Soret thermal diffusion coefficient
86  // functions [kg/m/s]
87  PtrList<Function2<scalar>> DTFuncs_;
88 
89  //- Generalised Fick's law diffusion coefficients field list.
90  // This is the diagonal of the mass diffusion coefficient matrix
92 
93  //- List of fields of the explicit part of the mass diffusion flux
94  // of the species
96 
97 
98  // Workspace for diffusion coefficient transformation
99 
100  //- Molecular weights of the species
101  scalarField W;
102 
103  //- List of mass-fraction field pointers
104  // for the internal mesh or a patch
106 
107  //- Matrix of mass diffusion coefficient field pointers
108  // for the internal mesh or a patch
110 
111  //- Mass-fractions at a cell or face
112  scalarField Y;
113 
114  //- Mole-fractions at a cell or face
115  scalarField X;
116 
117  //- Binary mass diffusion coefficient matrix at a cell or face
119 
120  //- Matrix form of the coefficients in the Maxwell-Stefan equation
121  // at a cell or face
122  LUscalarMatrix A;
123 
124  //- Matrix of composition coefficients at a cell or face
125  // used to transform the binary mass diffusion coefficients into
126  // the generalised Fick's law diffusion coefficients
128 
129  //- Inverse of A
130  scalarSquareMatrix invA;
131 
132  //- Matrix of the generalised Fick's law diffusion coefficients
133  // at a cell or face
135 
136 
137  // Private member functions
138 
139  //- Transform the Binary mass diffusion coefficient matrix DD into the
140  // matrix of the generalised Fick's law diffusion coefficients D
141  void transformDiffusionCoefficient();
142 
143  //- Calls transformDiffusionCoefficient() for each cell and patch face
144  // and maps the resulting D into DijPtrs
145  void transformDiffusionCoefficientFields();
146 
147  //- Calls transformDiffusionCoefficientFields() for the internal mesh
148  // and each patch and maps the resulting DijPtrs into Dii_ and Dij
149  void transform(List<PtrList<volScalarField>>& Dij);
150 
151 
152 public:
153 
154  typedef typename BasicThermophysicalTransportModel::alphaField
155  alphaField;
156 
157  typedef typename
160 
161  typedef typename BasicThermophysicalTransportModel::thermoModel
162  thermoModel;
163 
164 
165  // Constructors
166 
167  //- Construct from a momentum transport model and a thermo model
169  (
170  const word& type,
171  const momentumTransportModel& momentumTransport,
172  const thermoModel& thermo
173  );
174 
175 
176  //- Destructor
177  virtual ~MaxwellStefan()
178  {}
179 
180 
181  // Member Functions
182 
183  //- Read thermophysicalTransport dictionary
184  virtual bool read();
185 
186  //- Effective mass diffusion coefficient
187  // for a given specie mass-fraction [kg/m/s]
188  // This is the self-diffusion coefficient
189  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const;
190 
191  //- Effective mass diffusion coefficient
192  // for a given specie mass-fraction for patch [kg/m/s]
193  virtual tmp<scalarField> DEff
194  (
195  const volScalarField& Yi,
196  const label patchi
197  ) const;
198 
199  //- Return the heat flux [W/m^2]
200  virtual tmp<surfaceScalarField> q() const;
201 
202  //- Return the source term for the energy equation
203  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
204 
205  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
206  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
207 
208  //- Return the source term for the given specie mass-fraction equation
209  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
210 
211  //- Update the diffusion coefficients and flux corrections
212  virtual void correct();
213 };
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #ifdef NoRepository
223  #include "MaxwellStefan.C"
224 #endif
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual ~MaxwellStefan()
Destructor.
BasicThermophysicalTransportModel::thermoModel thermoModel
fluidReactionThermo & thermo
Definition: createFields.H:28
Class to perform the LU decomposition on a symmetric matrix.
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::momentumTransportModel momentumTransportModel
virtual void correct()
Update the diffusion coefficients and flux corrections.
BasicThermophysicalTransportModel::alphaField alphaField
compressibleMomentumTransportModel momentumTransportModel
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
MaxwellStefan(const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
Definition: MaxwellStefan.C:43
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: MaxwellStefan.C:94
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
thermo he()
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
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
A templated 2D square matrix of objects of <T>, where the n x n matrix dimension is known and used fo...
Definition: SquareMatrix.H:58
Namespace for OpenFOAM.