Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / solvers / multiphase / multiphaseInterFoam / multiphaseMixture / multiphaseMixture.H
blob8a2eff976c904c81c626523daf1d41997b7f2106
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::multiphaseMixture
27 Description
28     Incompressible multi-phase mixture with built in solution for the
29     phase fractions with interface compression for interface-capturing.
31     Derived from transportModel so that it can be unsed in conjunction with
32     the incompressible turbulence models.
34     Surface tension and contact-angle is handled for the interface
35     between each phase-pair.
37 SourceFiles
38     multiphaseMixture.C
40 \*---------------------------------------------------------------------------*/
42 #ifndef multiphaseMixture_H
43 #define multiphaseMixture_H
45 #include "incompressible/transportModel/transportModel.H"
46 #include "phase.H"
47 #include "PtrDictionary.H"
48 #include "volFields.H"
49 #include "surfaceFields.H"
50 #include "multivariateSurfaceInterpolationScheme.H"
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 namespace Foam
58 /*---------------------------------------------------------------------------*\
59                            Class multiphaseMixture Declaration
60 \*---------------------------------------------------------------------------*/
62 class multiphaseMixture
64     public transportModel
66 public:
68     class interfacePair
69     :
70         public Pair<word>
71     {
72     public:
74         class hash
75         :
76             public Hash<interfacePair>
77         {
78         public:
80             hash()
81             {}
83             label operator()(const interfacePair& key) const
84             {
85                 return word::hash()(key.first()) + word::hash()(key.second());
86             }
88             label operator()
89             (
90                 const interfacePair& key,
91                 const label tableSize
92             ) const
93             {
94                 return mag(operator()(key)) % tableSize;
95             }
96         };
99         // Constructors
101             interfacePair()
102             {}
104             interfacePair(const word& alpha1Name, const word& alpha2Name)
105             :
106                 Pair<word>(alpha1Name, alpha2Name)
107             {}
109             interfacePair(const phase& alpha1, const phase& alpha2)
110             :
111                 Pair<word>(alpha1.name(), alpha2.name())
112             {}
115         // Friend Operators
117             friend bool operator==
118             (
119                 const interfacePair& a,
120                 const interfacePair& b
121             )
122             {
123                 return
124                 (
125                     ((a.first() == b.first()) && (a.second() == b.second()))
126                  || ((a.first() == b.second()) && (a.second() == b.first()))
127                 );
128             }
130             friend bool operator!=
131             (
132                 const interfacePair& a,
133                 const interfacePair& b
134             )
135             {
136                 return (!(a == b));
137             }
138     };
141 private:
143     // Private data
145         //- Dictionary of phases
146         PtrDictionary<phase> phases_;
148         //- The phase chosen as reference, the one which is derived from
149         //  the others such thatr they sum to 1
150         phase& refPhase_;
152         const fvMesh& mesh_;
153         const volVectorField& U_;
154         const surfaceScalarField& phi_;
156         volScalarField nu_;
157         surfaceScalarField rhoPhi_;
159         volScalarField alphas_;
161         typedef HashTable<scalar, interfacePair, interfacePair::hash>
162             sigmaTable;
164         sigmaTable sigmas_;
165         dimensionSet dimSigma_;
167         //- Stabilisation for normalisation of the interface normal
168         const dimensionedScalar deltaN_;
170         //- Conversion factor for degrees into radians
171         static const scalar convertToRad;
173         //- Phase-fraction field table for multivariate discretisation
174         multivariateSurfaceInterpolationScheme<scalar>::fieldTable alphaTable_;
177     // Private member functions
179         void calcAlphas();
181         void solveAlphas
182         (
183             const label nAlphaCorr,
184             const bool cycleAlpha,
185             const scalar cAlpha
186         );
188         tmp<surfaceVectorField> nHatfv
189         (
190             const volScalarField& alpha1,
191             const volScalarField& alpha2
192         ) const;
194         tmp<surfaceScalarField> nHatf
195         (
196             const volScalarField& alpha1,
197             const volScalarField& alpha2
198         ) const;
200         void correctContactAngle
201         (
202             const phase& alpha1,
203             const phase& alpha2,
204             surfaceVectorField::GeometricBoundaryField& nHatb
205         ) const;
207         tmp<volScalarField> K(const phase& alpha1, const phase& alpha2) const;
210 public:
212     // Constructors
214         //- Construct from components
215         multiphaseMixture
216         (
217             const volVectorField& U,
218             const surfaceScalarField& phi
219         );
222     // Destructor
224         ~multiphaseMixture()
225         {}
228     // Member Functions
230         //- Return the phases
231         const PtrDictionary<phase>& phases() const
232         {
233             return phases_;
234         }
236         //- Return the velocity
237         const volVectorField& U() const
238         {
239             return U_;
240         }
242         //- Return the volumetric flux
243         const surfaceScalarField& phi() const
244         {
245             return phi_;
246         }
248         const surfaceScalarField& rhoPhi() const
249         {
250             return rhoPhi_;
251         }
253         //- Return the mixture density
254         tmp<volScalarField> rho() const;
256         //- Return the dynamic laminar viscosity
257         tmp<volScalarField> mu() const;
259         //- Return the face-interpolated dynamic laminar viscosity
260         tmp<surfaceScalarField> muf() const;
262         //- Return the kinematic laminar viscosity
263         const volScalarField& nu() const;
265         //- Return the face-interpolated dynamic laminar viscosity
266         tmp<surfaceScalarField> nuf() const;
268         tmp<surfaceScalarField> surfaceTensionForce() const;
270         //- Correct the mixture properties
271         void correct();
273         //- Read base transportProperties dictionary
274         bool read();
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 } // End namespace Foam
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 #endif
286 // ************************************************************************* //