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-2023 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  public UpdateableMeshObject<fvMesh>
77 {
78  // Private data
79 
80  // Model coefficients
81 
82  //- Array of specie binary mass diffusion coefficient functions
83  // [m^2/s]
85 
86  //- List of specie Soret thermal diffusion coefficient
87  // functions [kg/m/s]
88  PtrList<Function2<scalar>> DTFuncs_;
89 
90  //- Generalised Fick's law diffusion coefficients field list.
91  // This is the diagonal of the mass diffusion coefficient matrix
92  mutable PtrList<volScalarField> Dii_;
93 
94  //- List of fields of the explicit part of the mass diffusion flux
95  // of the species
96  mutable PtrList<surfaceScalarField> jexp_;
97 
98 
99  // Workspace for diffusion coefficient transformation
100 
101  //- Molecular weights of the species
102  mutable scalarField W;
103 
104  //- List of mass-fraction field pointers
105  // for the internal mesh or a patch
106  mutable List<const scalarField*> YPtrs;
107 
108  //- Matrix of mass diffusion coefficient field pointers
109  // for the internal mesh or a patch
110  mutable SquareMatrix<scalarField*> DijPtrs;
111 
112  //- Mass-fractions at a cell or face
113  mutable scalarField Y;
114 
115  //- Mole-fractions at a cell or face
116  mutable scalarField X;
117 
118  //- Binary mass diffusion coefficient matrix at a cell or face
119  mutable scalarSquareMatrix DD;
120 
121  //- Matrix form of the coefficients in the Maxwell-Stefan equation
122  // at a cell or face
123  mutable LUscalarMatrix A;
124 
125  //- Matrix of composition coefficients at a cell or face
126  // used to transform the binary mass diffusion coefficients into
127  // the generalised Fick's law diffusion coefficients
128  mutable scalarSquareMatrix B;
129 
130  //- Inverse of A
131  mutable scalarSquareMatrix invA;
132 
133  //- Matrix of the generalised Fick's law diffusion coefficients
134  // at a cell or face
135  mutable scalarSquareMatrix D;
136 
137 
138  // Private member functions
139 
140  //- Transform the Binary mass diffusion coefficient matrix DD into the
141  // matrix of the generalised Fick's law diffusion coefficients D
142  void transformDiffusionCoefficient() const;
143 
144  //- Calls transformDiffusionCoefficient() for each cell and patch face
145  // and maps the resulting D into DijPtrs
146  void transformDiffusionCoefficientFields() const;
147 
148  //- Calls transformDiffusionCoefficientFields() for the internal mesh
149  // and each patch and maps the resulting DijPtrs into Dii_ and Dij
150  void transform(List<PtrList<volScalarField>>& Dij) const;
151 
152  //- Update Dii_
153  void updateDii() const;
154 
155  //- Update Dii_ if not yet constructed and return
156  const PtrList<volScalarField>& Dii() const;
157 
158  //- Update jexp_ if not yet constructed and return
159  const PtrList<surfaceScalarField>& jexp() const;
160 
161 
162 public:
163 
164  typedef typename BasicThermophysicalTransportModel::alphaField
165  alphaField;
166 
167  typedef typename
170 
171  typedef typename BasicThermophysicalTransportModel::thermoModel
172  thermoModel;
173 
174 
175  // Constructors
176 
177  //- Construct from a momentum transport model and a thermo model
179  (
180  const word& type,
181  const momentumTransportModel& momentumTransport,
182  const thermoModel& thermo
183  );
184 
185 
186  //- Destructor
187  virtual ~MaxwellStefan()
188  {}
189 
190 
191  // Member Functions
192 
193  //- Read thermophysicalTransport dictionary
194  virtual bool read();
195 
196  //- Effective mass diffusion coefficient
197  // for a given specie mass-fraction [kg/m/s]
198  // This is the self-diffusion coefficient
199  virtual tmp<volScalarField> DEff(const volScalarField& Yi) const;
200 
201  //- Effective mass diffusion coefficient
202  // for a given specie mass-fraction for patch [kg/m/s]
203  virtual tmp<scalarField> DEff
204  (
205  const volScalarField& Yi,
206  const label patchi
207  ) const;
208 
209  //- Return the heat flux [W/m^2]
210  virtual tmp<surfaceScalarField> q() const;
211 
212  //- Return the source term for the energy equation
213  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const;
214 
215  //- Return the specie flux for the given specie mass-fraction [kg/m^2/s]
216  virtual tmp<surfaceScalarField> j(const volScalarField& Yi) const;
217 
218  //- Return the source term for the given specie mass-fraction equation
219  virtual tmp<fvScalarMatrix> divj(volScalarField& Yi) const;
220 
221  //- Update the diffusion coefficients and flux corrections
222  virtual void predict();
223 
224 
225  // Mesh changes
226 
227  //- Update for mesh motion
228  virtual bool movePoints();
229 
230  //- Update topology using the given map
231  virtual void topoChange(const polyTopoChangeMap& map);
232 
233  //- Update from another mesh using the given map
234  virtual void mapMesh(const polyMeshMap& map);
235 
236  //- Redistribute or update using the given distribution map
237  virtual void distribute(const polyDistributionMap& map);
238 };
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 } // End namespace Foam
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 #ifdef NoRepository
248  #include "MaxwellStefan.C"
249 #endif
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
Generic GeometricField class.
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: List.H:91
virtual bool movePoints()
Update for mesh motion.
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
BasicThermophysicalTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
virtual void predict()
Update the diffusion coefficients and flux corrections.
virtual void mapMesh(const polyMeshMap &map)
Update from another mesh using the given map.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
virtual ~MaxwellStefan()
Destructor.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
MaxwellStefan(const word &type, const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
BasicThermophysicalTransportModel::thermoModel thermoModel
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
virtual bool read()
Read thermophysicalTransport dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A templated 2D square matrix of objects of <T>, where the n x n matrix dimension is known and used fo...
Definition: SquareMatrix.H:61
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
compressibleMomentumTransportModel momentumTransportModel
Namespace for OpenFOAM.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()
fluidMulticomponentThermo & thermo
Definition: createFields.H:31