constantSaturationConditions.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 
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 
50 (
51  const dictionary& dict,
52  const phaseInterface& interface
53 )
54 :
55  saturationModel(dict, interface),
56  pSat_("pSat", dimPressure, dict),
57  Tsat_("Tsat", dimTemperature, dict)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
72 (
73  const volScalarField& T
74 ) const
75 {
76  return volScalarField::New
77  (
78  "pSat",
79  T.mesh(),
80  pSat_
81  );
82 }
83 
84 
87 (
88  const volScalarField& T
89 ) const
90 {
91  return volScalarField::New
92  (
93  "pSatPrime",
94  T.mesh(),
96  );
97 }
98 
99 
102 (
103  const volScalarField& T
104 ) const
105 {
106  return volScalarField::New
107  (
108  "lnPSat",
109  T.mesh(),
110  dimensionedScalar(dimless, log(pSat_.value()))
111  );
112 }
113 
114 
117 (
118  const volScalarField& p
119 ) const
120 {
121  return volScalarField::New
122  (
123  "Tsat",
124  p.mesh(),
125  Tsat_
126  );
127 }
128 
129 
130 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
const dimensionSet dimPressure
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivative w.r.t. temperature.
constantSaturationConditions(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimTemperature
Namespace for OpenFOAM.