47 #ifndef InjectionModel_H
48 #define InjectionModel_H
65 template<
class CloudType>
130 const scalar duration
154 bool errorOnNotFound =
true
161 typename CloudType::parcelType::trackingData& td,
178 typename parcelType::trackingData& td
184 const label parcelsAdded,
185 const scalar massAdded,
186 typename parcelType::trackingData& td
274 virtual scalar
timeEnd()
const = 0;
297 template<
class TrackCloudType>
300 TrackCloudType&
cloud,
301 typename CloudType::parcelType::trackingData& td
305 template<
class TrackCloudType>
308 TrackCloudType&
cloud,
309 typename CloudType::parcelType::trackingData& td
319 const label nParcels,
332 const label nParcels,
334 typename parcelType::trackingData& td,
355 #define makeInjectionModel(CloudType) \
357 typedef Foam::CloudType::momentumCloudType CloudType##momentumCloudType; \
359 defineNamedTemplateTypeNameAndDebug \
361 Foam::InjectionModel<CloudType##momentumCloudType>, \
367 defineTemplateRunTimeSelectionTable \
369 InjectionModel<CloudType##momentumCloudType>, \
375 #define makeInjectionModelType(SS, CloudType) \
377 typedef Foam::CloudType::momentumCloudType CloudType##momentumCloudType; \
379 defineNamedTemplateTypeNameAndDebug \
381 Foam::SS<CloudType##momentumCloudType>, \
385 Foam::InjectionModel<CloudType##momentumCloudType>:: \
386 adddictionaryConstructorToTable<Foam::SS<CloudType##momentumCloudType>>\
387 add##SS##CloudType##momentumCloudType##ConstructorToTable_;
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
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,...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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.
const word & modelName() const
Return const access to the name of the sub-model.
const word & modelType() const
Return const access to the sub-model type.
A class for handling words, derived from string.
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.
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.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Macros to ease declaration of run-time selection tables.