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