Removed unneeded lib dependency from mdInitialise
[foam-extend-3.2.git] / applications / solvers / multiphase / multiphaseInterFoam / multiphaseMixture / multiphaseMixture.H
blob9940742317fd6ba55e4d9a7a791700c0d5811068
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::multiphaseMixture
28 Description
29     Incompressible multi-phase mixture with built in solution for the
30     phase fractions with interface compression for interface-capturing.
32     Derived from transportModel so that it can be unsed in conjunction with
33     the incompressible turbulence models.
35     Surface tension and contact-angle is handled for the interface
36     between each phase-pair.
38 SourceFiles
39     multiphaseMixture.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef multiphaseMixture_H
44 #define multiphaseMixture_H
46 #include "incompressible/transportModel/transportModel.H"
47 #include "phase.H"
48 #include "PtrDictionary.H"
49 #include "volFields.H"
50 #include "surfaceFields.H"
51 #include "multivariateSurfaceInterpolationScheme.H"
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 namespace Foam
59 /*---------------------------------------------------------------------------*\
60                            Class multiphaseMixture Declaration
61 \*---------------------------------------------------------------------------*/
63 class multiphaseMixture
65     public transportModel
67 public:
69     class interfacePair
70     :
71         public Pair<word>
72     {
73     public:
75         class hash
76         :
77             public Hash<interfacePair>
78         {
79         public:
81             hash()
82             {}
84             label operator()(const interfacePair& key) const
85             {
86                 return word::hash()(key.first()) + word::hash()(key.second());
87             }
89             label operator()
90             (
91                 const interfacePair& key,
92                 const label tableSize
93             ) const
94             {
95                 return mag(operator()(key)) % tableSize;
96             }
97         };
100         // Constructors
102             interfacePair()
103             {}
105             interfacePair(const word& alpha1Name, const word& alpha2Name)
106             :
107                 Pair<word>(alpha1Name, alpha2Name)
108             {}
110             interfacePair(const phase& alpha1, const phase& alpha2)
111             :
112                 Pair<word>(alpha1.name(), alpha2.name())
113             {}
116         // Friend Operators
118             friend bool operator==
119             (
120                 const interfacePair& a,
121                 const interfacePair& b
122             )
123             {
124                 return
125                 (
126                     ((a.first() == b.first()) && (a.second() == b.second()))
127                  || ((a.first() == b.second()) && (a.second() == b.first()))
128                 );
129             }
131             friend bool operator!=
132             (
133                 const interfacePair& a,
134                 const interfacePair& b
135             )
136             {
137                 return (!(a == b));
138             }
139     };
142 private:
144     // Private data
146         //- Dictionary of phases
147         PtrDictionary<phase> phases_;
149         //- The phase chosen as reference, the one which is derived from
150         //  the others such thatr they sum to 1
151         phase& refPhase_;
153         const fvMesh& mesh_;
154         const volVectorField& U_;
155         const surfaceScalarField& phi_;
157         volScalarField nu_;
158         surfaceScalarField rhoPhi_;
160         volScalarField alphas_;
162         typedef HashTable<scalar, interfacePair, interfacePair::hash>
163             sigmaTable;
165         sigmaTable sigmas_;
166         dimensionSet dimSigma_;
168         //- Stabilisation for normalisation of the interface normal
169         const dimensionedScalar deltaN_;
171         //- Conversion factor for degrees into radians
172         static const scalar convertToRad;
174         //- Phase-fraction field table for multivariate discretisation
175         multivariateSurfaceInterpolationScheme<scalar>::fieldTable alphaTable_;
178     // Private member functions
180         void calcAlphas();
182         void solveAlphas
183         (
184             const label nAlphaCorr,
185             const bool cycleAlpha,
186             const scalar cAlpha
187         );
189         tmp<surfaceVectorField> nHatfv
190         (
191             const volScalarField& alpha1,
192             const volScalarField& alpha2
193         ) const;
195         tmp<surfaceScalarField> nHatf
196         (
197             const volScalarField& alpha1,
198             const volScalarField& alpha2
199         ) const;
201         void correctContactAngle
202         (
203             const phase& alpha1,
204             const phase& alpha2,
205             surfaceVectorField::GeometricBoundaryField& nHatb
206         ) const;
208         tmp<volScalarField> K(const phase& alpha1, const phase& alpha2) const;
211 public:
213     // Constructors
215         //- Construct from components
216         multiphaseMixture
217         (
218             const volVectorField& U,
219             const surfaceScalarField& phi
220         );
223     // Destructor
225         ~multiphaseMixture()
226         {}
229     // Member Functions
231         //- Return the phases
232         const PtrDictionary<phase>& phases() const
233         {
234             return phases_;
235         }
237         //- Return the velocity
238         const volVectorField& U() const
239         {
240             return U_;
241         }
243         //- Return the volumetric flux
244         const surfaceScalarField& phi() const
245         {
246             return phi_;
247         }
249         const surfaceScalarField& rhoPhi() const
250         {
251             return rhoPhi_;
252         }
254         //- Return the mixture density
255         tmp<volScalarField> rho() const;
257         //- Return the dynamic laminar viscosity
258         tmp<volScalarField> mu() const;
260         //- Return the face-interpolated dynamic laminar viscosity
261         tmp<surfaceScalarField> muf() const;
263         //- Return the kinematic laminar viscosity
264         const volScalarField& nu() const;
266         //- Return the face-interpolated dynamic laminar viscosity
267         tmp<surfaceScalarField> nuf() const;
269         tmp<surfaceScalarField> surfaceTensionForce() const;
271         //- Correct the mixture properties
272         void correct();
274         //- Read base transportProperties dictionary
275         bool read();
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 } // End namespace Foam
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 #endif
287 // ************************************************************************* //