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