injectionModel.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) 2021 OpenCFD Ltd.
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::areaSurfaceFilmModels::injectionModel
28 
29 Description
30  Base class for film injection models, handling mass transfer from the
31  film.
32 
33 SourceFiles
34  injectionModel.C
35  injectionModelNew.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef injectionModel_H
40 #define injectionModel_H
41 
42 #include "filmSubModelBase.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace regionModels
49 {
50 namespace areaSurfaceFilmModels
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class injectionModel Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class injectionModel
58 :
59  public filmSubModelBase
60 {
61  // Private Data
62 
63  //- Injected mass
64  scalar injectedMass_;
65 
66 
67  // Private Member Functions
68 
69  //- No copy construct
70  injectionModel(const injectionModel&) = delete;
71 
72  //- No copy assignment
73  void operator=(const injectionModel&) = delete;
74 
75 
76 protected:
77 
78  // Protected Member Functions
79 
80  //- Add to injected mass
81  void addToInjectedMass(const scalar dMass);
82 
83  //- Correct
84  void correct();
85 
86 
87 public:
88 
89  //- Runtime type information
90  TypeName("injectionModel");
91 
92 
93  // Declare runtime constructor selection table
94 
96  (
97  autoPtr,
99  dictionary,
100  (
102  const dictionary& dict
103  ),
104  (film, dict)
105  );
106 
107 
108  // Constructors
109 
110  //- Construct null
112 
113  //- Construct from type name, dictionary and surface film model
115  (
116  const word& modelType,
118  const dictionary& dict
119  );
120 
121 
122  // Selectors
123 
124  //- Return a reference to the selected injection model
126  (
128  const dictionary& dict,
129  const word& mdoelType
130  );
131 
132 
133  //- Destructor
134  virtual ~injectionModel();
135 
136 
137  // Member Functions
138 
139  //- Correct
140  virtual void correct
141  (
142  scalarField& availableMass,
143  scalarField& massToInject,
144  scalarField& diameterToInject
145  ) = 0;
146 
147  //- Return the total mass injected
148  virtual scalar injectedMassTotal() const;
149 
150  //- Accumulate the total mass injected for the patches into the
151  //- scalarField provided
152  virtual void patchInjectedMassTotals(scalar& patchMasses) const
153  {}
154 };
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace areaSurfaceFilmModels
160 } // End namespace regionModels
161 } // End namespace Foam
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 #endif
166 
167 // ************************************************************************* //
void addToInjectedMass(const scalar dMass)
Add to injected mass.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
TypeName("injectionModel")
Runtime type information.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:104
Base class for film injection models, handling mass transfer from the film.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const word & modelType() const
Return const access to the sub-model type.
Definition: subModelBase.C:116
virtual void patchInjectedMassTotals(scalar &patchMasses) const
Accumulate the total mass injected for the patches into the scalarField provided. ...
declareRunTimeSelectionTable(autoPtr, injectionModel, dictionary,(liquidFilmBase &film, const dictionary &dict),(film, dict))
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual scalar injectedMassTotal() const
Return the total mass injected.
Namespace for OpenFOAM.
static autoPtr< injectionModel > New(liquidFilmBase &film, const dictionary &dict, const word &mdoelType)
Return a reference to the selected injection model.
const liquidFilmBase & film() const
Return const access to the film surface film model.