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