InjectionModel.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-2024 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::InjectionModel
26 
27 Description
28  Templated injection model class.
29 
30  The injection model nominally describes the parcel:
31  - position
32  - diameter
33  - velocity
34  In this case, the fullyDescribed() flag should be set to 0 (false). When
35  the parcel is then added to the cloud, the remaining properties are
36  populated using values supplied in the constant properties.
37 
38  If, however, all of a parcel's properties are described in the model, the
39  fullyDescribed() flag should be set to 1 (true).
40 
41 SourceFiles
42  InjectionModel.C
43  InjectionModelNew.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef InjectionModel_H
48 #define InjectionModel_H
49 
50 #include "injectionModel.H"
51 #include "CloudSubModelBase.H"
52 #include "particle.H"
53 #include "Function1.H"
54 #include "runTimeSelectionTables.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class InjectionModel Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class CloudType>
66 class InjectionModel
67 :
68  public injectionModel,
69  public CloudSubModelBase<CloudType>
70 {
71 public:
72 
73  //- Convenience typedef for parcelType
74  typedef typename CloudType::parcelType parcelType;
75 
76 
77 protected:
78 
79  // Protected data
80 
81  // Global injection properties
82 
83  //- Start of injection [s]
84  scalar SOI_;
85 
86  //- Total mass injected to date [kg]
87  scalar massInjected_;
88 
89  //- Total number of parcels injected to date
91 
92  //- Mass deferred to be injected at a later time step
93  scalar massDeferred_;
94 
95  //- Number of parcels deferred to be injected at a later time step
96  scalar nParcelsDeferred_;
97 
98 
99  // Injection properties per Lagrangian time step
100 
101  //- Fixed nParticle to assign to parcels. Only valid if
102  // uniformParcelSize is nParticle.
103  scalar nParticleFixed_;
104 
105  //- Size uniform to all parcels
107 
108 
109  // Protected Member Functions
110 
111  //- Read the total mass value for instantaneous injections
112  scalar readMassTotal
113  (
114  const dictionary& dict,
116  );
117 
118  //- Read the duration for continuous injections
119  scalar readDuration
120  (
121  const dictionary& dict,
123  );
124 
125  //- Read the mass flow rate function for continuous injections
127  (
128  const dictionary& dict,
129  CloudType& owner,
130  const scalar duration
131  );
132 
133  //- Read the number of parcels injected per second for continuous
134  // injections
136  (
137  const dictionary& dict,
139  );
140 
141  //- Get the index of this injector
142  label index() const;
143 
144  //- Find the cell that contains the supplied position
145  // Will modify position slightly towards the owner cell centroid to
146  // ensure that it lies in a cell and not edge/face
147  bool findCellAtPosition
148  (
149  const point& position,
151  label& celli,
152  label& tetFacei,
153  label& tetPti,
154  bool errorOnNotFound = true
155  );
156 
157  //- Constrain a parcel's position appropriately to the geometric
158  // dimensions of the mesh
159  void constrainPosition
160  (
161  typename CloudType::parcelType::trackingData& td,
162  typename CloudType::parcelType& parcel
163  );
164 
165  //- Return the sampling moment to be used by the size distribution
166  label sizeSampleQ() const;
167 
168  //- Set number of particles to inject given parcel properties
170  (
171  PtrList<parcelType>& parcelPtrs,
172  const scalar mass
173  ) const;
174 
175  //- Pre injection hook
176  virtual void preInject
177  (
178  typename parcelType::trackingData& td
179  );
180 
181  //- Post injection hook
182  virtual void postInject
183  (
184  const label parcelsAdded,
185  const scalar massAdded,
186  typename parcelType::trackingData& td
187  );
188 
189 
190 public:
191 
192  //- Runtime type information
193  TypeName("injectionModel");
194 
195  //- Declare runtime constructor selection table
197  (
198  autoPtr,
200  dictionary,
201  (
202  const dictionary& dict,
203  CloudType& owner,
204  const word& modelType
205  ),
206  (dict, owner, modelType)
207  );
208 
209 
210  // Constructors
211 
212  //- Construct null from owner
214 
215  //- Construct from dictionary
217  (
218  const dictionary& dict,
219  CloudType& owner,
220  const word& modelName,
221  const word& modelType
222  );
223 
224  //- Construct copy
226 
227  //- Construct and return a clone
228  virtual autoPtr<InjectionModel<CloudType>> clone() const = 0;
229 
230 
231  //- Destructor
232  virtual ~InjectionModel();
233 
234 
235  // Selectors
236 
237  //- Selector with lookup from dictionary
239  (
240  const dictionary& dict,
242  );
243 
244  //- Selector with name and type
246  (
247  const dictionary& dict,
248  const word& modelName,
249  const word& modelType,
251  );
252 
253 
254  // Member Functions
255 
256  // Mapping
257 
258  //- Update mesh
259  virtual void topoChange();
260 
261 
262  // Global information
263 
264  //- Return the start-of-injection time
265  inline scalar timeStart() const;
266 
267  //- Return mass of parcels injected (cumulative)
268  inline scalar massInjected() const;
269 
270  //- Return number of parcels injected (cumulative)
271  inline label nParcelsInjected() const;
272 
273  //- Return the end-of-injection time
274  virtual scalar timeEnd() const = 0;
275 
276  //- Number of parcels to introduce relative to SOI
277  virtual scalar nParcelsToInject
278  (
279  const scalar time0,
280  const scalar time1
281  ) = 0;
282 
283  //- Parcel mass to introduce relative to SOI
284  virtual scalar massToInject
285  (
286  const scalar time0,
287  const scalar time1
288  ) = 0;
289 
290  //- Return the average injected parcel mass
291  scalar averageParcelMass();
292 
293 
294  // Per-injection event functions
295 
296  //- Main injection loop
297  template<class TrackCloudType>
298  void inject
299  (
300  TrackCloudType& cloud,
301  typename CloudType::parcelType::trackingData& td
302  );
303 
304  //- Main injection loop - steady-state
305  template<class TrackCloudType>
306  void injectSteadyState
307  (
308  TrackCloudType& cloud,
309  typename CloudType::parcelType::trackingData& td
310  );
311 
312 
313  // Injection geometry
314 
315  //- Set the injection position and owner cell, tetFace and tetPt
316  virtual void setPositionAndCell
317  (
318  const label parcelI,
319  const label nParcels,
320  const scalar time,
322  label& celli,
323  label& tetFacei,
324  label& tetPti,
325  label& facei
326  ) = 0;
327 
328  //- Set the parcel properties
329  virtual void setProperties
330  (
331  const label parcelI,
332  const label nParcels,
333  const scalar time,
334  typename parcelType::trackingData& td,
335  parcelType& parcel
336  ) = 0;
337 
338  //- Flag to identify whether model fully describes the parcel
339  virtual bool fullyDescribed() const = 0;
340 
341 
342  // I-O
343 
344  //- Write injection info to stream
345  virtual void info(Ostream& os);
346 };
347 
348 
349 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
350 
351 } // End namespace Foam
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 #define makeInjectionModel(CloudType) \
356  \
357  typedef Foam::CloudType::momentumCloudType CloudType##momentumCloudType; \
358  \
359  defineNamedTemplateTypeNameAndDebug \
360  ( \
361  Foam::InjectionModel<CloudType##momentumCloudType>, \
362  0 \
363  ); \
364  \
365  namespace Foam \
366  { \
367  defineTemplateRunTimeSelectionTable \
368  ( \
369  InjectionModel<CloudType##momentumCloudType>, \
370  dictionary \
371  ); \
372  }
373 
374 
375 #define makeInjectionModelType(SS, CloudType) \
376  \
377  typedef Foam::CloudType::momentumCloudType CloudType##momentumCloudType; \
378  \
379  defineNamedTemplateTypeNameAndDebug \
380  ( \
381  Foam::SS<CloudType##momentumCloudType>, \
382  0 \
383  ); \
384  \
385  Foam::InjectionModel<CloudType##momentumCloudType>:: \
386  adddictionaryConstructorToTable<Foam::SS<CloudType##momentumCloudType>>\
387  add##SS##CloudType##momentumCloudType##ConstructorToTable_;
388 
389 
390 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 
392 #include "InjectionModelI.H"
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 #ifdef NoRepository
397  #include "InjectionModel.C"
398 #endif
399 
400 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
401 
402 #endif
403 
404 // ************************************************************************* //
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
Templated injection model class.
virtual bool fullyDescribed() const =0
Flag to identify whether model fully describes the parcel.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)=0
Set the injection position and owner cell, tetFace and tetPt.
virtual void topoChange()
Update mesh.
virtual scalar timeEnd() const =0
Return the end-of-injection time.
virtual ~InjectionModel()
Destructor.
void inject(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td)
Main injection loop.
autoPtr< Function1< scalar > > readParcelsPerSecond(const dictionary &dict, CloudType &owner)
Read the number of parcels injected per second for continuous.
void constrainPosition(typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
Constrain a parcel's position appropriately to the geometric.
virtual void postInject(const label parcelsAdded, const scalar massAdded, typename parcelType::trackingData &td)
Post injection hook.
virtual void preInject(typename parcelType::trackingData &td)
Pre injection hook.
TypeName("injectionModel")
Runtime type information.
bool findCellAtPosition(const point &position, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
scalar readMassTotal(const dictionary &dict, CloudType &owner)
Read the total mass value for instantaneous injections.
virtual void info(Ostream &os)
Write injection info to stream.
label nParcelsInjected_
Total number of parcels injected to date.
scalar timeStart() const
Return the start-of-injection time.
void setNumberOfParticles(PtrList< parcelType > &parcelPtrs, const scalar mass) const
Set number of particles to inject given parcel properties.
autoPtr< Function1< scalar > > readMassFlowRate(const dictionary &dict, CloudType &owner, const scalar duration)
Read the mass flow rate function for continuous injections.
uniformParcelSize uniformParcelSize_
Size uniform to all parcels.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename parcelType::trackingData &td, parcelType &parcel)=0
Set the parcel properties.
scalar massDeferred_
Mass deferred to be injected at a later time step.
virtual scalar massToInject(const scalar time0, const scalar time1)=0
Parcel mass to introduce relative to SOI.
label sizeSampleQ() const
Return the sampling moment to be used by the size distribution.
static autoPtr< InjectionModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector with lookup from dictionary.
scalar averageParcelMass()
Return the average injected parcel mass.
virtual scalar nParcelsToInject(const scalar time0, const scalar time1)=0
Number of parcels to introduce relative to SOI.
declareRunTimeSelectionTable(autoPtr, InjectionModel, dictionary,(const dictionary &dict, CloudType &owner, const word &modelType),(dict, owner, modelType))
Declare runtime constructor selection table.
scalar nParticleFixed_
Fixed nParticle to assign to parcels. Only valid if.
label index() const
Get the index of this injector.
void injectSteadyState(TrackCloudType &cloud, typename CloudType::parcelType::trackingData &td)
Main injection loop - steady-state.
scalar readDuration(const dictionary &dict, CloudType &owner)
Read the duration for continuous injections.
scalar massInjected() const
Return mass of parcels injected (cumulative)
virtual autoPtr< InjectionModel< CloudType > > clone() const =0
Construct and return a clone.
InjectionModel(CloudType &owner)
Construct null from owner.
CloudType::parcelType parcelType
Convenience typedef for parcelType.
label nParcelsInjected() const
Return number of parcels injected (cumulative)
scalar SOI_
Start of injection [s].
scalar nParcelsDeferred_
Number of parcels deferred to be injected at a later time step.
scalar massInjected_
Total mass injected to date [kg].
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Non-templated base class for lagrangian injection models.
uniformParcelSize
Enumeration for the parcels' uniform size.
const dictionary & dict() const
Return const access to the cloud dictionary.
Definition: subModelBase.C:110
const word & modelName() const
Return const access to the name of the sub-model.
Definition: subModelBase.C:104
const word & modelType() const
Return const access to the sub-model type.
Definition: subModelBase.C:122
A class for handling words, derived from string.
Definition: word.H:62
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
Namespace for OpenFOAM.
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
Macros to ease declaration of run-time selection tables.