Merkle.C
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-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "Merkle.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace phaseChangeTwoPhaseMixtures
37{
38 defineTypeNameAndDebug(Merkle, 0);
39 addToRunTimeSelectionTable(phaseChangeTwoPhaseMixture, Merkle, components);
40}
41}
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46(
47 const volVectorField& U,
48 const surfaceScalarField& phi
49)
50:
51 phaseChangeTwoPhaseMixture(typeName, U, phi),
52
53 UInf_("UInf", dimVelocity, phaseChangeTwoPhaseMixtureCoeffs_),
54 tInf_("tInf", dimTime, phaseChangeTwoPhaseMixtureCoeffs_),
55 Cc_("Cc", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
56 Cv_("Cv", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
57
58 p0_(pSat().dimensions(), Zero),
59
60 mcCoeff_(Cc_/(0.5*sqr(UInf_)*tInf_)),
61 mvCoeff_(Cv_*rho1()/(0.5*sqr(UInf_)*tInf_*rho2()))
62{
63 correct();
64}
65
66
67// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68
71{
73
74 return Pair<tmp<volScalarField>>
75 (
76 mcCoeff_*max(p - pSat(), p0_),
77 mvCoeff_*min(p - pSat(), p0_)
78 );
79}
80
83{
84 const volScalarField& p = alpha1_.db().lookupObject<volScalarField>("p");
85 volScalarField limitedAlpha1(min(max(alpha1_, scalar(0)), scalar(1)));
86
87 return Pair<tmp<volScalarField>>
88 (
89 mcCoeff_*(1.0 - limitedAlpha1)*pos0(p - pSat()),
90 (-mvCoeff_)*limitedAlpha1*neg(p - pSat())
91 );
92}
93
94
96{}
97
98
100{
102 {
103 phaseChangeTwoPhaseMixtureCoeffs_ = optionalSubDict(type() + "Coeffs");
104
105 phaseChangeTwoPhaseMixtureCoeffs_.readEntry("UInf", UInf_);
106 phaseChangeTwoPhaseMixtureCoeffs_.readEntry("tInf", tInf_);
107 phaseChangeTwoPhaseMixtureCoeffs_.readEntry("Cc", Cc_);
108 phaseChangeTwoPhaseMixtureCoeffs_.readEntry("Cv", Cv_);
109
110 mcCoeff_ = Cc_/(0.5*sqr(UInf_)*tInf_);
111 mvCoeff_ = Cv_*rho1()/(0.5*sqr(UInf_)*tInf_*rho2());
112
113 return true;
114 }
115
116 return false;
117}
118
119
120// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
volScalarField & rho2
volScalarField & rho1
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
const Type & lookupObject(const word &name, const bool recursive=false) const
virtual bool read()=0
Read the transportProperties dictionary and update.
const dimensionedScalar & pSat() const
Return const-access to the saturation vapour pressure.
Merkle cavitation model.
Definition: Merkle.H:64
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
virtual void correct()
Correct the Merkle phaseChange model.
virtual bool read()
Read the transportProperties dictionary and update.
volScalarField alpha1_
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
thermo correct()
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimVelocity
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensionedScalar neg(const dimensionedScalar &ds)