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 | Copyright (C) 2015-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 "polynomial.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
35  defineTypeNameAndDebug(polynomial, 0);
36  addToRunTimeSelectionTable(saturationModel, polynomial, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 :
45  saturationModel(),
46  C_(dict.lookup("C<8>"))
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
51 
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
60 (
61  const volScalarField& T
62 ) const
63 {
65  return volScalarField::null();
66 }
67 
68 
71 (
72  const volScalarField& T
73 ) const
74 {
76  return volScalarField::null();
77 }
78 
79 
82 (
83  const volScalarField& T
84 ) const
85 {
87  return volScalarField::null();
88 }
89 
90 
93 (
94  const volScalarField& p
95 ) const
96 {
97  tmp<volScalarField> tTsat
98  (
99  new volScalarField
100  (
101  IOobject
102  (
103  "Tsat",
104  p.mesh().time().timeName(),
105  p.mesh(),
108  ),
109  p.mesh(),
111  )
112  );
113 
114  volScalarField& Tsat = tTsat.ref();
115 
116  forAll(Tsat,celli)
117  {
118  Tsat[celli] = C_.value(p[celli]);
119  }
120 
121  volScalarField::Boundary& TsatBf = Tsat.boundaryFieldRef();
122 
123  forAll(Tsat.boundaryField(), patchi)
124  {
125  scalarField& Tsatp = TsatBf[patchi];
126  const scalarField& pp = p.boundaryField()[patchi];
127 
128  forAll(Tsatp, facei)
129  {
130  Tsatp[facei] = C_.value(pp[facei]);
131  }
132  }
133 
134  return tTsat;
135 }
136 
137 
138 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual ~polynomial()
Destructor.
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Internal & ref()
Return a reference to the dimensioned internal field.
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
label patchi
polynomial(const dictionary &dict)
Construct from a dictionary.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
A class for managing temporary objects.
Definition: PtrList.H:54
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
Namespace for OpenFOAM.