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-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::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"
46 #include "runTimeSelectionTables.H"
47 #include "sidedPhaseInterface.H"
48 #include "SidedInterfacialModel.H"
49 
50 namespace Foam
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class interfaceCompositionModel Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 {
59  // Private data
60 
61  //- Interface
62  const sidedPhaseInterface interface_;
63 
64  //- Names of the transferring species
65  const hashedWordList species_;
66 
67  //- Lewis number
68  const dimensionedScalar Le_;
69 
70  //- Multi-component thermo model for this side of the interface
71  const rhoMulticomponentThermo& thermo_;
72 
73  //- General thermo model for the other side of the interface
74  const rhoThermo& otherThermo_;
75 
76 
77 public:
78 
79  //- Runtime type information
80  TypeName("interfaceCompositionModel");
81 
82 
83  // Declare runtime construction
84 
86  (
87  autoPtr,
89  dictionary,
90  (
91  const dictionary& dict,
93  ),
94  (dict, interface)
95  );
96 
97 
98  // Constructors
99 
100  //- Construct from a dictionary and an interface
102  (
103  const dictionary& dict,
105  );
106 
107 
108  //- Destructor
109  virtual ~interfaceCompositionModel();
110 
111 
112  // Selectors
113 
115  (
116  const dictionary& dict,
117  const phaseInterface& interface,
118  const bool outer=true
119  );
120 
121 
122  // Member Functions
123 
124  // Access
125 
126  //- Return the interface
127  inline const sidedPhaseInterface& interface() const;
128 
129  //- Return the transferring species names
130  inline const hashedWordList& species() const;
131 
132  //- Return the thermo
133  inline const rhoMulticomponentThermo& thermo() const;
134 
135  //- Return the composition
136  inline const basicSpecieMixture& composition() const;
137 
138  //- Return the other thermo
139  inline const rhoThermo& otherThermo() const;
140 
141  //- Return whether the other side has a multi-specie composition
142  inline bool otherHasComposition() const;
143 
144  //- Return the other composition
145  inline const basicSpecieMixture& otherComposition() const;
146 
147 
148  // Evaluation
149 
150  //- Interface mass fraction
151  virtual tmp<volScalarField> Yf
152  (
153  const word& speciesName,
154  const volScalarField& Tf
155  ) const = 0;
156 
157  //- The interface mass fraction derivative w.r.t. temperature
159  (
160  const word& speciesName,
161  const volScalarField& Tf
162  ) const = 0;
163 
164  //- Mass fraction difference between the interface and the field
166  (
167  const word& speciesName,
168  const volScalarField& Tf
169  ) const;
170 
171  //- Mass fraction difference between the interface and the field
172  // derivative w.r.t. temperature
174  (
175  const word& speciesName,
176  const volScalarField& Tf
177  ) const;
178 
179  //- Mass diffusivity
181  (
182  const word& speciesName
183  ) const;
184 
185 
186  //- Update the composition
187  virtual void update(const volScalarField& Tf) = 0;
188 };
189 
190 
191 /*---------------------------------------------------------------------------*\
192  Class sidedInterfaceCompositionModel Declaration
193 \*---------------------------------------------------------------------------*/
194 
196 :
197  public SidedInterfacialModel<interfaceCompositionModel>
198 {
199 public:
200 
201  // Constructors
202 
203  //- Inherit base class constructors
204  using
207 
208 
209  // Selectors
210 
212  (
213  const dictionary& dict,
215  );
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace Foam
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #endif
230 
231 // ************************************************************************* //
Generic GeometricField class.
const phaseInterface & interface() const
Access the interface.
SidedInterfacialModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A wordList with hashed indices for faster lookup by name.
Generic base class for interface composition models. These models describe the composition in phase 1...
interfaceCompositionModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
const rhoMulticomponentThermo & thermo() const
Return the thermo.
tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phaseInterface &interface),(dict, interface))
const basicSpecieMixture & otherComposition() const
Return the other composition.
const sidedPhaseInterface & interface() const
Return the interface.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
virtual void update(const volScalarField &Tf)=0
Update the composition.
tmp< volScalarField > dYfPrime(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
bool otherHasComposition() const
Return whether the other side has a multi-specie composition.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phaseInterface &interface, const bool outer=true)
const basicSpecieMixture & composition() const
Return the composition.
TypeName("interfaceCompositionModel")
Runtime type information.
const rhoThermo & otherThermo() const
Return the other thermo.
const hashedWordList & species() const
Return the transferring species names.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Base-class for multi-component fluid thermodynamic properties based on density.
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:55
static autoPtr< sidedInterfaceCompositionModel > New(const dictionary &dict, const phaseInterface &interface)
Class to represent a certain side of an interface between phases.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
void outer(FieldField< Field1, typename outerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Macros to ease declaration of run-time selection tables.
dictionary dict