InjectedParticleInjection.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "mathematicalConstants.H"
30 #include "bitSet.H"
31 #include "injectedParticleCloud.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class CloudType>
39 {
40  const injectedParticleCloud cloud(this->owner().mesh(), cloudName_);
41 
42  label nParticles = cloud.size();
43  List<scalar> time(nParticles);
44  List<vector> position(nParticles);
45  List<scalar> diameter(nParticles);
46  List<vector> U(nParticles);
47 
48  label particlei = 0;
49 
50  for (const injectedParticle& p : cloud)
51  {
52  time[particlei] = p.soi();
53  position[particlei] = p.position() + positionOffset_;
54  diameter[particlei] = p.d();
55  U[particlei] = p.U();
56 
57  ++particlei;
58  }
59 
60  // Combine all proc data
61  if (Pstream::parRun())
62  {
64  procTime[Pstream::myProcNo()].transfer(time);
65  Pstream::gatherList(procTime);
66  Pstream::scatterList(procTime);
67  time =
68  ListListOps::combine<List<scalar>>
69  (
70  procTime, accessOp<List<scalar>>()
71  );
72 
73  List<List<point>> procPosition(Pstream::nProcs());
74  procPosition[Pstream::myProcNo()].transfer(position);
75  Pstream::gatherList(procPosition);
76  Pstream::scatterList(procPosition);
77  position =
78  ListListOps::combine<List<point>>
79  (
80  procPosition, accessOp<List<point>>()
81  );
82 
84  procD[Pstream::myProcNo()].transfer(diameter);
85  Pstream::gatherList(procD);
86  Pstream::scatterList(procD);
87  diameter =
88  ListListOps::combine<List<scalar>>
89  (
90  procD, accessOp<List<scalar>>()
91  );
92 
94  procU[Pstream::myProcNo()].transfer(U);
95  Pstream::gatherList(procU);
96  Pstream::scatterList(procU);
97  U =
98  ListListOps::combine<List<vector>>
99  (
100  procU, accessOp<List<vector>>()
101  );
102  }
103 
104  nParticles = time.size();
105 
106  // Reset SOI according to user selection
107  scalar minTime = min(time);
108  forAll(time, i)
109  {
110  time[i] -= minTime;
111  }
112 
113  // Sort and renumber to ensure lists in ascending time
114  labelList sortedIndices(Foam::sortedOrder(time));
115 
116  time_ = UIndirectList<scalar>(time, sortedIndices);
117  position_ = UIndirectList<point>(position, sortedIndices);
118  diameter_ = UIndirectList<scalar>(diameter, sortedIndices);
119  U_ = UIndirectList<vector>(U, sortedIndices);
120 
121  // Pre-calculate injected particle volumes
122  List<scalar> volume(nParticles);
123  scalar sumVolume = 0;
124  forAll(volume, i)
125  {
126  scalar vol = pow3(diameter_[i])*mathematical::pi/16.0;
127  volume[i] = vol;
128  sumVolume += vol;
129  }
130  volume_.transfer(volume);
131 
132  // Set the volume of particles to inject
133  this->volumeTotal_ = sumVolume;
134 
135  // Provide some feedback
136  Info<< " Read " << nParticles << " particles" << endl;
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
142 template<class CloudType>
144 (
145  const dictionary& dict,
146  CloudType& owner,
147  const word& modelName
148 )
149 :
150  InjectionModel<CloudType>(dict, owner, modelName, typeName),
151  cloudName_(this->coeffDict().lookup("cloud")),
152  injectorCells_(),
153  injectorTetFaces_(),
154  injectorTetPts_(),
155  time_(this->template getModelProperty<scalarList>("time")),
156  position_(this->template getModelProperty<vectorList>("position")),
157  positionOffset_(this->coeffDict().lookup("positionOffset")),
158  diameter_(this->template getModelProperty<scalarList>("diameter")),
159  U_(this->template getModelProperty<vectorList>("U")),
160  volume_(this->template getModelProperty<scalarList>("volume")),
161  ignoreOutOfBounds_
162  (
163  this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
164  ),
165  currentParticlei_
166  (
167  this->template getModelProperty<label>
168  (
169  "currentParticlei",
170  -1
171  )
172  )
173 {
174  if (this->parcelBasis_ != InjectionModel<CloudType>::pbFixed)
175  {
177  << "Injector model: " << this->modelName()
178  << " Parcel basis must be set to fixed"
179  << exit(FatalError);
180  }
181 
182  if (!time_.size())
183  {
184  // Clean start
185  initialise();
186  }
187 
188  injectorCells_.setSize(position_.size());
189  injectorTetFaces_.setSize(position_.size());
190  injectorTetPts_.setSize(position_.size());
191 
192  updateMesh();
193 
194  this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
195 }
196 
197 
198 template<class CloudType>
200 (
202 )
203 :
205  cloudName_(im.cloudName_),
206  injectorCells_(im.injectorCells_),
207  injectorTetFaces_(im.injectorTetFaces_),
208  injectorTetPts_(im.injectorTetPts_),
209  time_(im.time_),
210  position_(im.position_),
211  positionOffset_(im.positionOffset_),
212  diameter_(im.diameter_),
213  U_(im.U_),
214  volume_(im.volume_),
215  ignoreOutOfBounds_(im.ignoreOutOfBounds_),
216  currentParticlei_(im.currentParticlei_)
217 {}
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
222 template<class CloudType>
224 {
225  label nRejected = 0;
226 
227  bitSet keep(position_.size(), true);
228 
229  forAll(position_, particlei)
230  {
231  if
232  (
233  !this->findCellAtPosition
234  (
235  injectorCells_[particlei],
236  injectorTetFaces_[particlei],
237  injectorTetPts_[particlei],
238  position_[particlei],
239  !ignoreOutOfBounds_
240  )
241  )
242  {
243  keep.unset(particlei);
244  nRejected++;
245  }
246  }
247 
248  if (nRejected > 0)
249  {
250  inplaceSubset(keep, time_);
251  inplaceSubset(keep, position_);
252  inplaceSubset(keep, diameter_);
253  inplaceSubset(keep, U_);
254  inplaceSubset(keep, volume_);
255  inplaceSubset(keep, injectorCells_);
256  inplaceSubset(keep, injectorTetFaces_);
257  inplaceSubset(keep, injectorTetPts_);
258 
259  Info<< " " << nRejected
260  << " particles ignored, out of bounds" << endl;
261  }
262 }
263 
264 
265 template<class CloudType>
267 {
268  return max(time_);
269 }
270 
271 
272 template<class CloudType>
274 (
275  const scalar time0,
276  const scalar time1
277 )
278 {
279  label nParticles = 0;
280  forAll(time_, particlei)
281  {
282  if ((time_[particlei] >= time0) && (time_[particlei] < time1))
283  {
284  ++nParticles;
285  }
286  }
287 
288  return nParticles;
289 }
290 
291 
292 template<class CloudType>
294 (
295  const scalar time0,
296  const scalar time1
297 )
298 {
299  scalar sumVolume = 0;
300  forAll(time_, particlei)
301  {
302  if ((time_[particlei] >= time0) && (time_[particlei] < time1))
303  {
304  sumVolume += volume_[particlei];
305  }
306  }
307 
308  return sumVolume;
309 }
310 
311 
312 template<class CloudType>
314 (
315  const label parceli,
316  const label nParcels,
317  const scalar time,
318  vector& position,
319  label& cellOwner,
320  label& tetFacei,
321  label& tetPti
322 )
323 {
324  // Note: optimisation - consume particle from lists to reduce storage
325  // as injection proceeds
326 
327  currentParticlei_++;
328 
329  position = position_[currentParticlei_];
330  cellOwner = injectorCells_[currentParticlei_];
331  tetFacei = injectorTetFaces_[currentParticlei_];
332  tetPti = injectorTetPts_[currentParticlei_];
333 }
334 
335 
336 template<class CloudType>
338 (
339  const label parceli,
340  const label,
341  const scalar,
342  typename CloudType::parcelType& parcel
343 )
344 {
345  // Set particle velocity
346  parcel.U() = U_[currentParticlei_];
347 
348  // Set particle diameter
349  parcel.d() = diameter_[currentParticlei_];
350 }
351 
352 
353 template<class CloudType>
355 {
356  return false;
357 }
358 
359 
360 template<class CloudType>
362 (
363  const label
364 )
365 {
366  return true;
367 }
368 
369 
370 template<class CloudType>
372 {
374 
375  if (this->writeTime())
376  {
377  this->setModelProperty("currentParticlei", currentParticlei_);
378  this->setModelProperty("time", time_);
379  this->setModelProperty("position", position_);
380  this->setModelProperty("diameter", diameter_);
381  this->setModelProperty("U", U_);
382  this->setModelProperty("volume", volume_);
383  }
384 }
385 
386 
387 // ************************************************************************* //
Foam::InjectedParticleInjection::positionOffset_
vector positionOffset_
Position offset to apply to input positions.
Definition: InjectedParticleInjection.H:102
Foam::HashTable::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::InjectedParticleInjection::volumeToInject
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
Definition: InjectedParticleInjection.C:294
mathematicalConstants.H
Foam::InjectionModel::info
virtual void info(Ostream &os)
Write injection info to stream.
Definition: InjectionModel.C:669
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::accessOp
Definition: UList.H:668
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::InjectedParticleInjection::volume_
scalarList volume_
List of volume per particle [m3].
Definition: InjectedParticleInjection.H:111
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::DSMCCloud::constProps
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:98
Foam::InjectedParticleInjection::currentParticlei_
label currentParticlei_
Index of current particle.
Definition: InjectedParticleInjection.H:117
Foam::InjectedParticleInjection::fullyDescribed
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
Definition: InjectedParticleInjection.C:354
Foam::InjectedParticleInjection::injectorCells_
labelList injectorCells_
List of cell label per injector.
Definition: InjectedParticleInjection.H:87
Foam::InjectedParticleInjection::injectorTetPts_
labelList injectorTetPts_
List of tetPt label per injector.
Definition: InjectedParticleInjection.H:93
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
Foam::InjectedParticleInjection::setPositionAndCell
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.
Definition: InjectedParticleInjection.C:314
Foam::injectedParticleCloud
Definition: injectedParticleCloud.H:53
Foam::InjectionModel
Templated injection model class.
Definition: InjectionModel.H:73
Foam::InjectedParticleInjection
Replays an set of particle data based on an injectedParticleCloud, using the assumption of one partic...
Definition: InjectedParticleInjection.H:75
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fieldTypes::volume
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
bitSet.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::InjectedParticleInjection::U_
vectorList U_
List of velocity per particle [m/s].
Definition: InjectedParticleInjection.H:108
Foam::InjectedParticleInjection::setProperties
virtual void setProperties(const label parceli, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
Definition: InjectedParticleInjection.C:338
Foam::InjectedParticleInjection::validInjection
virtual bool validInjection(const label parceli)
Return flag to identify whether or not injection of parcelI is.
Definition: InjectedParticleInjection.C:362
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
InjectedParticleInjection.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::injectedParticle
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Definition: injectedParticle.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::InjectedParticleInjection::diameter_
scalarList diameter_
List of diameter per particle [m].
Definition: InjectedParticleInjection.H:105
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::InjectedParticleInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: InjectedParticleInjection.C:266
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::InjectedParticleInjection::position_
vectorList position_
List of position per particle [m].
Definition: InjectedParticleInjection.H:99
Foam::InjectedParticleInjection::initialise
void initialise()
Initialise injectors.
Definition: InjectedParticleInjection.C:38
Foam::InjectedParticleInjection::info
void info(Ostream &os)
Write injection info to stream.
Definition: InjectedParticleInjection.C:371
Foam::inplaceSubset
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
Definition: ListOpsTemplates.C:583
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
U
U
Definition: pEqn.H:72
Foam::InjectedParticleInjection::InjectedParticleInjection
InjectedParticleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: InjectedParticleInjection.C:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< scalar >
Foam::InjectedParticleInjection::parcelsToInject
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
Definition: InjectedParticleInjection.C:274
Foam::InjectedParticleInjection::cloudName_
const word cloudName_
Name of cloud used to seed the new particles.
Definition: InjectedParticleInjection.H:84
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
injectedParticleCloud.H
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::InjectedParticleInjection::ignoreOutOfBounds_
Switch ignoreOutOfBounds_
Flag to suppress errors if particle injection site is out-of-bounds.
Definition: InjectedParticleInjection.H:114
Foam::InjectedParticleInjection::updateMesh
virtual void updateMesh()
Set injector locations when mesh is updated.
Definition: InjectedParticleInjection.C:223
Foam::InjectedParticleInjection::time_
scalarList time_
List of injection time per particle [s].
Definition: InjectedParticleInjection.H:96
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::InjectedParticleInjection::injectorTetFaces_
labelList injectorTetFaces_
List of tetFace label per injector.
Definition: InjectedParticleInjection.H:90