Kunz.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "Kunz.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace phaseChangeTwoPhaseMixtures
34 {
35  defineTypeNameAndDebug(Kunz, 0);
36  addToRunTimeSelectionTable(phaseChangeTwoPhaseMixture, Kunz, components);
37 }
38 }
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const volVectorField& U,
45  const surfaceScalarField& phi
46 )
47 :
48  phaseChangeTwoPhaseMixture(typeName, U, phi),
49 
50  UInf_("UInf", dimVelocity, phaseChangeTwoPhaseMixtureCoeffs_),
51  tInf_("tInf", dimTime, phaseChangeTwoPhaseMixtureCoeffs_),
52  Cc_("Cc", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
53  Cv_("Cv", dimless, phaseChangeTwoPhaseMixtureCoeffs_),
54 
55  p0_("0", pSat().dimensions(), 0.0),
56 
57  mcCoeff_(Cc_*rho2()/tInf_),
58  mvCoeff_(Cv_*rho2()/(0.5*rho1()*sqr(UInf_)*tInf_))
59 {
60  correct();
61 }
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
68 {
69  const volScalarField& p = alpha1_.db().lookupObject<volScalarField>("p");
70  volScalarField limitedAlpha1(min(max(alpha1_, scalar(0)), scalar(1)));
71 
72  return Pair<tmp<volScalarField>>
73  (
74  mcCoeff_*sqr(limitedAlpha1)
75  *max(p - pSat(), p0_)/max(p - pSat(), 0.01*pSat()),
76 
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_*sqr(limitedAlpha1)*(1.0 - limitedAlpha1)
90  *pos0(p - pSat())/max(p - pSat(), 0.01*pSat()),
91 
92  (-mvCoeff_)*limitedAlpha1*neg(p - pSat())
93  );
94 }
95 
96 
98 {}
99 
100 
102 {
104  {
105  phaseChangeTwoPhaseMixtureCoeffs_ = optionalSubDict(type() + "Coeffs");
106 
107  phaseChangeTwoPhaseMixtureCoeffs_.lookup("UInf") >> UInf_;
108  phaseChangeTwoPhaseMixtureCoeffs_.lookup("tInf") >> tInf_;
109  phaseChangeTwoPhaseMixtureCoeffs_.lookup("Cc") >> Cc_;
110  phaseChangeTwoPhaseMixtureCoeffs_.lookup("Cv") >> Cv_;
111 
112  mcCoeff_ = Cc_*rho2()/tInf_;
113  mvCoeff_ = Cv_*rho2()/(0.5*rho1()*sqr(UInf_)*tInf_);
114 
115  return true;
116  }
117  else
118  {
119  return false;
120  }
121 }
122 
123 
124 // ************************************************************************* //
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
virtual bool read()=0
Read the transportProperties dictionary and update.
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
surfaceScalarField & phi
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
dimensionedScalar pos0(const dimensionedScalar &ds)
Kunz(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionedScalar & rho2
Definition: createFields.H:40
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void correct()
Correct the Kunz phaseChange model.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
const dimensionedScalar & rho1
Definition: createFields.H:39
virtual bool read()
Read the transportProperties dictionary and update.
Namespace for OpenFOAM.
const dimensionSet dimVelocity