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) 2015-2018 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::interfaceCompositionModel
29
30Description
31 Generic base class for interface composition models. These models describe
32 the composition in phase 1 of the supplied pair at the interface with phase
33 2.
34
35SourceFiles
36 interfaceCompositionModel.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef interfaceCompositionModel_H
41#define interfaceCompositionModel_H
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45#include "volFields.H"
46#include "dictionary.H"
47#include "hashedWordList.H"
49
50namespace Foam
51{
52
53// Forward Declarations
54class phaseModel;
55class phasePair;
56
57/*---------------------------------------------------------------------------*\
58 Class interfaceCompositionModel Declaration
59\*---------------------------------------------------------------------------*/
62{
63protected:
64
65 // Protected Data
66
67 //- Phase pair
68 const phasePair& pair_;
69
70 //- Names of the transferring species
71 const hashedWordList speciesNames_;
72
73
74public:
75
76 //- Runtime type information
77 TypeName("interfaceCompositionModel");
78
79
80 // Declare runtime construction
83 (
84 autoPtr,
87 (
88 const dictionary& dict,
89 const phasePair& pair
90 ),
91 (dict, pair)
92 );
93
94
95 // Constructors
96
97 //- Construct from a dictionary and a phase pair
99 (
100 const dictionary& dict,
101 const phasePair& pair
102 );
103
104
105 //- Destructor
106 virtual ~interfaceCompositionModel() = default;
107
108
109 // Selectors
112 (
113 const dictionary& dict,
114 const phasePair& pair
115 );
116
117
118 // Member Functions
119
120 //- Update the composition
121 virtual void update(const volScalarField& Tf) = 0;
122
123 //- The transferring species names
124 const hashedWordList& species() const
125 {
126 return speciesNames_;
127 }
128
129 //- Returns whether the species is transported by the model and
130 //- provides the name of the diffused species
131 bool transports
132 (
133 word& speciesName
134 ) const;
135
136 //- Interface mass fraction
137 virtual tmp<volScalarField> Yf
138 (
139 const word& speciesName,
140 const volScalarField& Tf
141 ) const = 0;
142
143 //- The interface mass fraction derivative w.r.t. temperature
145 (
146 const word& speciesName,
147 const volScalarField& Tf
148 ) const = 0;
149
150 //- Mass fraction difference between the interface and the field
151 virtual tmp<volScalarField> dY
152 (
153 const word& speciesName,
154 const volScalarField& Tf
155 ) const = 0;
156
157 //- Mass diffusivity
158 virtual tmp<volScalarField> D
159 (
160 const word& speciesName
161 ) const = 0;
162
163 //- Latent heat
164 virtual tmp<volScalarField> L
165 (
166 const word& speciesName,
167 const volScalarField& Tf
168 ) const = 0;
169
170 //- Add latent heat flow rate to total
171 virtual void addMDotL
172 (
173 const volScalarField& K,
174 const volScalarField& Tf,
175 volScalarField& mDotL,
176 volScalarField& mDotLPrime
177 ) const = 0;
178};
179
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183} // End namespace Foam
184
185// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186
187#endif
188
189// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
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
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
virtual tmp< volScalarField > D(const word &speciesName) const =0
Mass diffusivity.
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat.
virtual void update(const volScalarField &Tf)=0
Update the composition.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
const hashedWordList speciesNames_
Names of the transferring species.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const =0
Mass fraction difference between the interface and the field.
virtual ~interfaceCompositionModel()=default
Destructor.
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
TypeName("interfaceCompositionModel")
Runtime type information.
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const =0
Add latent heat flow rate to total.
const hashedWordList & species() const
The transferring species names.
const phasePair & pair_
Phase pair.
const phasePair & pair() const
The phase pair.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
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.
virtual ~interfaceCompositionModel()=default
Destructor.
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
mesh update()
Namespace for OpenFOAM.
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
const dimensionedScalar & D
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
const vector L(dict.get< vector >("L"))