Loading...
Searching...
No Matches
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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::multiphaseSystem
29
30Description
31 Incompressible multi-phase mixture with built in solution for the
32 phase fractions with interface compression for interface-capturing.
33
34 Derived from transportModel so that it can be unused in conjunction with
35 the incompressible turbulence models.
36
37 Surface tension and contact-angle is handled for the interface
38 between each phase-pair.
39
40SourceFiles
41 multiphaseSystem.C
42
43\*---------------------------------------------------------------------------*/
44
45#ifndef Foam_multiphaseEuler_multiphaseSystem_H
46#define Foam_multiphaseEuler_multiphaseSystem_H
47
49#include "IOdictionary.H"
50#include "phaseModel.H"
51#include "PtrDictionary.H"
52#include "volFields.H"
53#include "surfaceFields.H"
54#include "dragModel.H"
55#include "HashPtrTable.H"
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59namespace Foam
60{
62/*---------------------------------------------------------------------------*\
63 Class multiphaseSystem Declaration
64\*---------------------------------------------------------------------------*/
65
67:
68 public IOdictionary,
69 public transportModel
70{
71public:
72
73 //- Name pair for the interface
74 class interfacePair
75 :
76 public Pair<word>
77 {
78 public:
79
80 // Ordered hashing (alias)
81 using hash = Pair<word>::hasher;
82
83 // Unordered hashing (alias)
85
87 // Constructors
89 interfacePair() = default;
90
91 interfacePair(const word& alpha1Name, const word& alpha2Name)
92 :
93 Pair<word>(alpha1Name, alpha2Name)
94 {}
95
97 :
99 {}
100
102 // Friend Operators
103
104 friend bool operator==
105 (
106 const interfacePair& a,
107 const interfacePair& b
108 )
109 {
110 return (0 != Pair<word>::compare(a, b));
111 }
112
113 friend bool operator!=
114 (
115 const interfacePair& a,
116 const interfacePair& b
117 )
118 {
119 return (!(a == b));
120 }
121 };
122
123
124 typedef
126 <
129
132
133
134private:
135
136 // Private data
137
138 //- Dictionary of phases
140
141 const fvMesh& mesh_;
142 const surfaceScalarField& phi_;
143
144 volScalarField alphas_;
145
147 scalarCoeffSymmTable;
148
150 scalarCoeffTable;
151
152 scalarCoeffSymmTable sigmas_;
153 dimensionSet dimSigma_;
154
155 scalarCoeffSymmTable cAlphas_;
156
157 scalarCoeffTable Cvms_;
158
160 interfaceDictTable;
161
162 dragModelTable dragModels_;
163
164 //- Stabilisation for normalisation of the interface normal
165 const dimensionedScalar deltaN_;
166
167
168 // Private member functions
169
170 void calcAlphas();
171
172 void solveAlphas();
173
175 (
176 const volScalarField& alpha1,
178 ) const;
179
181 (
182 const volScalarField& alpha1,
184 ) const;
185
186 void correctContactAngle
187 (
188 const phaseModel& alpha1,
189 const phaseModel& alpha2,
191 ) const;
192
194 (
195 const phaseModel& alpha1,
196 const phaseModel& alpha2
197 ) const;
198
199
200public:
201
202 // Constructors
203
204 //- Construct from components
206 (
207 const volVectorField& U,
209 );
210
211
212 //- Destructor
213 virtual ~multiphaseSystem() = default;
214
215
216 // Member Functions
217
218 //- Return the phases
219 const PtrDictionary<phaseModel>& phases() const
220 {
221 return phases_;
222 }
223
224 //- Return the phases
227 return phases_;
228 }
229
230 //- Return the mixture density
231 tmp<volScalarField> rho() const;
232
233 //- Return the mixture density for patch
234 tmp<scalarField> rho(const label patchi) const;
235
236 //- Return the mixture laminar viscosity
237 tmp<volScalarField> nu() const;
238
239 //- Return the laminar viscosity for patch
240 tmp<scalarField> nu(const label patchi) const;
241
242 //- Return the virtual-mass coefficient for the given phase
244
245 //- Return the virtual-mass source for the given phase
247
248 //- Return the table of drag models
249 const dragModelTable& dragModels() const
250 {
251 return dragModels_;
252 }
253
254 //- Return the drag coefficients for all of the interfaces
256
257 //- Return the sum of the drag coefficients for the given phase
259 (
260 const phaseModel& phase,
262 ) const;
263
265
266 //- Indicator of the proximity of the interface
267 // Field values are 1 near and 0 away for the interface.
269
270 //- Solve for the mixture phase-fractions
271 void solve();
273 //- Dummy correct
274 void correct()
275 {}
276
277 //- Read base transportProperties dictionary
278 bool read();
279};
280
281
282// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283
284} // End namespace Foam
285
286// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287
288#endif
289
290// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
const volScalarField & alpha1
const volScalarField & alpha2
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
static int compare(const Pair< T > &a, const Pair< T > &b)
Compare Pairs.
Definition PairI.H:24
Template dictionary class which manages the storage associated with it.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Name pair for the interface.
interfacePair(const word &alpha1Name, const word &alpha2Name)
interfacePair(const phaseModel &alpha1, const phaseModel &alpha2)
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
PtrDictionary< phaseModel > & phases()
Return the phases.
HashPtrTable< volScalarField, interfacePair, interfacePair::symmHash > dragCoeffFields
void correct()
Dummy correct.
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
const PtrDictionary< phaseModel > & phases() const
Return the phases.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
tmp< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient for the given phase.
const dragModelTable & dragModels() const
Return the table of drag models.
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
virtual ~multiphaseSystem()=default
Destructor.
void solve()
Solve for the mixture phase-fractions.
bool read()
Read base transportProperties dictionary.
HashPtrTable< multiphaseEuler::dragModel, interfacePair, interfacePair::symmHash > dragModelTable
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const surfaceScalarField & phi() const
Return the mixture flux.
tmp< volVectorField > U() const
Return the mixture velocity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & b
Symmetric hashing functor for Pair, hashes lower value first.
Definition Pair.H:202
Foam::surfaceFields.