PhaseTransferPhaseSystem.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018-2023 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::PhaseTransferPhaseSystem
26 
27 Description
28  Class which models non-thermally-coupled or weakly thermally coupled
29  mass transfers.
30 
31 SourceFiles
32  PhaseTransferPhaseSystem.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef PhaseTransferPhaseSystem_H
37 #define PhaseTransferPhaseSystem_H
38 
39 #include "phaseSystem.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 class blendedPhaseTransferModel;
47 
48 /*---------------------------------------------------------------------------*\
49  Class PhaseTransferPhaseSystem Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 template<class BasePhaseSystem>
54 :
55  public BasePhaseSystem
56 {
57 private:
58 
59  // Private typedefs
60 
61  typedef HashTable
62  <
67 
68 
69  // Private data
70 
71  // Sub Models
72 
73  //- Mass transfer models
74  phaseTransferModelTable phaseTransferModels_;
75 
76  //- Bulk mass transfer rates
78 
79  //- Bulk mass transfer pressure linearisation coefficients
80  phaseSystem::dmdtfTable d2mdtdpfs_;
81 
82  //- Specie mass transfer rates
83  phaseSystem::dmidtfTable dmidtfs_;
84 
85 
86  // Private Member Functions
87 
88  //- Return the summed mass transfers across each interface
89  autoPtr<phaseSystem::dmdtfTable> totalDmdtfs() const;
90 
91 
92 protected:
93 
94  // Protected member functions
95 
96  //- Add specie transfer terms which result from bulk mass transfers
97  void addDmdtYfs
98  (
99  const phaseSystem::dmdtfTable& dmdtfs,
101  ) const;
102 
103  //- Add specie transfer terms which result from specie mass transfers
104  void addDmidtYf
105  (
106  const phaseSystem::dmidtfTable& dmidtfs,
108  ) const;
109 
110 
111 public:
112 
113  // Constructors
114 
115  //- Construct from fvMesh
117 
118 
119  //- Destructor
120  virtual ~PhaseTransferPhaseSystem();
121 
122 
123  // Member Functions
124 
125  //- Return the mass transfer rate for an interface
126  virtual tmp<volScalarField> dmdtf(const phaseInterfaceKey& key) const;
127 
128  //- Return the mass transfer rates for each phase
129  virtual PtrList<volScalarField> dmdts() const;
130 
131  //- Return the mass transfer pressure linearisation coefficients for
132  // each phase
133  virtual PtrList<volScalarField> d2mdtdps() const;
134 
135  //- Return the momentum transfer matrices for the cell-based algorithm
137 
138  //- Return the momentum transfer matrices for the face-based algorithm
140 
141  //- Return the heat transfer matrices
143 
144  //- Return the specie transfer matrices
146  specieTransfer() const;
147 
148  //- Correct the mass transfer rates
149  virtual void correct();
150 
151  //- Read base phaseProperties dictionary
152  virtual bool read();
153 };
154 
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 } // End namespace Foam
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 #ifdef NoRepository
163  #include "PhaseTransferPhaseSystem.C"
164 #endif
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #endif
169 
170 // ************************************************************************* //
An STL-conforming hash table.
Definition: HashTable.H:127
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual void correct()
Correct the mass transfer rates.
void addDmidtYf(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from specie mass transfers.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual ~PhaseTransferPhaseSystem()
Destructor.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer pressure linearisation coefficients for.
virtual bool read()
Read base phaseProperties dictionary.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Word-pair based class used for keying interface models in hash tables.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.