All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
interfaceCompositionModel.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) 2015-2020 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::interfaceCompositionModel
26 
27 Description
28  Generic base class for interface composition models. These models describe
29  the composition in phase 1 of the supplied pair at the interface with phase
30  2.
31 
32 SourceFiles
33  interfaceCompositionModel.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef interfaceCompositionModel_H
38 #define interfaceCompositionModel_H
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 #include "volFields.H"
43 #include "dictionary.H"
44 #include "hashedWordList.H"
45 #include "rhoReactionThermo.H"
46 #include "runTimeSelectionTables.H"
47 
48 namespace Foam
49 {
50 
51 class phaseModel;
52 class phasePair;
53 
54 /*---------------------------------------------------------------------------*\
55  Class interfaceCompositionModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 {
60  // Private data
61 
62  //- Phase pair
63  const phasePair& pair_;
64 
65  //- Names of the transferring species
66  const hashedWordList species_;
67 
68  //- Lewis number
69  const dimensionedScalar Le_;
70 
71  //- Phase
72  const phaseModel& phase_;
73 
74  //- Other phase
75  const phaseModel& otherPhase_;
76 
77  //- Multi-component thermo model for this side of the interface
78  const rhoReactionThermo& thermo_;
79 
80  //- General thermo model for the other side of the interface
81  const rhoThermo& otherThermo_;
82 
83 
84 public:
85 
86  //- Runtime type information
87  TypeName("interfaceCompositionModel");
88 
89 
90  // Declare runtime construction
91 
93  (
94  autoPtr,
96  dictionary,
97  (
98  const dictionary& dict,
99  const phasePair& pair
100  ),
101  (dict, pair)
102  );
103 
104 
105  // Constructors
106 
107  //- Construct from a dictionary and a phase pair
109  (
110  const dictionary& dict,
111  const phasePair& pair
112  );
113 
114 
115  //- Destructor
116  virtual ~interfaceCompositionModel();
117 
118 
119  // Selectors
120 
122  (
123  const dictionary& dict,
124  const phasePair& pair
125  );
126 
127 
128  // Member Functions
129 
130  // Access
131 
132  //- Return the phase pair
133  inline const phasePair& pair() const;
134 
135  //- Return the transferring species names
136  inline const hashedWordList& species() const;
137 
138  //- Return the thermo
139  inline const rhoReactionThermo& thermo() const;
140 
141  //- Return the composition
142  inline const basicSpecieMixture& composition() const;
143 
144  //- Return the other thermo
145  inline const rhoThermo& otherThermo() const;
146 
147  //- Return whether the other side has a multi-specie composition
148  inline bool otherHasComposition() const;
149 
150  //- Return the other composition
151  inline const basicSpecieMixture& otherComposition() const;
152 
153 
154  // Evaluation
155 
156  //- Interface mass fraction
157  virtual tmp<volScalarField> Yf
158  (
159  const word& speciesName,
160  const volScalarField& Tf
161  ) const = 0;
162 
163  //- The interface mass fraction derivative w.r.t. temperature
165  (
166  const word& speciesName,
167  const volScalarField& Tf
168  ) const = 0;
169 
170  //- Mass fraction difference between the interface and the field
172  (
173  const word& speciesName,
174  const volScalarField& Tf
175  ) const;
176 
177  //- Mass fraction difference between the interface and the field
178  // derivative w.r.t. temperature
180  (
181  const word& speciesName,
182  const volScalarField& Tf
183  ) const;
184 
185  //- Mass diffusivity
187  (
188  const word& speciesName
189  ) const;
190 
191 
192  //- Update the composition
193  virtual void update(const volScalarField& Tf) = 0;
194 };
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace Foam
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 #endif
208 
209 // ************************************************************************* //
const rhoReactionThermo & thermo() const
Return the thermo.
dictionary dict
const hashedWordList & species() const
Return the transferring species names.
tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
tmp< volScalarField > dYfPrime(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
const basicSpecieMixture & otherComposition() const
Return the other composition.
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
const phasePair & pair() const
Return the phase pair.
virtual void update(const volScalarField &Tf)=0
Update the composition.
TypeName("interfaceCompositionModel")
Runtime type information.
Generic base class for interface composition models. These models describe the composition in phase 1...
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
A class for handling words, derived from string.
Definition: word.H:59
Base-class for multi-component fluid thermodynamic properties based on density.
virtual ~interfaceCompositionModel()
Destructor.
bool otherHasComposition() const
Return whether the other side has a multi-specie composition.
const rhoThermo & otherThermo() const
Return the other thermo.
const basicSpecieMixture & composition() const
Return the composition.
A wordList with hashed indices for faster lookup by name.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:52
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
Namespace for OpenFOAM.