constantFilmThermo.H
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) 2013-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 Class
27  Foam::regionModels::surfaceFilmModels::constantFilmThermo
28 
29 Description
30  Constant thermo model
31 
32 SourceFiles
33  constantFilmThermo.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef constantFilmThermo_H
38 #define constantFilmThermo_H
39 
40 #include "filmThermoModel.H"
41 #include "dimensionSet.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 namespace regionModels
48 {
49 namespace surfaceFilmModels
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class constantFilmThermo Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 :
58  public filmThermoModel
59 {
60 public:
61 
62  struct thermoData
63  {
65  scalar value_;
66  bool set_;
67 
68  thermoData()
69  :
70  name_("unknown"),
71  value_(0.0),
72  set_(false)
73  {}
74 
75  thermoData(const word& n)
76  :
77  name_(n),
78  value_(0.0),
79  set_(false)
80  {}
81  };
82 
83 
84 private:
85 
86  // Private data
87 
88  //- Specie name
89  word name_;
90 
91  //- Density [kg/m3]
92  mutable thermoData rho0_;
93 
94  //- Dynamic viscosity [Pa.s]
95  mutable thermoData mu0_;
96 
97  //- Surface tension [kg/s2]
98  mutable thermoData sigma0_;
99 
100  //- Specific heat capacity [J/kg/K]
101  mutable thermoData Cp0_;
102 
103  //- Thermal conductivity [W/m/K]
104  mutable thermoData kappa0_;
105 
106  //- Diffusivity [m2/s]
107  mutable thermoData D0_;
108 
109  //- Latent heat [J/kg]
110  mutable thermoData hl0_;
111 
112  //- Vapour pressure [Pa]
113  mutable thermoData pv0_;
114 
115  //- Molecular weight [kg/kmol]
116  mutable thermoData W0_;
117 
118  //- Boiling temperature [K]
119  mutable thermoData Tb0_;
120 
121 
122  // Private member functions
123 
124  //- Initialise thermoData object
125  void init(thermoData& td);
126 
127  //- No copy construct
128  constantFilmThermo(const constantFilmThermo&) = delete;
129 
130  //- No copy assignment
131  void operator=(const constantFilmThermo&) = delete;
132 
133 
134 public:
135 
136  //- Runtime type information
137  TypeName("constant");
138 
139 
140  // Constructors
141 
142  //- Construct from surface film model and dictionary
144  (
146  const dictionary& dict
147  );
148 
149 
150  //- Destructor
151  virtual ~constantFilmThermo();
152 
153 
154  // Member Functions
155 
156  //- Return the specie name
157  virtual const word& name() const;
158 
159 
160  // Elemental access
161 
162  //- Return density [kg/m3]
163  virtual scalar rho(const scalar p, const scalar T) const;
164 
165  //- Return dynamic viscosity [Pa.s]
166  virtual scalar mu(const scalar p, const scalar T) const;
167 
168  //- Return surface tension [kg/s2]
169  virtual scalar sigma(const scalar p, const scalar T) const;
170 
171  //- Return specific heat capacity [J/kg/K]
172  virtual scalar Cp(const scalar p, const scalar T) const;
173 
174  //- Return thermal conductivity [W/m/K]
175  virtual scalar kappa(const scalar p, const scalar T) const;
176 
177  //- Return diffusivity [m2/s]
178  virtual scalar D(const scalar p, const scalar T) const;
179 
180  //- Return latent heat [J/kg]
181  virtual scalar hl(const scalar p, const scalar T) const;
182 
183  //- Return vapour pressure [Pa]
184  virtual scalar pv(const scalar p, const scalar T) const;
185 
186  //- Return molecular weight [kg/kmol]
187  virtual scalar W() const;
188 
189  //- Return boiling temperature [K]
190  virtual scalar Tb(const scalar p) const;
191 
192 
193  // Field access
194 
195  //- Return density [kg/m3]
196  virtual tmp<volScalarField> rho() const;
197 
198  //- Return dynamic viscosity [Pa.s]
199  virtual tmp<volScalarField> mu() const;
200 
201  //- Return surface tension [kg/s2]
202  virtual tmp<volScalarField> sigma() const;
203 
204  //- Return specific heat capacity [J/kg/K]
205  virtual tmp<volScalarField> Cp() const;
206 
207  //- Return thermal conductivity [W/m/K]
208  virtual tmp<volScalarField> kappa() const;
209 };
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace surfaceFilmModels
215 } // End namespace regionModels
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::name_
word name_
Definition: constantFilmThermo.H:63
p
volScalarField & p
Definition: createFieldRefs.H:8
filmThermoModel.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData
Definition: constantFilmThermo.H:61
Foam::regionModels::surfaceFilmModels::constantFilmThermo::mu
virtual tmp< volScalarField > mu() const
Return dynamic viscosity [Pa.s].
Definition: constantFilmThermo.C:289
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::regionModels::surfaceFilmModels::constantFilmThermo::rho
virtual tmp< volScalarField > rho() const
Return density [kg/m3].
Definition: constantFilmThermo.C:262
Foam::regionModels::surfaceFilmModels::constantFilmThermo::Cp
virtual tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: constantFilmThermo.C:343
Foam::regionModels::surfaceFilmModels::filmSubModelBase::film
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Definition: filmSubModelBaseI.H:39
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::thermoData
thermoData()
Definition: constantFilmThermo.H:67
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::value_
scalar value_
Definition: constantFilmThermo.H:64
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::thermoData
thermoData(const word &n)
Definition: constantFilmThermo.H:74
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::subModelBase::dict
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:113
Foam::regionModels::surfaceFilmModels::constantFilmThermo::thermoData::set_
bool set_
Definition: constantFilmThermo.H:65
Foam::regionModels::surfaceFilmModels::constantFilmThermo::pv
virtual scalar pv(const scalar p, const scalar T) const
Return vapour pressure [Pa].
Definition: constantFilmThermo.C:223
Foam::regionModels::surfaceFilmModels::constantFilmThermo::D
virtual scalar D(const scalar p, const scalar T) const
Return diffusivity [m2/s].
Definition: constantFilmThermo.C:191
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
dimensionSet.H
Foam::regionModels::surfaceFilmModels::constantFilmThermo::W
virtual scalar W() const
Return molecular weight [kg/kmol].
Definition: constantFilmThermo.C:238
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::surfaceFilmModels::filmThermoModel
Base class for film thermo models.
Definition: filmThermoModel.H:56
Foam::regionModels::surfaceFilmModels::constantFilmThermo::sigma
virtual tmp< volScalarField > sigma() const
Return surface tension [kg/s2].
Definition: constantFilmThermo.C:316
Foam::regionModels::surfaceFilmModels::constantFilmThermo::hl
virtual scalar hl(const scalar p, const scalar T) const
Return latent heat [J/kg].
Definition: constantFilmThermo.C:207
Foam::regionModels::surfaceFilmModels::constantFilmThermo::Tb
virtual scalar Tb(const scalar p) const
Return boiling temperature [K].
Definition: constantFilmThermo.C:250
Foam::regionModels::surfaceFilmModels::constantFilmThermo
Constant thermo model.
Definition: constantFilmThermo.H:55
Foam::regionModels::surfaceFilmModels::constantFilmThermo::kappa
virtual tmp< volScalarField > kappa() const
Return thermal conductivity [W/m/K].
Definition: constantFilmThermo.C:370
Foam::regionModels::surfaceFilmModels::constantFilmThermo::TypeName
TypeName("constant")
Runtime type information.
Foam::regionModels::surfaceFilmModels::constantFilmThermo::~constantFilmThermo
virtual ~constantFilmThermo()
Destructor.
Definition: constantFilmThermo.C:98
Foam::regionModels::surfaceFilmModels::constantFilmThermo::name
virtual const word & name() const
Return the specie name.
Definition: constantFilmThermo.C:104