48 #ifndef InjectionModel_H
49 #define InjectionModel_H
66 template<
class CloudType>
137 const scalar duration
156 const point& position,
161 bool errorOnNotFound =
true
168 typename CloudType::parcelType::trackingData& td,
185 const label parcelsAdded,
186 const scalar massAdded
271 virtual scalar
timeEnd()
const = 0;
303 template<
class TrackCloudType>
306 TrackCloudType&
cloud,
307 typename CloudType::parcelType::trackingData& td
311 template<
class TrackCloudType>
314 TrackCloudType&
cloud,
315 typename CloudType::parcelType::trackingData& td
325 const label nParcels,
338 const label nParcels,
360 #define makeInjectionModel(CloudType) \
362 typedef Foam::CloudType::momentumCloudType momentumCloudType; \
363 defineNamedTemplateTypeNameAndDebug \
365 Foam::InjectionModel<momentumCloudType>, \
371 defineTemplateRunTimeSelectionTable \
373 InjectionModel<momentumCloudType>, \
379 #define makeInjectionModelType(SS, CloudType) \
381 typedef Foam::CloudType::momentumCloudType momentumCloudType; \
382 defineNamedTemplateTypeNameAndDebug(Foam::SS<momentumCloudType>, 0); \
384 Foam::InjectionModel<momentumCloudType>:: \
385 adddictionaryConstructorToTable<Foam::SS<momentumCloudType>> \
386 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.
scalar timeStep0_
Time at start of injection time step [s].
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.
void constrainPosition(typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
Constrain a parcel's position appropriately to the geometric.
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.
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.
uniformParcelSize uniformParcelSize_
Size uniform to all parcels.
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.
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.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, parcelType &parcel)=0
Set the parcel properties.
TimeFunction1< scalar > readParcelsPerSecond(const dictionary &dict, CloudType &owner)
Read the number of parcels injected per second for continuous.
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.
label parcelsAddedTotal() const
Return the total number parcels added.
scalar massInjected() const
Return mass of particles injected (cumulative)
virtual autoPtr< InjectionModel< CloudType > > clone() const =0
Construct and return a clone.
virtual label nParcelsToInject(const scalar time0, const scalar time1)=0
Number of parcels to introduce relative to SOI.
label nInjections_
Number of injections counter.
InjectionModel(CloudType &owner)
Construct null from owner.
TimeFunction1< scalar > readMassFlowRate(const dictionary &dict, CloudType &owner, const scalar duration)
Read the mass flow rate function for continuous injections.
label parcelsAddedTotal_
Running counter of total number of parcels added.
CloudType::parcelType parcelType
Convenience typedef for parcelType.
scalar time0_
Continuous phase time at start of injection time step [s].
scalar SOI_
Start of injection [s].
label nInjections() const
Return the number of injections.
void postInjectCheck(const label parcelsAdded, const scalar massAdded)
Post injection checks.
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...
A cloud is a collection of lagrangian particles.
A list of keyword definitions, which are a keyword followed by any number of values (e....
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.
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.