HeatAndMassTransferPhaseSystem.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-2016 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::HeatAndMassTransferPhaseSystem
26 
27 Description
28  Base class to support interfacial heat and mass transfer between a number
29  of phases.
30 
31 SourceFiles
32  HeatAndMassTransferPhaseSystem.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef HeatAndMassTransferPhaseSystem_H
37 #define HeatAndMassTransferPhaseSystem_H
38 
39 #include "phaseSystem.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 template<class modelType>
47 class BlendedInterfacialModel;
48 
49 class blendingMethod;
50 class heatTransferModel;
51 class massTransferModel;
52 
53 /*---------------------------------------------------------------------------*\
54  Class HeatAndMassTransferPhaseSystem Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 template<class BasePhaseSystem>
59 :
60  public BasePhaseSystem
61 {
62 protected:
63 
64  // Protected typedefs
65 
66  typedef HashTable
67  <
68  HashTable
69  <
71  >,
75 
76  typedef HashTable
77  <
78  HashTable
79  <
81  >,
82  phasePairKey,
85 
86 
87  // Protected data
88 
89  //- Mass transfer rate
91  dmdt_;
92 
93  //- Explicit part of the mass transfer rate
96 
97  //- Interface temperatures
99 
100  // Sub Models
101 
102  //- Heat transfer models
103  heatTransferModelTable heatTransferModels_;
104 
105  //- Mass transfer models
106  massTransferModelTable massTransferModels_;
107 
108 
109 public:
110 
111  // Constructors
112 
113  //- Construct from fvMesh
115 
116 
117  //- Destructor
119 
120 
121  // Member Functions
122 
123  //- Return true if there is mass transfer for phase
124  virtual bool transfersMass(const phaseModel& phase) const;
125 
126  //- Return the interfacial mass flow rate
127  virtual tmp<volScalarField> dmdt(const phasePairKey& key) const;
128 
129  //- Return the total interfacial mass transfer rate for phase
130  virtual tmp<volScalarField> dmdt(const phaseModel& phase) const;
131 
132  //- Return the momentum transfer matrices
134  momentumTransfer() const;
135 
136  //- Return the heat transfer matrices
138  heatTransfer() const;
139 
140  //- Return the mass transfer matrices
142  massTransfer() const = 0;
143 
144  //- Correct the thermodynamics
145  virtual void correctThermo() = 0;
146 
147  //- Read base phaseProperties dictionary
148  virtual bool read();
149 };
150 
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 } // End namespace Foam
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #ifdef NoRepository
160 #endif
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 #endif
165 
166 // ************************************************************************* //
Base class to support interfacial heat and mass transfer between a number of phases.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
HashTable< HashTable< autoPtr< BlendedInterfacialModel< massTransferModel > > >, phasePairKey, phasePairKey::hash > massTransferModelTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > Tf_
Interface temperatures.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdt_
Mass transfer rate.
virtual bool transfersMass(const phaseModel &phase) const
Return true if there is mass transfer for phase.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtExplicit_
Explicit part of the mass transfer rate.
massTransferModelTable massTransferModels_
Mass transfer models.
HashTable< HashTable< autoPtr< BlendedInterfacialModel< heatTransferModel > > >, phasePairKey, phasePairKey::hash > heatTransferModelTable
virtual ~HeatAndMassTransferPhaseSystem()
Destructor.
An STL-conforming hash table.
Definition: HashTable.H:62
virtual void correctThermo()=0
Correct the thermodynamics.
HeatAndMassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const =0
Return the mass transfer matrices.
heatTransferModelTable heatTransferModels_
Heat transfer models.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer() const
Return the momentum transfer matrices.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the interfacial mass flow rate.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
virtual bool read()
Read base phaseProperties dictionary.
Namespace for OpenFOAM.