This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/multiphase/interPhaseChangeFoam/phaseChangeTwoPhaseMixtures/Merkle/Merkle.H

129 lines
3.6 KiB
C++
Raw Normal View History

/*---------------------------------------------------------------------------*\
========Merkle= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.1
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
General Public License for more details.
2013-12-11 16:09:41 +00:00
You should have received a copy of the GNU General Public License
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::phaseChangeTwoPhaseMixtures::Merkle
Description
Merkle cavitation model.
Reference:
@verbatim
C. L. Merkle, J. Feng, and P. E. O. Buelow,
"Computational modeling of the dynamics of sheet cavitation",
in Proceedings Third International Symposium on Cavitation
Grenoble, France 1998.
@verbatim
SourceFiles
Merkle.C
\*--------------------------------------------------------------------*/
#ifndef Merkle_H
#define Merkle_H
#include "phaseChangeTwoPhaseMixture.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace phaseChangeTwoPhaseMixtures
{
/*--------------------------------------------------------------------*\
Class Merkle
\*--------------------------------------------------------------------*/
class Merkle
:
public phaseChangeTwoPhaseMixture
{
// Private data
dimensionedScalar UInf_;
dimensionedScalar tInf_;
dimensionedScalar Cc_;
dimensionedScalar Cv_;
dimensionedScalar p0_;
dimensionedScalar mcCoeff_;
dimensionedScalar mvCoeff_;
public:
//- Runtime type information
TypeName("Merkle");
// Constructors
//- construct from components
Merkle
(
const volVectorField& U,
const surfaceScalarField& phi,
const word& alpha1Name = "alpha1"
);
// Destructor
virtual ~Merkle()
{}
// Member Functions
//- Return the mass condensation and vaporisation rates as a
// coefficient to multiply (1 - alphal) for the condensation rate
// and a coefficient to multiply alphal for the vaporisation rate
virtual Pair<tmp<volScalarField> > mDotAlphal() const;
//- Return the mass condensation and vaporisation rates as an
// explicit term for the condensation rate and a coefficient to
// multiply (p - pSat) for the vaporisation rate
virtual Pair<tmp<volScalarField> > mDotP() const;
//- Correct the Merkle phaseChange model
virtual void correct();
//- Read the transportProperties dictionary and update
virtual bool read();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace phaseChangeTwoPhaseMixtures
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //