SchnerrSauer.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 | Copyright (C) 2011-2016 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::phaseChangeTwoPhaseMixtures::SchnerrSauer
26 
27 Description
28  SchnerrSauer cavitation model.
29 
30  Reference:
31  \verbatim
32  Schnerr, G. H., And Sauer, J.,
33  "Physical and Numerical Modeling of Unsteady Cavitation Dynamics",
34  Proc. 4th International Conference on Multiphase Flow,
35  New Orleans, U.S.A., 2001.
36  \endverbatim
37 
38 SourceFiles
39  SchnerrSauer.C
40 
41 \*--------------------------------------------------------------------*/
42 
43 #ifndef SchnerrSauer_H
44 #define SchnerrSauer_H
45 
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 namespace phaseChangeTwoPhaseMixtures
53 {
54 
55 /*--------------------------------------------------------------------*\
56  Class SchnerrSauer
57 \*--------------------------------------------------------------------*/
58 
59 class SchnerrSauer
60 :
62 {
63  // Private data
64 
65  //- Bubble number density
67 
68  //- Nucleation site diameter
69  dimensionedScalar dNuc_;
70 
71  //- Condensation rate coefficient
73 
74  //- Vapourisation rate coefficient
76 
78 
79  //- Nucleation site volume-fraction
80  dimensionedScalar alphaNuc() const;
81 
82  //- Reciprocal bubble radius
83  tmp<volScalarField>rRb(const volScalarField& limitedAlpha1) const;
84 
85  //- Part of the condensation and vapourisation rates
86  tmp<volScalarField> pCoeff(const volScalarField& p) const;
87 
88 
89 public:
90 
91  //- Runtime type information
92  TypeName("SchnerrSauer");
93 
94 
95  // Constructors
96 
97  //- Construct from components
99  (
100  const volVectorField& U,
101  const surfaceScalarField& phi
102  );
103 
104 
105  //- Destructor
106  virtual ~SchnerrSauer()
107  {}
108 
109 
110  // Member Functions
111 
112  //- Return the mass condensation and vaporisation rates as a
113  // coefficient to multiply (1 - alphal) for the condensation rate
114  // and a coefficient to multiply alphal for the vaporisation rate
115  virtual Pair<tmp<volScalarField>> mDotAlphal() const;
116 
117  //- Return the mass condensation and vaporisation rates as coefficients
118  // to multiply (p - pSat)
119  virtual Pair<tmp<volScalarField>> mDotP() const;
120 
121  //- Correct the SchnerrSauer phaseChange model
122  virtual void correct();
123 
124  //- Read the transportProperties dictionary and update
125  virtual bool read();
126 };
127 
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 } // End namespace phaseChangeTwoPhaseMixtures
132 } // End namespace Foam
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 #endif
137 
138 // ************************************************************************* //
surfaceScalarField & phi
virtual void correct()
Correct the SchnerrSauer phaseChange model.
TypeName("SchnerrSauer")
Runtime type information.
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
virtual bool read()
Read the transportProperties dictionary and update.
const volVectorField & U() const
Return const-access to the mixture velocity.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:54
SchnerrSauer(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
SchnerrSauer cavitation model.
Definition: SchnerrSauer.H:58
Namespace for OpenFOAM.