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 | Copyright (C) 2013-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::multiphaseSystem
26 
27 Description
28  Class which solves the volume fraction equations for two phases.
29 
30 SourceFiles
31  multiphaseSystem.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef multiphaseSystem_H
36 #define multiphaseSystem_H
37 
38 #include "phaseSystem.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 class dragModel;
46 class virtualMassModel;
47 
48 /*---------------------------------------------------------------------------*\
49  Class multiphaseSystem Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 class multiphaseSystem
53 :
54  public phaseSystem
55 {
56  // Private data
57 
58  volScalarField alphas_;
59 
60  typedef HashTable<scalar, phasePairKey, phasePairKey::hash>
61  cAlphaTable;
62 
63  cAlphaTable cAlphas_;
64 
65  //- Stabilisation for normalisation of the interface normal
66  const dimensionedScalar deltaN_;
67 
68  //- Conversion factor for degrees into radians
69  static const scalar convertToRad;
70 
71 
72  // Private member functions
73 
74  void calcAlphas();
75 
76  void solveAlphas();
77 
78  tmp<surfaceVectorField> nHatfv
79  (
80  const volScalarField& alpha1,
81  const volScalarField& alpha2
82  ) const;
83 
84  tmp<surfaceScalarField> nHatf
85  (
86  const volScalarField& alpha1,
87  const volScalarField& alpha2
88  ) const;
89 
90  void correctContactAngle
91  (
92  const phaseModel& alpha1,
93  const phaseModel& alpha2,
94  surfaceVectorField::Boundary& nHatb
95  ) const;
96 
97  tmp<volScalarField> K
98  (
99  const phaseModel& alpha1,
100  const phaseModel& alpha2
101  ) const;
102 
103 
104  //- Return the drag coefficient for phase pair
105  virtual tmp<volScalarField> Kd(const phasePairKey& key) const = 0;
106 
107  //- Return the face drag coefficient for phase pair
108  virtual tmp<surfaceScalarField> Kdf(const phasePairKey& key) const = 0;
109 
110  //- Return the virtual mass coefficient for phase pair
111  virtual tmp<volScalarField> Vm(const phasePairKey& key) const = 0;
112 
113  //- Return the face virtual mass coefficient for phase pair
114  virtual tmp<surfaceScalarField> Vmf(const phasePairKey& key) const = 0;
115 
116  //- Return the turbulent diffusivity for phase pair
117  // Multiplies the phase-fraction gradient
118  virtual tmp<volScalarField> D(const phasePairKey& key) const = 0;
119 
120  //- Return the interfacial mass flow rate for phase pair
121  virtual tmp<volScalarField> dmdt(const phasePairKey& key) const = 0;
122 
123 
124 public:
125 
126  //- Runtime type information
127  TypeName("multiphaseSystem");
128 
129  // Declare runtime construction
130 
132  (
133  autoPtr,
135  dictionary,
136  (
137  const fvMesh& mesh
138  ),
139  (mesh)
140  );
141 
142 
143  // Constructors
144 
145  //- Construct from fvMesh
146  multiphaseSystem(const fvMesh&);
147 
148 
149  //- Destructor
150  virtual ~multiphaseSystem();
151 
152 
153  // Selectors
154 
155  static autoPtr<multiphaseSystem> New
156  (
157  const fvMesh& mesh
158  );
159 
160 
161  // Member Functions
162 
163  //- Return the drag coefficient for all phase-pairs
164  virtual const phaseSystem::KdTable& Kds() const = 0;
165 
166  //- Return the drag coefficient for phase
167  virtual tmp<volScalarField> Kd(const phaseModel& phase) const = 0;
168 
169  //- Return the combined force (lift + wall-lubrication) for phase pair
170  virtual autoPtr<PtrList<Foam::volVectorField>> Fs() const = 0;
171 
172  //- Return the turbulent dispersion force on faces for phase pair
173  virtual autoPtr<PtrList<Foam::surfaceScalarField>> phiDs
174  (
175  const PtrList<volScalarField>& rAUs
176  ) const = 0;
177 
178  //- Return true if there is mass transfer for phase
179  virtual bool transfersMass(const phaseModel& phase) const = 0;
180 
181  //- Return the total interfacial mass transfer rate for phase
182  virtual tmp<volScalarField> dmdt(const phaseModel& phase) const = 0;
183 
184  //- Return the momentum transfer matrices
185  virtual autoPtr<momentumTransferTable> momentumTransfer() const = 0;
186 
187  //- Return the heat transfer matrices
188  virtual autoPtr<heatTransferTable> heatTransfer() const = 0;
189 
190  //- Return the mass transfer matrices
191  virtual autoPtr<massTransferTable> massTransfer() const = 0;
192 
193  tmp<surfaceScalarField> surfaceTension(const phaseModel& phase) const;
194 
195  //- Indicator of the proximity of the interface
196  // Field values are 1 near and 0 away for the interface.
197  tmp<volScalarField> nearInterface() const;
198 
199  //- Solve for the phase fractions
200  virtual void solve();
201 };
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 } // End namespace Foam
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 #include "multiphaseSystemI.H"
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 #endif
215 
216 // ************************************************************************* //
virtual autoPtr< PtrList< Foam::surfaceScalarField > > phiDs(const PtrList< volScalarField > &rAUs) const =0
Return the turbulent dispersion force on faces for phase pair.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
Definition: phaseSystem.H:82
static autoPtr< multiphaseSystem > New(const fvMesh &mesh)
virtual autoPtr< PtrList< Foam::volVectorField > > Fs() const =0
Return the combined force (lift + wall-lubrication) for phase pair.
virtual autoPtr< heatTransferTable > heatTransfer() const =0
Return the heat transfer matrices.
virtual autoPtr< massTransferTable > massTransfer() const =0
Return the mass transfer matrices.
virtual autoPtr< momentumTransferTable > momentumTransfer() const =0
Return the momentum transfer matrices.
CGAL::Exact_predicates_exact_constructions_kernel K
declareRunTimeSelectionTable(autoPtr, multiphaseSystem, dictionary,(const fvMesh &mesh),(mesh))
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:113
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual const phaseSystem::KdTable & Kds() const =0
Return the drag coefficient for all phase-pairs.
const fvMesh & mesh() const
Constant access the mesh.
Definition: phaseSystemI.H:28
PtrList< volScalarField > rAUs(fluid.phases().size())
alpha2
Definition: alphaEqn.H:112
void solve()
Solve for the mixture phase-fractions.
volScalarField & alpha1
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool transfersMass(const phaseModel &phase) const =0
Return true if there is mass transfer for phase.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
TypeName("multiphaseSystem")
Runtime type information.
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
virtual ~multiphaseSystem()
Destructor.
Namespace for OpenFOAM.