liquidFilmThermo.H
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) 2013-2019 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 Class
25  Foam::regionModels::surfaceFilmModels::liquidFilmThermo
26 
27 Description
28  Liquid thermo model
29 
30 SourceFiles
31  liquidFilmThermo.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef liquidFilmThermo_H
36 #define liquidFilmThermo_H
37 
38 #include "filmThermoModel.H"
39 #include "liquidProperties.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace regionModels
46 {
47 namespace surfaceFilmModels
48 {
49 
50 // Forward class declarations
51 class thermoSingleLayer;
52 
53 /*---------------------------------------------------------------------------*\
54  Class liquidFilmThermo Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class liquidFilmThermo
58 :
59  public filmThermoModel
60 {
61 protected:
62 
63  // Protected data
64 
65  //- Liquid name
66  word name_;
67 
68  //- Pointer to the liquid properties
70 
71  //- Flag to indicate that model owns the liquid object
72  bool ownLiquid_;
73 
74  //- Flag to indicate that reference values of p and T should be used
76 
77  //- Reference pressure [pa]
78  scalar pRef_;
79 
80  //- Reference temperature [K]
81  scalar TRef_;
82 
83 
84  // Protected member functions
85 
86  //- Return a reference to a thermo film
87  const thermoSingleLayer& thermoFilm() const;
88 
89  //- Initialise the liquid pointer
90  void initLiquid(const dictionary& dict);
91 
92 
93 public:
94 
95  //- Runtime type information
96  TypeName("liquid");
97 
98 
99  // Constructors
100 
101  //- Construct from surface film model and dictionary
103  (
105  const dictionary& dict
106  );
107 
108  //- Disallow default bitwise copy construction
109  liquidFilmThermo(const liquidFilmThermo&) = delete;
110 
111 
112  //- Destructor
113  virtual ~liquidFilmThermo();
114 
115 
116  // Member Functions
117 
118  //- Return the specie name
119  virtual const word& name() const;
120 
121 
122  // Elemental access
123 
124  //- Return density [kg/m^3]
125  virtual scalar rho(const scalar p, const scalar T) const;
126 
127  //- Return dynamic viscosity [Pa.s]
128  virtual scalar mu(const scalar p, const scalar T) const;
129 
130  //- Return surface tension [kg/s^2]
131  virtual scalar sigma(const scalar p, const scalar T) const;
132 
133  //- Return specific heat capacity [J/kg/K]
134  virtual scalar Cp(const scalar p, const scalar T) const;
135 
136  //- Return thermal conductivity [W/m/K]
137  virtual scalar kappa(const scalar p, const scalar T) const;
138 
139  //- Return diffusivity [m^2/s]
140  virtual scalar D(const scalar p, const scalar T) const;
141 
142  //- Return latent heat [J/kg]
143  virtual scalar hl(const scalar p, const scalar T) const;
144 
145  //- Return vapour pressure [Pa]
146  virtual scalar pv(const scalar p, const scalar T) const;
147 
148  //- Return molecular weight [kg/kmol]
149  virtual scalar W() const;
150 
151  //- Return boiling temperature [K]
152  virtual scalar Tb(const scalar p) const;
153 
154 
155  // Field access
156 
157  //- Return density [kg/m^3]
158  virtual tmp<volScalarField> rho() const;
159 
160  //- Return dynamic viscosity [Pa.s]
161  virtual tmp<volScalarField> mu() const;
162 
163  //- Return surface tension [kg/s^2]
164  virtual tmp<volScalarField> sigma() const;
165 
166  //- Return specific heat capacity [J/kg/K]
167  virtual tmp<volScalarField> Cp() const;
168 
169  //- Return thermal conductivity [W/m/K]
170  virtual tmp<volScalarField> kappa() const;
171 
172 
173  // Member Operators
174 
175  //- Disallow default bitwise assignment
176  void operator=(const liquidFilmThermo&) = delete;
177 };
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace surfaceFilmModels
183 } // End namespace regionModels
184 } // End namespace Foam
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #endif
189 
190 // ************************************************************************* //
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m^2/s].
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
const liquidProperties * liquidPtr_
Pointer to the liquid properties.
virtual scalar W() const
Return molecular weight [kg/kmol].
bool useReferenceValues_
Flag to indicate that reference values of p and T should be used.
void initLiquid(const dictionary &dict)
Initialise the liquid pointer.
A class for handling words, derived from string.
Definition: word.H:59
const thermoSingleLayer & thermoFilm() const
Return a reference to a thermo film.
virtual tmp< volScalarField > rho() const
Return density [kg/m^3].
The thermophysical properties of a liquid.
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s^2].
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
virtual const word & name() const
Return the specie name.
TypeName("liquid")
Runtime type information.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
void operator=(const liquidFilmThermo &)=delete
Disallow default bitwise assignment.
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
liquidFilmThermo(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model and dictionary.
bool ownLiquid_
Flag to indicate that model owns the liquid object.