constantSaturationConditions.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 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 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
35  defineTypeNameAndDebug(constantSaturationConditions, 0);
37  (
38  saturationModel,
39  constantSaturationConditions,
40  dictionary
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 constantSaturationConditions(const dictionary& dict)
50 :
51  saturationModel(),
52  pSat_("pSat", dimPressure, dict),
53  Tsat_("Tsat", dimTemperature, dict)
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
68 (
69  const volScalarField& T
70 ) const
71 {
72  return tmp<volScalarField>
73  (
74  new volScalarField
75  (
76  IOobject
77  (
78  "pSat",
79  T.mesh().time().timeName(),
80  T.mesh(),
83  false
84  ),
85  T.mesh(),
86  pSat_
87  )
88  );
89 }
90 
91 
94 (
95  const volScalarField& T
96 ) const
97 {
98  return tmp<volScalarField>
99  (
100  new volScalarField
101  (
102  IOobject
103  (
104  "pSatPrime",
105  T.mesh().time().timeName(),
106  T.mesh(),
109  false
110  ),
111  T.mesh(),
113  )
114  );
115 }
116 
117 
120 (
121  const volScalarField& T
122 ) const
123 {
124  return tmp<volScalarField>
125  (
126  new volScalarField
127  (
128  IOobject
129  (
130  "lnPSat",
131  T.mesh().time().timeName(),
132  T.mesh(),
135  false
136  ),
137  T.mesh(),
138  dimensionedScalar("lnPSat", dimless, log(pSat_.value()))
139  )
140  );
141 }
142 
143 
146 (
147  const volScalarField& p
148 ) const
149 {
150  return tmp<volScalarField>
151  (
152  new volScalarField
153  (
154  IOobject
155  (
156  "Tsat",
157  p.mesh().time().timeName(),
158  p.mesh(),
161  false
162  ),
163  p.mesh(),
164  Tsat_
165  )
166  );
167 }
168 
169 
170 // ************************************************************************* //
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
constantSaturationConditions(const dictionary &dict)
Construct from a dictionary.
const dimensionSet dimPressure
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:54
Namespace for OpenFOAM.