KocamustafaogullariIshiiDepartureFrequency.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) 2019-2020 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 
31 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace wallBoilingModels
38 {
39 namespace departureFrequencyModels
40 {
41  defineTypeNameAndDebug(KocamustafaogullariIshiiDepartureFrequency, 0);
43  (
44  departureFrequencyModel,
45  KocamustafaogullariIshiiDepartureFrequency,
46  dictionary
47  );
48 }
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
58 :
59  departureFrequencyModel(),
60  Cf_(dict.lookupOrDefault<scalar>("Cf", 1.18))
61 {}
62 
63 
67 (
68  const KocamustafaogullariIshiiDepartureFrequency& model
69 )
70 :
72  Cf_(model.Cf_)
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
77 
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
89 (
90  const phaseModel& liquid,
91  const phaseModel& vapor,
92  const label patchi,
93  const scalarField& Tl,
94  const scalarField& Tsatw,
95  const scalarField& L,
96  const scalarField& dDep
97 ) const
98 {
99  // Gravitational acceleration
101  liquid.mesh().lookupObject<uniformDimensionedVectorField>("g");
102 
103  const scalarField rhoLiquid(liquid.thermo().rho(patchi));
104  const scalarField rhoVapor(min(vapor.thermo().rho(patchi), rhoLiquid));
105 
106  const tmp<volScalarField> tsigma
107  (
108  liquid.fluid().sigma(phasePairKey(liquid.name(), vapor.name()))
109  );
110 
111  const volScalarField& sigma = tsigma();
112  const fvPatchScalarField& sigmaw = sigma.boundaryField()[patchi];
113 
114  return (Cf_/dDep)*pow025
115  (
116  sigmaw*mag(g.value())*(rhoLiquid - rhoVapor)/sqr(rhoLiquid)
117  );
118 }
119 
120 
123 {
125  writeKeyword(os, "Cf") << Cf_ << token::END_STATEMENT << nl;
126 }
127 
128 
129 // ************************************************************************* //
dictionary dict
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
UniformDimensionedField< vector > uniformDimensionedVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow025(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual void write(Ostream &os) const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label patchi
virtual tmp< scalarField > fDeparture(const phaseModel &liquid, const phaseModel &vapor, const label patchi, const scalarField &Tl, const scalarField &Tsatw, const scalarField &L, const scalarField &dDep) const
Calculate and return the bubble departure frequency.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
KocamustafaogullariIshiiDepartureFrequency(const dictionary &dict)
Construct from a dictionary.
Namespace for OpenFOAM.