surfaceFilmModel.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 | Copyright (C) 2011-2017 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::surfaceFilmModel
26 
27 Description
28  Base class for surface film models
29 
30 SourceFiles
31  surfaceFilmModelI.H
32  surfaceFilmModel.C
33  surfaceFilmModelNew.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef surfaceFilmModel_H
38 #define surfaceFilmModel_H
39 
40 #include "singleLayerRegion.H"
41 
42 #include "dimensionedVector.H"
43 #include "runTimeSelectionTables.H"
44 #include "volFieldsFwd.H"
45 #include "DimensionedField.H"
46 #include "labelList.H"
47 #include "NamedEnum.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace regionModels
54 {
55 namespace surfaceFilmModels
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class surfaceFilmModel Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 class surfaceFilmModel
63 :
64  public singleLayerRegion
65 {
66 
67 private:
68 
69  // Private Member Functions
70 
71  //- Disallow default bitwise copy construct
73 
74  //- Disallow default bitwise assignment
75  void operator=(const surfaceFilmModel&);
76 
77 
78 protected:
79 
80  // Protected data
81 
82  //- Acceleration due to gravity [m/s2]
83  const dimensionedVector& g_;
84 
85 
86  // Protected member functions
87 
88  //- Read control parameters from dictionary
89  virtual bool read();
90 
91 
92 public:
93 
94  //- Runtime type information
95  TypeName("surfaceFilmModel");
96 
97  //- Reference temperature for enthalpy
98  static const dimensionedScalar Tref;
99 
100 
101  // Declare runtime constructor selection table
102 
104  (
105  autoPtr,
107  mesh,
108  (
109  const word& modelType,
110  const fvMesh& mesh,
111  const dimensionedVector& g,
112  const word& regionType
113  ),
114  (modelType, mesh, g, regionType)
115  );
116 
117  // Constructors
118 
119  //- Construct from type name, mesh and gravity vector
121  (
122  const word& modelType,
123  const fvMesh& mesh,
124  const dimensionedVector& g,
125  const word& regionType
126  );
127 
128 
129  // Selectors
130 
131  //- Return a reference to the selected surface film model
133  (
134  const fvMesh& mesh,
135  const dimensionedVector& g,
136  const word& regionType = "surfaceFilm"
137  );
138 
139 
140  //- Destructor
141  virtual ~surfaceFilmModel();
142 
143 
144  // Member Functions
145 
146  // Access
147 
148  //- Return the accleration due to gravity
149  inline const dimensionedVector& g() const;
150 
151  //- External hook to add sources to the film
152  virtual void addSources
153  (
154  const label patchi,
155  const label facei,
156  const scalar massSource,
157  const vector& momentumSource,
158  const scalar pressureSource,
159  const scalar energySource
160  ) = 0;
161 
162 
163  // Solution parameters
164 
165  //- Courant number evaluation
166  virtual scalar CourantNumber() const;
167 
168 
169  // Fields
170 
171  //- Return the film thickness [m]
172  virtual const volScalarField& delta() const = 0;
173 
174  //- Return the film coverage, 1 = covered, 0 = uncovered / []
175  virtual const volScalarField& alpha() const = 0;
176 
177  //- Return the film velocity [m/s]
178  virtual const volVectorField& U() const = 0;
179 
180  //- Return the film surface velocity [m/s]
181  virtual const volVectorField& Us() const = 0;
182 
183  //- Return the film wall velocity [m/s]
184  virtual const volVectorField& Uw() const = 0;
185 
186  //- Return the film density [kg/m3]
187  virtual const volScalarField& rho() const = 0;
188 
189  //- Return the film mean temperature [K]
190  virtual const volScalarField& T() const = 0;
191 
192  //- Return the film surface temperature [K]
193  virtual const volScalarField& Ts() const = 0;
194 
195  //- Return the film wall temperature [K]
196  virtual const volScalarField& Tw() const = 0;
197 
198  //- Return the film surface temperature [J/kg]
199  virtual const volScalarField& hs() const = 0;
200 
201  //- Return the film specific heat capacity [J/kg/K]
202  virtual const volScalarField& Cp() const = 0;
203 
204  //- Return the film thermal conductivity [W/m/K]
205  virtual const volScalarField& kappa() const = 0;
206 
207  //- Return the film surface tension [N/m]
208  virtual const volScalarField& sigma() const = 0;
209 
210 
211  // Transfer fields - to the primary region
212 
213  //- Return mass transfer source - Eulerian phase only
214  virtual tmp<volScalarField> primaryMassTrans() const = 0;
215 
216  //- Return the film mass available for transfer
217  virtual const volScalarField& cloudMassTrans() const = 0;
218 
219  //- Return the parcel diameters originating from film
220  virtual const volScalarField& cloudDiameterTrans() const = 0;
221 
222 
223  // Source fields
224 
225  // Mapped into primary region
226 
227  //- Return total mass source - Eulerian phase only
228  virtual tmp<volScalarField::Internal> Srho() const;
229 
230  //- Return mass source for specie i - Eulerian phase only
232  (
233  const label i
234  ) const;
235 
236  //- Return enthalpy source - Eulerian phase only
237  virtual tmp<volScalarField::Internal> Sh() const;
238 };
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 } // End namespace surfaceFilmModels
244 } // End namespace regionModels
245 } // End namespace Foam
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #include "surfaceFilmModelI.H"
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual const volScalarField & Tw() const =0
Return the film wall temperature [K].
const dimensionedVector & g() const
Return the accleration due to gravity.
TypeName("surfaceFilmModel")
Runtime type information.
virtual const volScalarField & kappa() const =0
Return the film thermal conductivity [W/m/K].
virtual const volScalarField & T() const =0
Return the film mean temperature [K].
declareRunTimeSelectionTable(autoPtr, surfaceFilmModel, mesh,(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType),(modelType, mesh, g, regionType))
const dimensionedVector & g_
Acceleration due to gravity [m/s2].
dynamicFvMesh & mesh
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
A class for handling words, derived from string.
Definition: word.H:59
virtual const volScalarField & Ts() const =0
Return the film surface temperature [K].
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
virtual const volScalarField & hs() const =0
Return the film surface temperature [J/kg].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
virtual bool read()
Read control parameters from dictionary.
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volScalarField & Cp() const =0
Return the film specific heat capacity [J/kg/K].
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource)=0
External hook to add sources to the film.
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
virtual scalar CourantNumber() const
Courant number evaluation.
static autoPtr< surfaceFilmModel > New(const fvMesh &mesh, const dimensionedVector &g, const word &regionType="surfaceFilm")
Return a reference to the selected surface film model.
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Macros to ease declaration of run-time selection tables.
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > primaryMassTrans() const =0
Return mass transfer source - Eulerian phase only.
virtual const volVectorField & Uw() const =0
Return the film wall velocity [m/s].
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
Namespace for OpenFOAM.
static const dimensionedScalar Tref
Reference temperature for enthalpy.
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].