CellZoneInjection.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "CellZoneInjection.H"
27 #include "mathematicalConstants.H"
29 #include "globalIndex.H"
30 #include "Pstream.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const labelList& cellZoneCells
38 )
39 {
40  const fvMesh& mesh = this->owner().mesh();
41 
42  Random& rnd = this->owner().rndGen();
43 
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);
49 
50  scalar newParticlesTotal = 0;
51  label addParticlesTotal = 0;
52 
53  forAll(cellZoneCells, i)
54  {
55  const label celli = cellZoneCells[i];
56 
57  // Calc number of particles to add
58  const scalar newParticles = mesh.V()[celli]*numberDensity_;
59  newParticlesTotal += newParticles;
60  label addParticles = floor(newParticles);
61  addParticlesTotal += addParticles;
62 
63  const scalar diff = newParticlesTotal - addParticlesTotal;
64  if (diff > 1)
65  {
66  label corr = floor(diff);
67  addParticles += corr;
68  addParticlesTotal += corr;
69  }
70 
71  // Construct cell tet indices
72  const List<tetIndices> cellTetIs =
73  polyMeshTetDecomposition::cellTetIndices(mesh, celli);
74 
75  // Construct cell tet volume fractions
76  scalarField cTetVFrac(cellTetIs.size());
77  cTetVFrac[0] = cellTetIs[0].tet(mesh).mag();
78  for (label tetI = 1; tetI < cellTetIs.size(); tetI++)
79  {
80  cTetVFrac[tetI] =
81  cTetVFrac[tetI-1] + cellTetIs[tetI].tet(mesh).mag();
82  }
83  cTetVFrac /= cTetVFrac.last();
84 
85  // Set new particle position and cellId
86  for (label pI = 0; pI < addParticles; pI++)
87  {
88  const scalar volFrac = rnd.sample01<scalar>();
89  label tetI = 0;
90  forAll(cTetVFrac, vfI)
91  {
92  if (cTetVFrac[vfI] > volFrac)
93  {
94  tetI = vfI;
95  break;
96  }
97  }
98 
99  injectorCoordinates.append(barycentric01(rnd));
100  injectorCells.append(celli);
101  injectorTetFaces.append(cellTetIs[tetI].face());
102  injectorTetPts.append(cellTetIs[tetI].tetPt());
103  }
104  }
105 
106  // Accumulate into global lists. Set the coordinates and topology for local
107  // particles, and leave remote ones with invalid data.
108  globalIndex globalPositions(injectorCoordinates.size());
109 
110  List<barycentric> allInjectorCoordinates
111  (
112  globalPositions.size(),
113  barycentric::uniform(NaN)
114  );
115  List<label> allInjectorCells(globalPositions.size(), -1);
116  List<label> allInjectorTetFaces(globalPositions.size(), -1);
117  List<label> allInjectorTetPts(globalPositions.size(), -1);
118 
119  SubList<barycentric>
120  (
121  allInjectorCoordinates,
122  globalPositions.localSize(Pstream::myProcNo()),
123  globalPositions.offset(Pstream::myProcNo())
124  ) = injectorCoordinates;
125  SubList<label>
126  (
127  allInjectorCells,
128  globalPositions.localSize(Pstream::myProcNo()),
129  globalPositions.offset(Pstream::myProcNo())
130  ) = injectorCells;
131  SubList<label>
132  (
133  allInjectorTetFaces,
134  globalPositions.localSize(Pstream::myProcNo()),
135  globalPositions.offset(Pstream::myProcNo())
136  ) = injectorTetFaces;
137  SubList<label>
138  (
139  allInjectorTetPts,
140  globalPositions.localSize(Pstream::myProcNo()),
141  globalPositions.offset(Pstream::myProcNo())
142  ) = injectorTetPts;
143 
144  // Transfer data
145  injectorCoordinates_.transfer(allInjectorCoordinates);
146  injectorCells_.transfer(allInjectorCells);
147  injectorTetFaces_.transfer(allInjectorTetFaces);
148  injectorTetPts_.transfer(allInjectorTetPts);
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
153 
154 template<class CloudType>
156 (
157  const dictionary& dict,
158  CloudType& owner,
159  const word& modelName
160 )
161 :
162  InjectionModel<CloudType>(dict, owner, modelName, typeName),
163  cellZoneName_(this->coeffDict().lookup("cellZone")),
164  massTotal_(this->readMassTotal(dict, owner)),
165  numberDensity_(this->coeffDict().template lookup<scalar>("numberDensity")),
166  injectorCoordinates_(),
167  injectorCells_(),
168  injectorTetFaces_(),
169  injectorTetPts_(),
170  diameters_(),
171  U0_(this->coeffDict().lookup("U0")),
172  sizeDistribution_
173  (
175  (
176  this->coeffDict().subDict("sizeDistribution"),
177  owner.rndGen(),
178  this->sizeSampleQ()
179  )
180  )
181 {
182  topoChange();
183 }
184 
185 
186 template<class CloudType>
188 (
190 )
191 :
193  cellZoneName_(im.cellZoneName_),
194  massTotal_(im.massTotal_),
195  numberDensity_(im.numberDensity_),
196  injectorCoordinates_(im.injectorCoordinates_),
197  injectorCells_(im.injectorCells_),
198  injectorTetFaces_(im.injectorTetFaces_),
199  injectorTetPts_(im.injectorTetPts_),
200  diameters_(im.diameters_),
201  U0_(im.U0_),
202  sizeDistribution_(im.sizeDistribution_().clone().ptr())
203 {}
204 
205 
206 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
207 
208 template<class CloudType>
210 {}
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
215 template<class CloudType>
217 {
218  // Set/cache the injector cells
219  const fvMesh& mesh = this->owner().mesh();
220  const label zoneI = mesh.cellZones().findZoneID(cellZoneName_);
221 
222  if (zoneI < 0)
223  {
225  << "Unknown cell zone name: " << cellZoneName_
226  << ". Valid cell zones are: " << mesh.cellZones().names()
227  << nl << exit(FatalError);
228  }
229 
230  const labelList& cellZoneCells = mesh.cellZones()[zoneI];
231  const label nCells = cellZoneCells.size();
232  const scalar nCellsTotal = returnReduce(nCells, sumOp<label>());
233  const scalar VCells = sum(scalarField(mesh.V(), cellZoneCells));
234  const scalar VCellsTotal = returnReduce(VCells, sumOp<scalar>());
235  Info<< " cell zone size = " << nCellsTotal << endl;
236  Info<< " cell zone volume = " << VCellsTotal << endl;
237 
238  if ((nCellsTotal == 0) || (VCellsTotal*numberDensity_ < 1))
239  {
241  << "Number of particles to be added to cellZone " << cellZoneName_
242  << " is zero" << endl;
243  }
244  else
245  {
246  setPositions(cellZoneCells);
247 
248  Info<< " number density = " << numberDensity_ << nl
249  << " number of particles = " << injectorCoordinates_.size()
250  << endl;
251 
252  // Construct parcel diameters
253  diameters_.setSize(injectorCoordinates_.size());
254  forAll(diameters_, i)
255  {
256  diameters_[i] = sizeDistribution_->sample();
257  }
258  }
259 }
260 
261 
262 template<class CloudType>
264 {
265  // Not used
266  return this->SOI_;
267 }
268 
269 
270 template<class CloudType>
272 (
273  const scalar time0,
274  const scalar time1
275 )
276 {
277  // All parcels introduced at SOI
278  if (0 >= time0 && 0 < time1)
279  {
280  return injectorCoordinates_.size();
281  }
282  else
283  {
284  return 0;
285  }
286 }
287 
288 
289 template<class CloudType>
291 (
292  const scalar time0,
293  const scalar time1
294 )
295 {
296  // All parcels introduced at SOI
297  if (0 >= time0 && 0 < time1)
298  {
299  return massTotal_;
300  }
301  else
302  {
303  return 0;
304  }
305 }
306 
307 
308 template<class CloudType>
310 (
311  const label parcelI,
312  const label nParcels,
313  const scalar time,
314  barycentric& coordinates,
315  label& celli,
316  label& tetFacei,
317  label& tetPti,
318  label& facei
319 )
320 {
321  coordinates = injectorCoordinates_[parcelI];
322  celli = injectorCells_[parcelI];
323  tetFacei = injectorTetFaces_[parcelI];
324  tetPti = injectorTetPts_[parcelI];
325 }
326 
327 
328 template<class CloudType>
330 (
331  const label parcelI,
332  const label,
333  const scalar,
334  typename CloudType::parcelType& parcel
335 )
336 {
337  // set particle velocity
338  parcel.U() = U0_;
339 
340  // set particle diameter
341  parcel.d() = diameters_[parcelI];
342 }
343 
344 
345 template<class CloudType>
347 {
348  return false;
349 }
350 
351 
352 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Injection positions specified by a particle number density within a cell set.
virtual void topoChange()
Set injector locations when mesh is updated.
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 void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
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.
Definition: DSMCCloud.H:79
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Templated injection model class.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
wordList names() const
Return a list of zone names.
Definition: MeshZones.C:256
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:451
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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)
static Type NaN()
Return a primitive with all components set to NaN.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
barycentric barycentric01(Random &rndGen)
Generate a random barycentric coordinate within the unit tetrahedron.
Definition: barycentric.C:31
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
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.
Definition: triad.C:403
T clone(const T &t)
Definition: List.H:55
error FatalError
static const char nl
Definition: Ostream.H:260
dictionary dict
Random rndGen(label(0))