34 template<
class CloudType>
40 const fvMesh& mesh = this->owner().mesh();
42 randomGenerator&
rndGen = this->owner().rndGen();
44 const label nCells = cellZoneCells.size();
45 DynamicList<barycentric> injectorCoordinates(nCells);
46 DynamicList<label> injectorCells(nCells);
47 DynamicList<label> injectorTetFaces(nCells);
48 DynamicList<label> injectorTetPts(nCells);
50 scalar newParticlesTotal = 0;
51 label addParticlesTotal = 0;
55 const label celli = cellZoneCells[i];
58 const scalar newParticles = mesh.V()[celli]*numberDensity_;
59 newParticlesTotal += newParticles;
60 label addParticles = floor(newParticles);
61 addParticlesTotal += addParticles;
63 const scalar
diff = newParticlesTotal - addParticlesTotal;
68 addParticlesTotal += corr;
72 const List<tetIndices> cellTetIs =
73 polyMeshTetDecomposition::cellTetIndices(mesh, celli);
77 cTetVFrac[0] = cellTetIs[0].tet(mesh).mag();
78 for (
label tetI = 1; tetI < cellTetIs.size(); tetI++)
81 cTetVFrac[tetI-1] + cellTetIs[tetI].tet(mesh).mag();
83 cTetVFrac /= cTetVFrac.last();
86 for (
label pI = 0; pI < addParticles; pI++)
88 const scalar volFrac =
rndGen.sample01<scalar>();
92 if (cTetVFrac[vfI] > volFrac)
100 injectorCells.append(celli);
101 injectorTetFaces.append(cellTetIs[tetI].face());
102 injectorTetPts.append(cellTetIs[tetI].tetPt());
108 globalIndex globalPositions(injectorCoordinates.size());
110 List<barycentric> allInjectorCoordinates
112 globalPositions.size(),
113 barycentric::uniform(NaN)
115 List<label> allInjectorCells(globalPositions.size(), -1);
116 List<label> allInjectorTetFaces(globalPositions.size(), -1);
117 List<label> allInjectorTetPts(globalPositions.size(), -1);
121 allInjectorCoordinates,
122 globalPositions.localSize(Pstream::myProcNo()),
123 globalPositions.offset(Pstream::myProcNo())
124 ) = injectorCoordinates;
128 globalPositions.localSize(Pstream::myProcNo()),
129 globalPositions.offset(Pstream::myProcNo())
134 globalPositions.localSize(Pstream::myProcNo()),
135 globalPositions.offset(Pstream::myProcNo())
136 ) = injectorTetFaces;
140 globalPositions.localSize(Pstream::myProcNo()),
141 globalPositions.offset(Pstream::myProcNo())
145 injectorCoordinates_.transfer(allInjectorCoordinates);
146 injectorCells_.transfer(allInjectorCells);
147 injectorTetFaces_.transfer(allInjectorTetFaces);
148 injectorTetPts_.transfer(allInjectorTetPts);
154 template<
class CloudType>
159 const word& modelName
163 cellZoneName_(this->coeffDict().lookup(
"cellZone")),
164 massTotal_(this->readMassTotal(
dict, owner)),
165 numberDensity_(this->coeffDict().template lookup<scalar>(
"numberDensity")),
166 injectorCoordinates_(),
171 U0_(this->coeffDict().lookup(
"U0")),
177 this->coeffDict().subDict(
"sizeDistribution"),
179 owner.
rndGen().generator()
187 template<
class CloudType>
194 cellZoneName_(im.cellZoneName_),
195 massTotal_(im.massTotal_),
196 numberDensity_(im.numberDensity_),
197 injectorCoordinates_(im.injectorCoordinates_),
198 injectorCells_(im.injectorCells_),
199 injectorTetFaces_(im.injectorTetFaces_),
200 injectorTetPts_(im.injectorTetPts_),
201 diameters_(im.diameters_),
203 sizeDistribution_(im.sizeDistribution_, false)
209 template<
class CloudType>
216 template<
class CloudType>
226 <<
"Unknown cell zone name: " << cellZoneName_
232 const label nCells = cellZoneCells.
size();
236 Info<<
" cell zone size = " << nCellsTotal <<
endl;
237 Info<<
" cell zone volume = " << VCellsTotal <<
endl;
239 if ((nCellsTotal == 0) || (VCellsTotal*numberDensity_ < 1))
242 <<
"Number of particles to be added to cellZone " << cellZoneName_
243 <<
" is zero" <<
endl;
247 setPositions(cellZoneCells);
249 Info<<
" number density = " << numberDensity_ <<
nl
250 <<
" number of particles = " << injectorCoordinates_.size()
254 diameters_.setSize(injectorCoordinates_.size());
257 diameters_[i] = sizeDistribution_->sample();
263 template<
class CloudType>
271 template<
class CloudType>
279 if (0 >= time0 && 0 < time1)
281 return injectorCoordinates_.size();
290 template<
class CloudType>
298 if (0 >= time0 && 0 < time1)
309 template<
class CloudType>
313 const label nParcels,
322 coordinates = injectorCoordinates_[parcelI];
323 celli = injectorCells_[parcelI];
324 tetFacei = injectorTetFaces_[parcelI];
325 tetPti = injectorTetPts_[parcelI];
329 template<
class CloudType>
335 typename CloudType::parcelType::trackingData& td,
343 parcel.d() = diameters_[parcelI];
347 template<
class CloudType>
#define forAll(list, i)
Loop across all elements in list.
Injection positions specified by a particle number density within a cell set.
virtual void topoChange()
Set injector locations when mesh is updated.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType::trackingData &td, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual void setPositionAndCell(const label parcelI, const label nParcels, const scalar time, barycentric &coordinates, label &celli, label &tetFacei, label &tetPti, label &facei)
Set the injection position and owner cell, tetFace and tetPt.
virtual ~CellZoneInjection()
Destructor.
virtual label nParcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
CellZoneInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
virtual scalar timeEnd() const
Return the end-of-injection time.
virtual scalar massToInject(const scalar time0, const scalar time1)
Parcel mass to introduce relative to SOI.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
wordList toc() const
Return the table of contents.
Templated injection model class.
void size(const label)
Override size to be inconsistent with allocated storage.
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Base class for statistical distributions.
Mesh data needed to do the Finite Volume discretisation.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const polyMesh & mesh() const
Return reference to polyMesh.
const cellZoneList & cellZones() const
Return cell zones.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField scalarField(fieldObject, mesh)
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
barycentric barycentric01(randomGenerator &rndGen)
Generate a random barycentric coordinate within the unit tetrahedron.
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
randomGenerator rndGen(653213)