interfaceCompositionModel.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::interfaceCompositionModel
28
29Description
30 Generic base class for interface models. Mass transfer models are
31 interface models between two thermos.
32 Abstract class for mass transfer functions
33
34SourceFiles
35 interfaceCompositionModel.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef interfaceCompositionModel_H
40#define interfaceCompositionModel_H
41
42#include "volFields.H"
43#include "dictionary.H"
44#include "Enum.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50namespace Foam
51{
52
53// Forward Declarations
54class phasePair;
55
56namespace multiphaseInter
57{
58
59// Forward Declarations
60class phaseModel;
61
62/*---------------------------------------------------------------------------*\
63 Class interfaceCompositionModel Declaration
64\*---------------------------------------------------------------------------*/
67{
68public:
69
70 // Public Types
71
72 //- Enumeration for variable based mass transfer models
73 enum modelVariable
74 {
75 T, /* temperature based */
76 P, /* pressure based */
77 Y, /* mass fraction based */
78 alpha /* alpha source */
79 };
80
81protected:
82
83 // Protected Data
84
85 //- Selection names for the modelVariable
87
88 //- Enumeration for the model variable
90
91 //- Add volume change in pEq
93
94 //- Phase pair
95 const phasePair& pair_;
96
97 //- Names of the transferring specie
99
100 //- Reference to mesh
101 const fvMesh& mesh_;
102
103
104public:
105
106 //- Runtime type information
107 TypeName("interfaceCompositionModel");
108
109
110 // Declare runtime construction
113 (
114 autoPtr,
117 (
118 const dictionary& dict,
119 const phasePair& pair
120 ),
121 (dict, pair)
122 );
123
124
125 // Constructors
126
127 //- Construct from a dictionary and a phase pair
129 (
130 const dictionary& dict,
131 const phasePair& pair
132 );
133
134
135 //- Destructor
136 virtual ~interfaceCompositionModel() = default;
137
138
139 // Selectors
140
142 (
143 const dictionary& dict,
144 const phasePair& pair
145 );
146
147
148 // Member Functions
149
150 //- Return the transferring species name
151 const word transferSpecie() const;
152
153 //- The phase pair
154 const phasePair& pair() const;
155
156 //- Return the multiphaseInterSystem this interface belongs to
157 const multiphaseInterSystem& fluid() const;
158
159 //- Interface mass fraction
160 virtual tmp<volScalarField> Yf
161 (
162 const word& speciesName,
163 const volScalarField& Tf
164 ) const = 0;
165
166 //- Mass fraction difference between the interface and the field
167 virtual tmp<volScalarField> dY
168 (
169 const word& speciesName,
170 const volScalarField& Tf
171 ) const = 0;
172
173 //- Specie mass diffusivity for pure mixture
175 (
176 const word& speciesName
177 ) const = 0;
178
179 //- Specie mass diffusivity for specie in a multicomponent
181 (
182 const word& speciesName
183 ) const = 0;
184
185 //- Latent heat (delta Hc)
186 virtual tmp<volScalarField> L
187 (
188 const word& speciesName,
189 const volScalarField& Tf
190 ) const = 0;
191
192 //- Explicit full mass transfer
194 (
195 const volScalarField& field
196 ) = 0;
197
198 //- Implicit mass transfer
200 (
201 label modelVariable,
202 const volScalarField& field
203 ) = 0;
204
205 //- Explicit mass transfer
207 (
208 label modelVariable,
209 const volScalarField& field
210 ) = 0;
211
212 //- Reference value
213 virtual const dimensionedScalar& Tactivate() const noexcept = 0;
214
215 //- Add/subtract alpha*div(U) as a source term
216 //- for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2)
217 virtual bool includeDivU() const noexcept;
218
219 //- Add volume change in pEq
220 bool includeVolChange();
221
222 //- Returns the variable on which the model is based
223 const word& variable() const;
224};
225
226
227// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228
229} // End namespace Foam
230} // End multiphaseInter
231
232// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233
234#endif
235
236// ************************************************************************* //
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
static const Enum< modelVariable > modelVariableNames_
Selection names for the modelVariable.
virtual const dimensionedScalar & Tactivate() const noexcept=0
Reference value.
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)=0
Implicit mass transfer.
const word transferSpecie() const
Return the transferring species name.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat (delta Hc)
const phasePair & pair() const
The phase pair.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
virtual tmp< volScalarField > Dfrom(const word &speciesName) const =0
Specie mass diffusivity for pure mixture.
const multiphaseInterSystem & fluid() const
Return the multiphaseInterSystem this interface belongs to.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const =0
Mass fraction difference between the interface and the field.
modelVariable
Enumeration for variable based mass transfer models.
virtual ~interfaceCompositionModel()=default
Destructor.
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
const word & variable() const
Returns the variable on which the model is based.
TypeName("interfaceCompositionModel")
Runtime type information.
word speciesName_
Names of the transferring specie.
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)=0
Explicit mass transfer.
virtual tmp< volScalarField > Dto(const word &speciesName) const =0
Specie mass diffusivity for specie in a multicomponent.
virtual tmp< volScalarField > Kexp(const volScalarField &field)=0
Explicit full mass transfer.
modelVariable modelVariable_
Enumeration for the model variable.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
Namespace for OpenFOAM.
const direction noexcept
Definition: Scalar.H:223
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
const vector L(dict.get< vector >("L"))