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) 2013-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::multiphaseSystem
28
29Description
30 Class which solves the volume fraction equations for two phases.
31
32SourceFiles
33 multiphaseSystem.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef multiphaseSystem_H
38#define multiphaseSystem_H
39
40#include "phaseSystem.H"
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44namespace Foam
45{
46
47// Forward Declarations
48class dragModel;
49class virtualMassModel;
50
51/*---------------------------------------------------------------------------*\
52 Class multiphaseSystem Declaration
53\*---------------------------------------------------------------------------*/
54
55class multiphaseSystem
56:
57 public phaseSystem
58{
59 // Private Typedefs
60
61 typedef HashTable<scalar, phasePairKey, phasePairKey::hash>
62 cAlphaTable;
63
64
65 // Private Data
66
67 //- The indexed phase-fraction field;
68 // E.g., 1 for water, 2 for air, 3 for oil, etc...
69 volScalarField alphas_;
70
71 //-
72 cAlphaTable cAlphas_;
73
74 //- Stabilisation for normalisation of the interface normal
75 const dimensionedScalar deltaN_;
76
77
78 // Private Member Functions
79
80 void calcAlphas();
81
82 void solveAlphas();
83
84 tmp<surfaceVectorField> nHatfv
85 (
88 ) const;
89
90 tmp<surfaceScalarField> nHatf
91 (
94 ) const;
95
96 void correctContactAngle
97 (
98 const phaseModel& alpha1,
99 const phaseModel& alpha2,
101 ) const;
102
103 tmp<volScalarField> K
104 (
105 const phaseModel& alpha1,
106 const phaseModel& alpha2
107 ) const;
108
109 //- Return the drag coefficient for phase pair
110 virtual tmp<volScalarField> Kd(const phasePairKey& key) const = 0;
111
112 //- Return the virtual mass coefficient for phase pair
113 virtual tmp<volScalarField> Vm(const phasePairKey& key) const = 0;
114
115
116protected:
117
118 // Protected data
119
120 //- Flag to indicate that returned lists of fields are "complete"; i.e.,
121 // that an absence of force is returned as a constructed list of zeros,
122 // rather than a null pointer
123 static const bool fillFields_ = false;
124
125
126public:
127
128 //- Runtime type information
129 TypeName("multiphaseSystem");
130
131 // Declare runtime construction
134 (
135 autoPtr,
138 (
139 const fvMesh& mesh
140 ),
141 (mesh)
142 );
143
144
145 // Constructors
146
147 //- Construct from fvMesh
148 explicit multiphaseSystem(const fvMesh& mesh);
149
150
151 //- Destructor
152 virtual ~multiphaseSystem() = default;
153
154
155 // Selectors
156
158 (
159 const fvMesh& mesh
160 );
161
162
163 // Member Functions
164
165 using phaseSystem::sigma;
166 using phaseSystem::dmdts;
167
168 //- Return the surface tension force
170
171 //- Indicator of the proximity of the interface
172 // Field values are 1 near and 0 away for the interface.
174
175 //- Solve for the phase fractions
176 virtual void solve();
177};
178
179
180// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181
182} // End namespace Foam
183
184// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185
186#endif
187
188// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
const volScalarField & alpha1
const volScalarField & alpha2
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Incompressible multi-phase mixture with built in solution for the phase fractions with interface comp...
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
declareRunTimeSelectionTable(autoPtr, multiphaseSystem, dictionary,(const fvMesh &mesh),(mesh))
static autoPtr< multiphaseSystem > New(const fvMesh &mesh)
static const bool fillFields_
Flag to indicate that returned lists of fields are "complete"; i.e.,.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
Return the surface tension force.
virtual void solve()
Solve for the phase fractions.
virtual ~multiphaseSystem()=default
Destructor.
void solve()
Solve for the mixture phase-fractions.
TypeName("multiphaseSystem")
Runtime type information.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: phaseSystem.C:355
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:30
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73