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