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 | Copyright (C) 2011-2016 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 "Merkle.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace phaseChangeTwoPhaseMixtures
34 {
35  defineTypeNameAndDebug(Merkle, 0);
36  addToRunTimeSelectionTable(phaseChangeTwoPhaseMixture, Merkle, 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_(phaseChangeTwoPhaseMixtureCoeffs_.lookup("UInf")),
51  tInf_(phaseChangeTwoPhaseMixtureCoeffs_.lookup("tInf")),
52  Cc_(phaseChangeTwoPhaseMixtureCoeffs_.lookup("Cc")),
53  Cv_(phaseChangeTwoPhaseMixtureCoeffs_.lookup("Cv")),
54 
55  p0_("0", pSat().dimensions(), 0.0),
56 
57  mcCoeff_(Cc_/(0.5*sqr(UInf_)*tInf_)),
58  mvCoeff_(Cv_*rho1()/(0.5*sqr(UInf_)*tInf_*rho2()))
59 {
60  correct();
61 }
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
68 {
70 
71  return Pair<tmp<volScalarField>>
72  (
73  mcCoeff_*max(p - pSat(), p0_),
74  mvCoeff_*min(p - pSat(), p0_)
75  );
76 }
77 
80 {
82  volScalarField limitedAlpha1(min(max(alpha1_, scalar(0)), scalar(1)));
83 
84  return Pair<tmp<volScalarField>>
85  (
86  mcCoeff_*(1.0 - limitedAlpha1)*pos(p - pSat()),
87  (-mvCoeff_)*limitedAlpha1*neg(p - pSat())
88  );
89 }
90 
91 
93 {}
94 
95 
97 {
99  {
101 
102  phaseChangeTwoPhaseMixtureCoeffs_.lookup("UInf") >> UInf_;
103  phaseChangeTwoPhaseMixtureCoeffs_.lookup("tInf") >> tInf_;
106 
107  mcCoeff_ = Cc_/(0.5*sqr(UInf_)*tInf_);
108  mvCoeff_ = Cv_*rho1()/(0.5*sqr(UInf_)*tInf_*rho2());
109 
110  return true;
111  }
112  else
113  {
114  return false;
115  }
116 }
117 
118 
119 // ************************************************************************* //
surfaceScalarField & phi
U
Definition: pEqn.H:83
virtual bool read()=0
Read the transportProperties dictionary and update.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Merkle(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
rho1
Definition: pEqn.H:114
const dimensionedScalar & rho2() const
Return const-access to phase2 density.
const dimensionedScalar & pSat() const
Return const-access to the saturation vapour pressure.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
dimensionedScalar pos(const dimensionedScalar &ds)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
virtual void correct()
Correct the Merkle phaseChange model.
stressControl lookup("compactNormalStress") >> compactNormalStress
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
rho2
Definition: pEqn.H:115
volScalarField alpha1_
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read()
Read the transportProperties dictionary and update.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
const dimensionedScalar & rho1() const
Return const-access to phase1 density.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void correct()
Correct the Kunz phaseChange model.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451