phaseSystem.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::phaseSystem
26 
27 Description
28  Class to represent a system of phases and model interfacial transfers
29  between them.
30 
31 SourceFiles
32  phaseSystem.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef phaseSystem_H
37 #define phaseSystem_H
38 
39 #include "IOdictionary.H"
40 
41 #include "phaseModel.H"
42 #include "phasePair.H"
43 #include "orderedPhasePair.H"
44 #include "HashPtrTable.H"
45 #include "PtrListDictionary.H"
46 
47 #include "IOMRFZoneList.H"
48 #include "fvIOoptionList.H"
49 
50 #include "volFields.H"
51 #include "surfaceFields.H"
52 #include "fvMatricesFwd.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 class blendingMethod;
60 template <class modelType> class BlendedInterfacialModel;
61 class surfaceTensionModel;
62 class aspectRatioModel;
63 
64 /*---------------------------------------------------------------------------*\
65  Class phaseSystem Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class phaseSystem
69 :
70  public IOdictionary
71 {
72 public:
73 
74  // Public typedefs
75 
76  typedef
78  <
82  >
83  KdTable;
84 
85  typedef
87  <
89  phasePairKey,
91  >
92  VmTable;
93 
94  typedef
96  <
98  word,
100  >
102 
103  typedef
105  <
107  word,
109  >
111 
112  typedef
114  <
116  word,
118  >
122 
123 
124 protected:
125 
126  // Protected typedefs
127 
128  typedef
130  dictTable;
131 
132  typedef
135 
136  typedef
139 
140  typedef
141  HashTable
142  <
144  phasePairKey,
146  >
148 
149  typedef
150  HashTable
151  <
153  phasePairKey,
155  >
157 
158 
159  // Protected data
160 
161  //- Reference to the mesh
162  const fvMesh& mesh_;
163 
164  //- Phase models
165  phaseModelList phaseModels_;
166 
167  //- Phase pairs
168  phasePairTable phasePairs_;
169 
170  //- Total volumetric flux
172 
173  //- Rate of change of pressure
174  volScalarField dpdt_;
175 
176  //- Optional MRF zones
178 
179  //- Optional FV-options
181 
182  //- Blending methods
183  blendingMethodTable blendingMethods_;
184 
185 
186  // Sub Models
187 
188  //- Surface tension models
190 
191  //- Aspect ratio models
193 
194 
195  // Protected member functions
196 
197  //- Calculate and return the mixture flux
199  (
200  const phaseModelList& phaseModels
201  ) const;
202 
203  //- Generate pairs
204  void generatePairs
205  (
206  const dictTable& modelDicts
207  );
208 
209  //- Generate pairs and sub-model tables
210  template<class modelType>
211  void createSubModels
212  (
213  const dictTable& modelDicts,
214  HashTable
215  <
217  phasePairKey,
219  >& models
220  );
221 
222  //- Generate pairs and sub-model tables
223  template<class modelType>
225  (
226  const word& modelName,
227  HashTable
228  <
230  phasePairKey,
232  >& models
233  );
234 
235  //- Generate pairs and blended sub-model tables
236  template<class modelType>
238  (
239  const word& modelName,
240  HashTable
241  <
243  phasePairKey,
245  >& models
246  );
247 
248  //- Generate pairs and per-phase sub-model tables
249  template<class modelType>
251  (
252  const word& modelName,
253  HashTable
254  <
256  phasePairKey,
258  >& models
259  );
260 
261 
262 public:
263 
264  //- Runtime type information
265  TypeName("phaseSystem");
266 
267  //- Default name of the phase properties dictionary
268  static const word propertiesName;
269 
270 
271  // Constructors
272 
273  //- Construct from fvMesh
274  phaseSystem(const fvMesh& mesh);
275 
276 
277  //- Destructor
278  virtual ~phaseSystem();
279 
280 
281  // Member Functions
282 
283  //- Constant access the mesh
284  inline const fvMesh& mesh() const;
285 
286  //- Constant access the phase models
287  inline const phaseModelList& phases() const;
288 
289  //- Access the phase models
290  inline phaseModelList& phases();
291 
292  //- Constant access the phase pairs
293  inline const phasePairTable& phasePairs() const;
294 
295  //- Return the mixture density
296  tmp<volScalarField> rho() const;
297 
298  //- Return the mixture velocity
299  tmp<volVectorField> U() const;
300 
301  //- Constant access the mixture flux
302  inline const surfaceScalarField& phi() const;
303 
304  //- Access the mixture flux
305  inline surfaceScalarField& phi();
306 
307  //- Constant access the rate of change of the pressure
308  inline const volScalarField& dpdt() const;
309 
310  //- Access the rate of change of the pressure
311  inline volScalarField& dpdt();
312 
313  //- Return the surface tension coefficient
314  tmp<volScalarField> sigma(const phasePairKey& key) const;
315 
316  //- Return MRF zones
317  inline const IOMRFZoneList& MRF() const;
318 
319  //- Optional FV-options
320  inline fv::IOoptionList& fvOptions() const;
321 
322  //- Access a sub model between a phase pair
323  template <class modelType>
324  const modelType& lookupSubModel(const phasePair& key) const;
325 
326  //- Access a sub model between two phases
327  template <class modelType>
328  const modelType& lookupSubModel
329  (
330  const phaseModel& dispersed,
331  const phaseModel& continuous
332  ) const;
333 
334  //- Solve for the phase fractions
335  virtual void solve();
336 
337  //- Correct the fluid properties other than the thermo and turbulence
338  virtual void correct();
339 
340  //- Correct the kinematics
341  virtual void correctKinematics();
342 
343  //- Correct the thermodynamics
344  virtual void correctThermo();
345 
346  //- Correct the turbulence
347  virtual void correctTurbulence();
348 
349  //- Correct the energy transport e.g. alphat
350  virtual void correctEnergyTransport();
351 
352  //- Read base phaseProperties dictionary
353  virtual bool read();
354 };
355 
356 
357 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358 
359 } // End namespace Foam
360 
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362 
363 #include "phaseSystemI.H"
364 
365 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
366 
367 #ifdef NoRepository
368 # include "phaseSystemTemplates.C"
369 #endif
370 
371 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
372 
373 #endif
374 
375 // ************************************************************************* //
const phasePairTable & phasePairs() const
Constant access the phase pairs.
Definition: phaseSystemI.H:49
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition: phaseSystem.H:129
virtual void solve()
Solve for the phase fractions.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::surfaceFields.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
Definition: phaseSystem.H:82
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient.
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:133
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
void generatePairs(const dictTable &modelDicts)
Generate pairs.
const surfaceScalarField & phi() const
Constant access the mixture flux.
Definition: phaseSystemI.H:55
An STL-conforming hash table.
Definition: HashTable.H:61
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:188
Hashing function class, shared by all the derived classes.
Definition: string.H:90
IOMRFZoneList MRF_
Optional MRF zones.
Definition: phaseSystem.H:176
A class for handling words, derived from string.
Definition: word.H:59
virtual void correctKinematics()
Correct the kinematics.
HashTable< autoPtr< aspectRatioModel >, phasePairKey, phasePairKey::hash > aspectRatioModelTable
Definition: phaseSystem.H:155
volScalarField dpdt_
Rate of change of pressure.
Definition: phaseSystem.H:173
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:167
tmp< volVectorField > U() const
Return the mixture velocity.
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:79
Namespace for OpenFOAM.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
TypeName("phaseSystem")
Runtime type information.
aspectRatioModelTable aspectRatioModels_
Aspect ratio models.
Definition: phaseSystem.H:191
virtual void correctTurbulence()
Correct the turbulence.
fv::IOoptionList & fvOptions() const
Optional FV-options.
Definition: phaseSystemI.H:85
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:182
const phaseModelList & phases() const
Constant access the phase models.
Definition: phaseSystemI.H:35
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashTable< autoPtr< blendingMethod >, word, word::hash > blendingMethodTable
Definition: phaseSystem.H:137
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:267
HashPtrTable< fvScalarMatrix, word, string::hash > heatTransferTable
Definition: phaseSystem.H:109
const volScalarField & dpdt() const
Constant access the rate of change of the pressure.
Definition: phaseSystemI.H:67
HashPtrTable< fvScalarMatrix, word, string::hash > massTransferTable
Definition: phaseSystem.H:118
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
PtrListDictionary< phaseModel > phaseModelList
Definition: phaseSystem.H:120
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
Definition: phaseSystem.H:91
Forward declarations of fvMatrix specializations.
fv::IOoptionList fvOptions_
Optional FV-options.
Definition: phaseSystem.H:179
HashPtrTable< fvVectorMatrix, word, string::hash > momentumTransferTable
Definition: phaseSystem.H:100
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
tmp< volScalarField > rho() const
Return the mixture density.
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:164
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:170
virtual bool read()
Read base phaseProperties dictionary.
virtual ~phaseSystem()
Destructor.
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
virtual void correctThermo()
Correct the thermodynamics.
HashTable< autoPtr< surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
Definition: phaseSystem.H:146
virtual void correct()
Correct the fluid properties other than the thermo and turbulence.
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:161
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries...
Definition: IOMRFZoneList.H:66
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
const fvMesh & mesh() const
Constant access the mesh.
Definition: phaseSystemI.H:28
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
A class for managing temporary objects.
Definition: PtrList.H:118
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50