Merkle.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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::Merkle
26 
27 Description
28  Merkle cavitation model.
29 
30  Reference:
31  \verbatim
32  C. L. Merkle, J. Feng, and P. E. O. Buelow,
33  "Computational modeling of the dynamics of sheet cavitation",
34  in Proceedings Third International Symposium on Cavitation
35  Grenoble, France 1998.
36  \endverbatim
37 
38 SourceFiles
39  Merkle.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef Merkle_H
44 #define Merkle_H
45 
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 namespace phaseChangeTwoPhaseMixtures
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class Merkle
57 \*---------------------------------------------------------------------------*/
58 
59 class Merkle
60 :
62 {
63  // Private data
64 
65  dimensionedScalar UInf_;
66  dimensionedScalar tInf_;
69 
71 
72  dimensionedScalar mcCoeff_;
73  dimensionedScalar mvCoeff_;
74 
75 
76 public:
77 
78  //- Runtime type information
79  TypeName("Merkle");
80 
81 
82  // Constructors
83 
84  //- Construct from components
85  Merkle
86  (
87  const volVectorField& U,
88  const surfaceScalarField& phi
89  );
90 
91 
92  //- Destructor
93  virtual ~Merkle()
94  {}
95 
96 
97  // Member Functions
98 
99  //- Return the mass condensation and vaporisation rates as a
100  // coefficient to multiply (1 - alphal) for the condensation rate
101  // and a coefficient to multiply alphal for the vaporisation rate
102  virtual Pair<tmp<volScalarField>> mDotAlphal() const;
103 
104  //- Return the mass condensation and vaporisation rates as coefficients
105  // to multiply (p - pSat)
106  virtual Pair<tmp<volScalarField>> mDotP() const;
107 
108  //- Correct the Merkle phaseChange model
109  virtual void correct();
110 
111  //- Read the transportProperties dictionary and update
112  virtual bool read();
113 };
114 
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 } // End namespace phaseChangeTwoPhaseMixtures
119 } // End namespace Foam
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 #endif
124 
125 // ************************************************************************* //
surfaceScalarField & phi
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
virtual ~Merkle()
Destructor.
Definition: Merkle.H:92
Merkle(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const volVectorField & U() const
Return const-access to the mixture velocity.
TypeName("Merkle")
Runtime type information.
virtual void correct()
Correct the Merkle phaseChange model.
Merkle cavitation model.
Definition: Merkle.H:58
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual bool read()
Read the transportProperties dictionary and update.
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
Namespace for OpenFOAM.