multiphaseSystem.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::multiphaseSystem
28 
29 Description
30  Class which solves the volume fraction equations for multiple phases
31 
32 SourceFiles
33  multiphaseSystem.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef multiphaseSystem_H
38 #define multiphaseSystem_H
39 
40 #include "multiphaseInterSystem.H"
41 #include "UPtrList.H"
42 #include "phasePairKey.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace multiphaseInter
49 {
50 /*---------------------------------------------------------------------------*\
51  Class multiphaseSystem Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class multiphaseSystem
55 :
57 {
58 protected:
59 
61 
64 
66 
67 
68  // Protected data
69 
70  //- Unallocated phase list
72 
73  //- Table for compression factors between phases
75 
76  //- Maximum volumen rate change
78 
79  //- Compression fluxed for phases
81 
82  //- Su phase source terms
84 
85  //- Sp phase source terms
86  SuSpTable Sp_;
87 
88 
89  // Protected members
90 
91  //- Calculate Sp and Su
92  void calculateSuSp();
93 
94  //- Solve alphas
95  void solveAlphas();
96 
97 
98 public:
99 
100  //- Runtime type information
101  TypeName("multiphaseSystem");
102 
103  // Declare runtime construction
104 
106  (
107  autoPtr,
109  dictionary,
110  (
111  const fvMesh& mesh
112  ),
113  (mesh)
114  );
115 
116 
117  // Constructors
118 
119  //- Construct from fvMesh
120  multiphaseSystem(const fvMesh&);
121 
122 
123  //- Destructor
124  virtual ~multiphaseSystem() = default;
125 
126 
127  // Selectors
128 
129  static autoPtr<multiphaseSystem> New(const fvMesh& mesh);
130 
131 
132  // Member Functions
133 
134  //- Solve for the phase fractions
135  virtual void solve();
136 
137 
138  // Access
139 
140  //- Return phases
141  const UPtrList<phaseModel>& phases() const;
142 
143  //- Return phases
145 
146  //- Constant access phase model i
147  const phaseModel& phase(const label i) const;
148 
149  //- Access phase model i
150  phaseModel& phase(const label i);
151 
152  //- Access to ddtAlphaMax
154 
155  //- Maximum diffusion number
156  scalar maxDiffNo() const;
157 
158  //- Access to compression fluxes for phaes
159  const compressionFluxTable& limitedPhiAlphas() const;
160 
161  //- Access Su
162  SuSpTable& Su();
163 
164  //- Access Sp
165  SuSpTable& Sp();
166 
167  //- Read thermophysical properties dictionary
168  virtual bool read();
169 };
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace Foam
174 } // End namespace multiphaseInter
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 #endif
179 
180 // ************************************************************************* //
const fvMesh & mesh() const
Return mesh.
UPtrList< phaseModel > phases_
Unallocated phase list.
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
virtual bool read()
Read thermophysical properties dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual ~multiphaseSystem()=default
Destructor.
HashTable< surfaceScalarField > compressionFluxTable
declareRunTimeSelectionTable(autoPtr, multiphaseSystem, dictionary,(const fvMesh &mesh),(mesh))
virtual void solve()
Solve for the phase fractions.
multiphaseSystem(const fvMesh &)
Construct from fvMesh.
const phaseModel & phase(const label i) const
Constant access phase model i.
dimensionedScalar ddtAlphaMax() const
Access to ddtAlphaMax.
TypeName("multiphaseSystem")
Runtime type information.
const compressionFluxTable & limitedPhiAlphas() const
Access to compression fluxes for phaes.
dimensionedScalar ddtAlphaMax_
Maximum volumen rate change.
scalarTable cAlphas_
Table for compression factors between phases.
SuSpTable Su_
Su phase source terms.
compressionFluxTable limitedPhiAlphas_
Compression fluxed for phases.
SuSpTable Sp_
Sp phase source terms.
HashTable< volScalarField::Internal > SuSpTable
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
scalar maxDiffNo() const
Maximum diffusion number.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const UPtrList< phaseModel > & phases() const
Return phases.
static autoPtr< multiphaseSystem > New(const fvMesh &mesh)
void calculateSuSp()
Calculate Sp and Su.
Namespace for OpenFOAM.