Loading...
Searching...
No Matches
polynomial.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) 2015-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "polynomial.H"
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace saturationModels
36{
39}
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
46(
47 const dictionary& dict,
48 const objectRegistry& db
49)
50:
51 saturationModel(db),
52 C_(dict.lookup("C<8>"))
53{}
54
55
56// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57
59{}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
66(
67 const volScalarField& T
68) const
77(
78 const volScalarField& T
79) const
88(
89 const volScalarField& T
90) const
93 return volScalarField::null();
94}
95
96
99(
100 const volScalarField& p
101) const
102{
104 (
106 (
107 "Tsat",
108 p.mesh(),
110 )
111 );
112
113 volScalarField& Tsat = tTsat.ref();
114
115 forAll(Tsat, celli)
116 {
117 Tsat[celli] = C_.value(p[celli]);
118 }
119
121
122 forAll(Tsat.boundaryField(), patchi)
123 {
124 scalarField& Tsatp = TsatBf[patchi];
125 const scalarField& pp = p.boundaryField()[patchi];
126
127 forAll(Tsatp, facei)
128 {
129 Tsatp[facei] = C_.value(pp[facei]);
130 }
131 }
132
133 return tTsat;
134}
135
136
137// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Registry of regIOobjects.
Lookup type of boundary radiation properties.
Definition lookup.H:60
Polynomial equation for the saturation vapour temperature in terms of the vapour pressure (in Pa).
Definition polynomial.H:62
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
Definition polynomial.C:81
virtual ~polynomial()
Destructor.
Definition polynomial.C:51
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
Definition polynomial.C:92
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
Definition polynomial.C:70
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
Definition polynomial.C:59
polynomial(const dictionary &dict, const objectRegistry &db)
Construct from a dictionary.
Definition polynomial.C:39
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299