38template<
class CloudType>
48 injectorCells_(positionAxis_.size()),
49 injectorTetFaces_(positionAxis_.size()),
50 injectorTetPts_(positionAxis_.size()),
51 duration_(this->
coeffDict().getScalar(
"duration")),
54 this->
coeffDict().getScalar(
"parcelsPerInjector")
100 injectorOrder_(
identity(positionAxis_.size())),
106 tanVec1_.setSize(positionAxis_.size());
107 tanVec2_.setSize(positionAxis_.size());
112 flowRateProfile_->userTimeToTime(time);
113 Umag_->userTimeToTime(time);
114 thetaInner_->userTimeToTime(time);
115 thetaOuter_->userTimeToTime(time);
121 vector& axis = positionAxis_[i].second();
125 scalar magTangent = 0.0;
128 while (magTangent < SMALL)
132 tangent = v - (v & axis)*axis;
133 magTangent =
mag(tangent);
136 tanVec1_[i] = tangent/magTangent;
137 tanVec2_[i] = axis^tanVec1_[i];
141 this->
volumeTotal_ = flowRateProfile_->integrate(0.0, duration_);
145template<
class CloudType>
152 positionAxis_(im.positionAxis_),
153 injectorCells_(im.injectorCells_),
154 injectorTetFaces_(im.injectorTetFaces_),
155 injectorTetPts_(im.injectorTetPts_),
156 duration_(im.duration_),
157 parcelsPerInjector_(im.parcelsPerInjector_),
158 flowRateProfile_(im.flowRateProfile_.clone()),
159 Umag_(im.Umag_.clone()),
160 thetaInner_(im.thetaInner_.clone()),
161 thetaOuter_(im.thetaOuter_.clone()),
162 sizeDistribution_(im.sizeDistribution_.clone()),
163 nInjected_(im.nInjected_),
164 injectorOrder_(im.injectorOrder_),
165 tanVec1_(im.tanVec1_),
166 tanVec2_(im.tanVec2_)
172template<
class CloudType>
175 bitSet reject(positionAxis_.size());
182 !this->findCellAtPosition
185 injectorTetFaces_[i],
187 positionAxis_[i].first(),
188 !this->ignoreOutOfBounds_
197 const label nRejected = reject.count();
208 <<
" positions rejected, out of bounds" <<
endl;
213template<
class CloudType>
216 return this->
SOI_ + duration_;
220template<
class CloudType>
227 if ((time0 >= 0.0) && (time0 < duration_))
229 const scalar targetVolume = flowRateProfile_->integrate(0, time1);
231 const scalar volumeFraction = targetVolume/this->volumeTotal_;
233 const label targetParcels =
234 ceil(positionAxis_.size()*parcelsPerInjector_*volumeFraction);
243template<
class CloudType>
250 if ((time0 >= 0.0) && (time0 < duration_))
252 return flowRateProfile_->integrate(time0, time1);
259template<
class CloudType>
274 const label i = injectorOrder_[parcelI % positionAxis_.size()];
275 position = positionAxis_[i].first();
276 cellOwner = injectorCells_[i];
277 tetFacei = injectorTetFaces_[i];
278 tetPti = injectorTetPts_[i];
282template<
class CloudType>
293 const label i = injectorOrder_[parcelI % positionAxis_.size()];
296 scalar t = time - this->SOI_;
297 scalar ti = thetaInner_->value(t);
298 scalar to = thetaOuter_->value(t);
302 scalar dcorr =
cos(coneAngle);
306 vector dirVec = dcorr*positionAxis_[i].second();
311 parcel.
U() = Umag_->value(t)*dirVec;
314 parcel.d() = sizeDistribution_().sample();
321template<
class CloudType>
328template<
class CloudType>
const CloudType & owner() const
Return const access to the owner cloud.
Multi-point cone injection model.
ConeInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual autoPtr< InjectionModel< CloudType > > clone() const
Construct and return a clone.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, vector &position, label &cellOwner, label &tetFacei, label &tetPti)
Set the injection position and owner cell, tetFace and tetPt.
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Random & rndGen()
Return reference to the random object.
const vector & U() const
Return const access to velocity.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Templated injection model class.
static autoPtr< InjectionModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector with lookup from dictionary.
bool ignoreOutOfBounds_
Optional flag to indicate that injections attempted outside the mesh should be ignored.
label parcelsAddedTotal() const
Return the total number parcels added.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
InjectionModel(CloudType &owner)
Construct null from owner.
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
scalar SOI_
Start of injection [s].
Inter-processor communications stream.
Type position(const Type &start, const Type &end)
Return a sample on the interval [start,end].
Type sample01()
Return a sample whose components lie in the range [0,1].
void shuffle(UList< Type > &values)
Shuffle the values in the list.
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
unsigned int count(const bool on=true) const
Count number of bits set.
void flip()
Invert all bits in the addressable region.
void set(const bitSet &bitset)
Set specified bits from another bitset.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Lookup type of boundary radiation properties.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
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.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
constexpr scalar twoPi(2 *M_PI)
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
static constexpr const zero Zero
Global zero (0).
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Unit conversion functions.