exponential.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) 2021 OpenCFD Ltd
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "exponential.H"
31 #include "phasePairKey.H"
32 #include "phaseSystem.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace wallBoilingModels
39 {
40 namespace nucleateFluxModels
41 {
44  (
48  );
49 }
50 }
51 }
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 Foam::wallBoilingModels::nucleateFluxModels::exponential::exponential
56 (
57  const dictionary& dict
58 )
59 :
61  a_(dict.getOrDefault<scalar>("a", 6309)),
62  b_(dict.getOrDefault<scalar>("b", 2.52))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
70 (
71  const phaseModel& liquid,
72  const phaseModel& vapor,
73  const label patchi,
74  const scalarField& Tl,
75  const scalarField& Tsatw,
76  const scalarField& L
77 ) const
78 {
79  const fvPatchScalarField& Tw =
80  liquid.thermo().T().boundaryField()[patchi];
81 
82  return a_*pow(max((Tw-Tsatw), scalar(0)), b_);
83 }
84 
85 
87 (
88  Ostream& os
89 ) const
90 {
92  os.writeEntry("a", a_);
93  os.writeEntry("b", b_);
94 }
95 
96 
97 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
L
const vector L(dict.get< vector >("L"))
Foam::wallBoilingModels::nucleateFluxModel
Base class for nucleation flux models.
Definition: nucleateFluxModel.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::wallBoilingModels::nucleateFluxModels::exponential::qNucleate
virtual tmp< scalarField > qNucleate(const phaseModel &liquid, const phaseModel &vapor, const label patchi, const scalarField &Tl, const scalarField &Tsatw, const scalarField &L) const
Calculate and return the nucleation-site density.
Definition: exponential.C:70
Foam::wallBoilingModels::nucleateFluxModels::defineTypeNameAndDebug
defineTypeNameAndDebug(exponential, 0)
Foam::liquid
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:54
Foam::Field< scalar >
Foam::wallBoilingModels::nucleateFluxModels::exponential::write
virtual void write(Ostream &os) const
Write.
Definition: exponential.C:87
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
uniformDimensionedFields.H
Foam::wallBoilingModels::nucleateFluxModels::exponential
Nucleate flux sub-cooling correlation.
Definition: exponential.H:115
Foam::wallBoilingModels::nucleateFluxModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(nucleateFluxModel, exponential, dictionary)
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
exponential.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56