polynomialTemperature.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) 2015-2022 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 "polynomialTemperature.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
37  (
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict
51 )
52 :
54  C_(dict.lookup("C<8>"))
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
68 (
70 ) const
71 {
73  (
75  (
76  "Tsat",
77  p.mesh(),
79  )
80  );
81 
82  volScalarField::Internal& Tsat = tTsat.ref();
83 
84  forAll(Tsat, celli)
85  {
86  Tsat[celli] = C_.value(p[celli]);
87  }
88 
89  return tTsat;
90 }
91 
92 
95 (
96  const volScalarField& p
97 ) const
98 {
100  (
102  (
103  "Tsat",
104  p.mesh(),
106  )
107  );
108 
109  volScalarField& Tsat = tTsat.ref();
110 
111  forAll(Tsat, celli)
112  {
113  Tsat[celli] = C_.value(p[celli]);
114  }
115 
117 
118  forAll(Tsat.boundaryField(), patchi)
119  {
120  scalarField& Tsatp = TsatBf[patchi];
121  const scalarField& pp = p.boundaryField()[patchi];
122 
123  forAll(Tsatp, facei)
124  {
125  Tsatp[facei] = C_.value(pp[facei]);
126  }
127  }
128 
129  return tTsat;
130 }
131 
132 
133 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Polynomial equation for the saturation vapour temperature in terms of the vapour pressure (in Pa).
virtual tmp< volScalarField::Internal > Tsat(const volScalarField::Internal &p) const
Saturation temperature for volScalarField::Internal.
polynomialTemperature(const dictionary &dict)
Construct from a dictionary.
Model to describe the dependence of saturation pressure on temperature, and vice versa.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
addToRunTimeSelectionTable(saturationPressureModel, Antoine, dictionary)
defineTypeNameAndDebug(Antoine, 0)
Namespace for OpenFOAM.
const dimensionSet dimTemperature
dictionary dict
volScalarField & p