PopulationBalancePhaseSystem.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) 2017-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::PopulationBalancePhaseSystem
26 
27 Description
28  Class which provides population balance functionality. Stores the mass
29  transfer rates resulting from coalescence, breakup or drift across
30  representative phases that collectively define a dispersed phase.
31 
32 See also
33  Foam::diameterModels::populationBalanceModel
34 
35 SourceFiles
36  PopulationBalancePhaseSystem.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef PopulationBalancePhaseSystem_H
41 #define PopulationBalancePhaseSystem_H
42 
43 #include "phaseSystem.H"
44 #include "populationBalanceModel.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class PopulationBalancePhaseSystem Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class BasePhaseSystem>
57 :
58  public BasePhaseSystem
59 {
60  // Private data
61 
62  //- Population balances
64 
65 
66 public:
67 
68  // Constructors
69 
70  //- Construct from fvMesh
72 
73 
74  //- Destructor
76 
77 
78  // Member Functions
79 
80  //- Return the mass transfer rate for an interface
81  virtual tmp<volScalarField> dmdtf(const phaseInterfaceKey& key) const;
82 
83  //- Return the mass transfer rates for each phase
84  virtual PtrList<volScalarField> dmdts() const;
85 
86  //- Return the momentum transfer matrices for the cell-based algorithm
88 
89  //- Return the momentum transfer matrices for the face-based algorithm
91 
92  //- Return the heat transfer matrices
94 
95  //- Return the specie transfer matrices
97  specieTransfer() const;
98 
99  //- Solve all population balance equations
100  virtual void solve
101  (
102  const PtrList<volScalarField>& rAUs,
103  const PtrList<surfaceScalarField>& rAUfs
104  );
105 
106  //- Correct derived properties
107  virtual void correct();
108 
109  //- Read base phaseProperties dictionary
110  virtual bool read();
111 };
112 
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 } // End namespace Foam
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 #ifdef NoRepository
122 #endif
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 #endif
127 
128 // ************************************************************************* //
Class which provides population balance functionality. Stores the mass transfer rates resulting from ...
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 derived properties.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve all population balance equations.
PopulationBalancePhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
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.