surfaceFilmRegionModel.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) 2011-2018 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::surfaceFilmRegionModel
26 
27 Description
28  Base class for surface film models
29 
30 SourceFiles
31  surfaceFilmRegionModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef surfaceFilmRegionModel_H
36 #define surfaceFilmRegionModel_H
37 
38 #include "surfaceFilmModel.H"
39 #include "singleLayerRegion.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace regionModels
46 {
47 namespace surfaceFilmModels
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class surfaceFilmRegionModel Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 :
56  public surfaceFilmModel,
57  public singleLayerRegion
58 {
59  // Private Member Functions
60 
61  //- Disallow default bitwise copy construct
63 
64  //- Disallow default bitwise assignment
65  void operator=(const surfaceFilmRegionModel&);
66 
67 
68 protected:
69 
70  // Protected data
71 
72  //- Acceleration due to gravity [m/s2]
73  const dimensionedVector& g_;
74 
75 
76  // Protected member functions
77 
78  //- Read control parameters from dictionary
79  virtual bool read();
80 
81 
82 public:
83 
84  //- Runtime type information
85  TypeName("surfaceFilmRegionModel");
86 
87 
88  // Constructors
89 
90  //- Construct from type name, mesh and gravity vector
92  (
93  const word& modelType,
94  const fvMesh& mesh,
95  const dimensionedVector& g,
96  const word& regionType
97  );
98 
99 
100  //- Destructor
101  virtual ~surfaceFilmRegionModel();
102 
103 
104  // Member Functions
105 
106  // Access
107 
108  //- Return the accleration due to gravity
109  inline const dimensionedVector& g() const;
110 
111  //- External hook to add sources to the film
112  virtual void addSources
113  (
114  const label patchi,
115  const label facei,
116  const scalar massSource,
117  const vector& momentumSource,
118  const scalar pressureSource,
119  const scalar energySource
120  ) = 0;
121 
122 
123  // Fields
124 
125  //- Return the film thickness [m]
126  virtual const volScalarField& delta() const = 0;
127 
128  //- Return the film coverage, 1 = covered, 0 = uncovered / []
129  virtual const volScalarField& alpha() const = 0;
130 
131  //- Return the film velocity [m/s]
132  virtual const volVectorField& U() const = 0;
133 
134  //- Return the film surface velocity [m/s]
135  virtual const volVectorField& Us() const = 0;
136 
137  //- Return the film wall velocity [m/s]
138  virtual const volVectorField& Uw() const = 0;
139 
140  //- Return the film density [kg/m3]
141  virtual const volScalarField& rho() const = 0;
142 
143  //- Return the film mean temperature [K]
144  virtual const volScalarField& T() const = 0;
145 
146  //- Return the film surface temperature [K]
147  virtual const volScalarField& Ts() const = 0;
148 
149  //- Return the film wall temperature [K]
150  virtual const volScalarField& Tw() const = 0;
151 
152  //- Return the film surface temperature [J/kg]
153  virtual const volScalarField& hs() const = 0;
154 
155  //- Return the film specific heat capacity [J/kg/K]
156  virtual const volScalarField& Cp() const = 0;
157 
158  //- Return the film thermal conductivity [W/m/K]
159  virtual const volScalarField& kappa() const = 0;
160 
161  //- Return the film surface tension [N/m]
162  virtual const volScalarField& sigma() const = 0;
163 
164 
165  // Transfer fields - to the primary region
166 
167  //- Return mass transfer source - Eulerian phase only
168  virtual tmp<volScalarField> primaryMassTrans() const = 0;
169 
170  //- Return the film mass available for transfer
171  virtual const volScalarField& cloudMassTrans() const = 0;
172 
173  //- Return the parcel diameters originating from film
174  virtual const volScalarField& cloudDiameterTrans() const = 0;
175 
176 
177  // Evolution
178 
179  //- Main driver routing to evolve the region - calls other evolves
180  virtual void evolve();
181 };
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace surfaceFilmModels
187 } // End namespace regionModels
188 } // End namespace Foam
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #include "surfaceFilmRegionModelI.H"
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 #endif
197 
198 // ************************************************************************* //
virtual bool read()
Read control parameters from dictionary.
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.
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 volVectorField & Us() const =0
Return the film surface velocity [m/s].
Base class for surface film models.
virtual const volScalarField & hs() const =0
Return the film surface temperature [J/kg].
Base class for single layer region models.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
TypeName("surfaceFilmRegionModel")
Runtime type information.
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
virtual const volVectorField & Uw() const =0
Return the film wall velocity [m/s].
virtual tmp< volScalarField > primaryMassTrans() const =0
Return mass transfer source - Eulerian phase only.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
virtual const volScalarField & Tw() const =0
Return the film wall temperature [K].
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
const dimensionedVector & g() const
Return the accleration due to gravity.
virtual const volScalarField & Ts() const =0
Return the film surface temperature [K].
virtual const volScalarField & kappa() const =0
Return the film thermal conductivity [W/m/K].
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const volScalarField & Cp() const =0
Return the film specific heat capacity [J/kg/K].
const dimensionedVector & g_
Acceleration due to gravity [m/s2].
A class for managing temporary objects.
Definition: PtrList.H:53
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
Namespace for OpenFOAM.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
virtual const volScalarField & T() const =0
Return the film mean temperature [K].