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 | Copyright (C) 2015 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 "runTimeSelectionTables.H"
46 
47 namespace Foam
48 {
49 
50 class phaseModel;
51 class phasePair;
52 
53 /*---------------------------------------------------------------------------*\
54  Class interfaceCompositionModel Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 {
59 protected:
60 
61  // Protected data
62 
63  //- Phase pair
64  const phasePair& pair_;
65 
66  //- Names of the transferring species
68 
69 
70 public:
71 
72  //- Runtime type information
73  TypeName("interfaceCompositionModel");
74 
75 
76  // Declare runtime construction
77 
79  (
80  autoPtr,
82  dictionary,
83  (
84  const dictionary& dict,
85  const phasePair& pair
86  ),
87  (dict, pair)
88  );
89 
90 
91  // Constructors
92 
93  //- Construct from a dictionary and a phase pair
95  (
96  const dictionary& dict,
97  const phasePair& pair
98  );
99 
100 
101  //- Destructor
102  virtual ~interfaceCompositionModel();
103 
104 
105  // Selectors
106 
108  (
109  const dictionary& dict,
110  const phasePair& pair
111  );
112 
113 
114  // Member Functions
115 
116  //- Update the composition
117  virtual void update(const volScalarField& Tf) = 0;
118 
119  //- Return the transferring species names
120  const hashedWordList& species() const;
121 
122  //- Returns whether the species is transported by the model and
123  // provides the name of the diffused species
124  bool transports
125  (
126  word& speciesName
127  ) const;
128 
129  //- Interface mass fraction
130  virtual tmp<volScalarField> Yf
131  (
132  const word& speciesName,
133  const volScalarField& Tf
134  ) const = 0;
135 
136  //- The interface mass fraction derivative w.r.t. temperature
138  (
139  const word& speciesName,
140  const volScalarField& Tf
141  ) const = 0;
142 
143  //- Mass fraction difference between the interface and the field
144  virtual tmp<volScalarField> dY
145  (
146  const word& speciesName,
147  const volScalarField& Tf
148  ) const = 0;
149 
150  //- Mass diffusivity
151  virtual tmp<volScalarField> D
152  (
153  const word& speciesName
154  ) const = 0;
155 
156  //- Latent heat
157  virtual tmp<volScalarField> L
158  (
159  const word& speciesName,
160  const volScalarField& Tf
161  ) const = 0;
162 
163  //- Add latent heat flow rate to total
164  virtual void addMDotL
165  (
166  const volScalarField& K,
167  const volScalarField& Tf,
168  volScalarField& mDotL,
169  volScalarField& mDotLPrime
170  ) const = 0;
171 };
172 
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 } // End namespace Foam
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 #endif
181 
182 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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.
bool transports(word &speciesName) const
Returns whether the species is transported by the model and.
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.
const hashedWordList & species() const
Return the transferring species names.
CGAL::Exact_predicates_exact_constructions_kernel K
virtual void update(const volScalarField &Tf)=0
Update the composition.
TypeName("interfaceCompositionModel")
Runtime type information.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat.
const phasePair & pair_
Phase pair.
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
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const =0
Add latent heat flow rate to total.
virtual ~interfaceCompositionModel()
Destructor.
const hashedWordList speciesNames_
Names of the transferring species.
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.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:54
virtual tmp< volScalarField > D(const word &speciesName) const =0
Mass diffusivity.
Namespace for OpenFOAM.