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 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::interfaceCompositionModel
28 
29 Description
30  Generic base class for interface composition models. These models describe
31  the composition in phase 1 of the supplied pair at the interface with phase
32  2.
33 
34 SourceFiles
35  interfaceCompositionModel.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef interfaceCompositionModel_H
40 #define interfaceCompositionModel_H
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 #include "volFields.H"
45 #include "dictionary.H"
46 #include "hashedWordList.H"
47 #include "runTimeSelectionTables.H"
48 
49 namespace Foam
50 {
51 
52 class phaseModel;
53 class phasePair;
54 
55 /*---------------------------------------------------------------------------*\
56  Class interfaceCompositionModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 {
61 protected:
62 
63  // Protected data
64 
65  //- Phase pair
66  const phasePair& pair_;
67 
68  //- Names of the transferring species
70 
71 
72 public:
73 
74  //- Runtime type information
75  TypeName("interfaceCompositionModel");
76 
77 
78  // Declare runtime construction
79 
81  (
82  autoPtr,
84  dictionary,
85  (
86  const dictionary& dict,
87  const phasePair& pair
88  ),
89  (dict, pair)
90  );
91 
92 
93  // Constructors
94 
95  //- Construct from a dictionary and a phase pair
97  (
98  const dictionary& dict,
99  const phasePair& pair
100  );
101 
102 
103  //- Destructor
104  virtual ~interfaceCompositionModel() = default;
105 
106 
107  // Selectors
108 
110  (
111  const dictionary& dict,
112  const phasePair& pair
113  );
114 
115 
116  // Member Functions
117 
118  //- Update the composition
119  virtual void update(const volScalarField& Tf) = 0;
120 
121  //- Return the transferring species names
122  const hashedWordList& species() const;
123 
124  //- Returns whether the species is transported by the model and
125  //- provides the name of the diffused species
126  bool transports
127  (
128  word& speciesName
129  ) const;
130 
131  //- Interface mass fraction
132  virtual tmp<volScalarField> Yf
133  (
134  const word& speciesName,
135  const volScalarField& Tf
136  ) const = 0;
137 
138  //- The interface mass fraction derivative w.r.t. temperature
140  (
141  const word& speciesName,
142  const volScalarField& Tf
143  ) const = 0;
144 
145  //- Mass fraction difference between the interface and the field
146  virtual tmp<volScalarField> dY
147  (
148  const word& speciesName,
149  const volScalarField& Tf
150  ) const = 0;
151 
152  //- Mass diffusivity
153  virtual tmp<volScalarField> D
154  (
155  const word& speciesName
156  ) const = 0;
157 
158  //- Latent heat
159  virtual tmp<volScalarField> L
160  (
161  const word& speciesName,
162  const volScalarField& Tf
163  ) const = 0;
164 
165  //- Add latent heat flow rate to total
166  virtual void addMDotL
167  (
168  const volScalarField& K,
169  const volScalarField& Tf,
170  volScalarField& mDotL,
171  volScalarField& mDotLPrime
172  ) const = 0;
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace Foam
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 #endif
183 
184 // ************************************************************************* //
volFields.H
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:51
Foam::interfaceCompositionModel::New
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
Definition: newInterfaceCompositionModel.C:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::interfaceCompositionModel::~interfaceCompositionModel
virtual ~interfaceCompositionModel()=default
Destructor.
Foam::tmp< volScalarField >
Foam::interfaceCompositionModel::D
virtual tmp< volScalarField > D(const word &speciesName) const =0
Mass diffusivity.
Foam::interfaceCompositionModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
Foam::interfaceCompositionModel::dY
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const =0
Mass fraction difference between the interface and the field.
Foam::interfaceCompositionModel
Generic base class for interface composition models. These models describe the composition in phase 1...
Definition: interfaceCompositionModel.H:58
Foam::interfaceCompositionModel::pair
const phasePair & pair() const
Return pair.
Foam::interfaceCompositionModel::species
const hashedWordList & species() const
Return the transferring species names.
Definition: interfaceCompositionModel.C:56
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::interfaceCompositionModel::L
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat.
Foam::interfaceCompositionModel::TypeName
TypeName("interfaceCompositionModel")
Runtime type information.
Foam::hashedWordList
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
Definition: hashedWordList.H:54
Foam::interfaceCompositionModel::YfPrime
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
Foam::interfaceCompositionModel::pair_
const phasePair & pair_
Phase pair.
Definition: interfaceCompositionModel.H:65
Foam::interfaceCompositionModel::update
virtual void update(const volScalarField &Tf)=0
Update the composition.
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::interfaceCompositionModel::transports
bool transports(word &speciesName) const
Definition: interfaceCompositionModel.C:62
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
hashedWordList.H
Foam::interfaceCompositionModel::speciesNames_
const hashedWordList speciesNames_
Names of the transferring species.
Definition: interfaceCompositionModel.H:68
Foam::interfaceCompositionModel::interfaceCompositionModel
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
Definition: interfaceCompositionModel.C:44
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::interfaceCompositionModel::Yf
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
Foam::interfaceCompositionModel::addMDotL
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const =0
Add latent heat flow rate to total.
dictionary.H
Foam::GeometricField< scalar, fvPatchField, volMesh >