Go to the documentation of this file.
34 template<
class CloudType>
42 const label nCells = cellZoneCells.
size();
50 scalar newParticlesTotal = 0.0;
51 label addParticlesTotal = 0;
55 const label cellI = cellZoneCells[i];
58 const scalar newParticles = V[cellI]*numberDensity_;
59 newParticlesTotal += newParticles;
60 label addParticles = floor(newParticles);
61 addParticlesTotal += addParticles;
63 const scalar
diff = newParticlesTotal - addParticlesTotal;
68 addParticlesTotal += corr;
73 polyMeshTetDecomposition::cellTetIndices(
mesh, cellI);
77 for (
label tetI = 1; tetI < cellTetIs.
size() - 1; tetI++)
80 cTetVFrac[tetI-1] + cellTetIs[tetI].tet(
mesh).mag()/V[cellI];
82 cTetVFrac.last() = 1.0;
85 for (
label pI = 0; pI < addParticles; pI++)
87 const scalar volFrac = rnd.
sample01<scalar>();
91 if (cTetVFrac[vfI] > volFrac)
97 positions.
append(cellTetIs[tetI].tet(
mesh).randomPoint(rnd));
99 injectorCells.
append(cellI);
100 injectorTetFaces.
append(cellTetIs[tetI].
face());
101 injectorTetPts.
append(cellTetIs[tetI].tetPt());
116 globalPositions.
localSize(Pstream::myProcNo()),
117 globalPositions.
offset(Pstream::myProcNo())
121 Pstream::listCombineScatter(allPositions);
127 globalPositions.
localSize(Pstream::myProcNo()),
128 globalPositions.
offset(Pstream::myProcNo())
133 globalPositions.
localSize(Pstream::myProcNo()),
134 globalPositions.
offset(Pstream::myProcNo())
135 ).
assign(injectorTetFaces);
139 globalPositions.
localSize(Pstream::myProcNo()),
140 globalPositions.
offset(Pstream::myProcNo())
144 positions_.transfer(allPositions);
145 injectorCells_.transfer(allInjectorCells);
146 injectorTetFaces_.transfer(allInjectorTetFaces);
147 injectorTetPts_.transfer(allInjectorTetPts);
163 template<
class CloudType>
168 const word& modelName
172 cellZoneName_(this->coeffDict().
lookup(
"cellZone")),
179 U0_(this->coeffDict().
lookup(
"U0")),
184 this->coeffDict().subDict(
"sizeDistribution"), owner.
rndGen()
192 template<
class CloudType>
213 template<
class CloudType>
220 template<
class CloudType>
225 const label zoneI =
mesh.cellZones().findZoneID(cellZoneName_);
230 <<
"Unknown cell zone name: " << cellZoneName_
231 <<
". Valid cell zones are: " <<
mesh.cellZones().names()
236 const label nCells = cellZoneCells.
size();
240 Info<<
" cell zone size = " << nCellsTotal <<
endl;
241 Info<<
" cell zone volume = " << VCellsTotal <<
endl;
243 if ((nCellsTotal == 0) || (VCellsTotal*numberDensity_ < 1))
246 <<
"Number of particles to be added to cellZone " << cellZoneName_
247 <<
" is zero" <<
endl;
251 setPositions(cellZoneCells);
253 Info<<
" number density = " << numberDensity_ <<
nl
254 <<
" number of particles = " << positions_.size() <<
endl;
257 diameters_.setSize(positions_.size());
260 diameters_[i] = sizeDistribution_->sample();
269 template<
class CloudType>
277 template<
class CloudType>
284 if ((0.0 >= time0) && (0.0 < time1))
286 return positions_.size();
295 template<
class CloudType>
303 if ((0.0 >= time0) && (0.0 < time1))
305 return this->volumeTotal_;
314 template<
class CloudType>
318 const label nParcels,
326 position = positions_[parcelI];
327 cellOwner = injectorCells_[parcelI];
328 tetFaceI = injectorTetFaces_[parcelI];
329 tetPtI = injectorTetPts_[parcelI];
333 template<
class CloudType>
346 parcel.d() = diameters_[parcelI];
350 template<
class CloudType>
357 template<
class CloudType>
labelList injectorTetFaces_
List of tetFace labels corresponding to injector positions.
Type sample01()
Return a sample whose components lie in the range 0-1.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Injection positions specified by a particle number density within a cell set.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
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.
#define forAll(list, i)
Loop across all elements in list.
void writeOBJ(Ostream &os, label &vertI, const tetPoints &tet)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A List obtained as a section of another List.
Templated injection model class.
const word cellZoneName_
Name of cell zone.
label localSize() const
My local size.
labelList injectorCells_
List of cell labels corresponding to injector positions.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Random & rndGen()
Return refernce to the random object.
const vector U0_
Initial parcel velocity.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
dimensionedScalar pow3(const dimensionedScalar &ds)
const fvMesh & mesh() const
Return refernce to the mesh.
scalar timeEnd() const
Return the end-of-injection time.
void assign(const UList< T > &)
Assign elements to those from UList.
label offset(const label procI) const
Start of procI data.
Templated base class for dsmc cloud.
label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const autoPtr< distributionModels::distributionModel > sizeDistribution_
Parcel size distribution model.
virtual ~CellZoneInjection()
Destructor.
Mesh data needed to do the Finite Volume discretisation.
CellZoneInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
virtual void updateMesh()
Set injector locations when mesh is updated.
label size() const
Global sum of localSizes.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
volScalarField scalarField(fieldObject, mesh)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
A face is a list of labels corresponding to mesh vertices.
void size(const label)
Override size to be inconsistent with allocated storage.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
#define WarningInFunction
Report a warning using Foam::Warning.
scalarList diameters_
Field of parcel diameters.
labelList injectorTetPts_
List of tetPt labels corresponding to injector positions.
const scalar numberDensity_
Number density.
void setPositions(const labelList &cellZoneCells)
Set the parcel injection positions.
stressControl lookup("compactNormalStress") >> compactNormalStress
List< vector > positions_
Field of parcel positions.