multiphaseSystem.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) 2013-2018 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:
57 
58  // Private typedefs
59 
60  typedef HashTable<scalar, phasePairKey, phasePairKey::hash>
61  cAlphaTable;
62 
63 
64  // Private data
65 
66  //- The indexed phase-fraction field; e.g., 1 for water, 2 for air, 3
67  // for oil, etc...
68  volScalarField alphas_;
69 
70  //-
71  cAlphaTable cAlphas_;
72 
73  //- Stabilisation for normalisation of the interface normal
74  const dimensionedScalar deltaN_;
75 
76  //- Conversion factor for degrees into radians
77  static const scalar convertToRad;
78 
79 
80  // Private member functions
81 
82  void calcAlphas();
83 
84  void solveAlphas();
85 
86  tmp<surfaceVectorField> nHatfv
87  (
88  const volScalarField& alpha1,
89  const volScalarField& alpha2
90  ) const;
91 
92  tmp<surfaceScalarField> nHatf
93  (
94  const volScalarField& alpha1,
95  const volScalarField& alpha2
96  ) const;
97 
98  void correctContactAngle
99  (
100  const phaseModel& alpha1,
101  const phaseModel& alpha2,
102  surfaceVectorField::Boundary& nHatb
103  ) const;
104 
105  tmp<volScalarField> K
106  (
107  const phaseModel& alpha1,
108  const phaseModel& alpha2
109  ) const;
110 
111  //- Return the drag coefficient for phase pair
112  virtual tmp<volScalarField> Kd(const phasePairKey& key) const = 0;
113 
114  //- Return the virtual mass coefficient for phase pair
115  virtual tmp<volScalarField> Vm(const phasePairKey& key) const = 0;
116 
117 
118 protected:
119 
120  // Protected data
121 
122  //- Flag to indicate that returned lists of fields are "complete"; i.e.,
123  // that an absence of force is returned as a constructed list of zeros,
124  // rather than a null pointer
125  static const bool fillFields_ = false;
126 
127 
128 public:
129 
130  //- Runtime type information
131  TypeName("multiphaseSystem");
132 
133  // Declare runtime construction
134 
136  (
137  autoPtr,
139  dictionary,
140  (
141  const fvMesh& mesh
142  ),
143  (mesh)
144  );
145 
146 
147  // Constructors
148 
149  //- Construct from fvMesh
150  multiphaseSystem(const fvMesh&);
151 
152 
153  //- Destructor
154  virtual ~multiphaseSystem();
155 
156 
157  // Selectors
158 
160  (
161  const fvMesh& mesh
162  );
163 
164 
165  // Member Functions
166 
167  using phaseSystem::sigma;
168  using phaseSystem::dmdts;
169 
170  //- Return the surface tension force
172 
173  //- Indicator of the proximity of the interface
174  // Field values are 1 near and 0 away for the interface.
176 
177  //- Solve for the phase fractions
178  virtual void solve();
179 };
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace Foam
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #include "multiphaseSystemI.H"
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #endif
193 
194 // ************************************************************************* //
static autoPtr< multiphaseSystem > New(const fvMesh &mesh)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
CGAL::Exact_predicates_exact_constructions_kernel K
declareRunTimeSelectionTable(autoPtr, multiphaseSystem, dictionary,(const fvMesh &mesh),(mesh))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
const volScalarField & alpha1
alpha2
Definition: alphaEqn.H:115
void solve()
Solve for the mixture phase-fractions.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
TypeName("multiphaseSystem")
Runtime type information.
A class for managing temporary objects.
Definition: PtrList.H:53
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Incompressible multi-phase mixture with built in solution for the phase fractions with interface comp...
virtual ~multiphaseSystem()
Destructor.
Namespace for OpenFOAM.
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.