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 "SortableList.H"
32 #include "injectedParticleCloud.H"
33 
34 using namespace Foam::constant;
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class CloudType>
40 {
41  const injectedParticleCloud cloud(this->owner().mesh(), cloudName_);
42 
43  label nParticles = cloud.size();
44  List<scalar> time(nParticles);
45  List<vector> position(nParticles);
46  List<scalar> diameter(nParticles);
47  List<vector> U(nParticles);
48 
49  label particlei = 0;
50 
51  for (const injectedParticle& p : cloud)
52  {
53  time[particlei] = p.soi();
54  position[particlei] = p.position() + positionOffset_;
55  diameter[particlei] = p.d();
56  U[particlei] = p.U();
57 
58  ++particlei;
59  }
60 
61  // Combine all proc data
62  if (Pstream::parRun())
63  {
65  procTime[Pstream::myProcNo()].transfer(time);
66  Pstream::gatherList(procTime);
67  Pstream::scatterList(procTime);
68  time =
69  ListListOps::combine<List<scalar>>
70  (
71  procTime, accessOp<List<scalar>>()
72  );
73 
74  List<List<point>> procPosition(Pstream::nProcs());
75  procPosition[Pstream::myProcNo()].transfer(position);
76  Pstream::gatherList(procPosition);
77  Pstream::scatterList(procPosition);
78  position =
79  ListListOps::combine<List<point>>
80  (
81  procPosition, accessOp<List<point>>()
82  );
83 
85  procD[Pstream::myProcNo()].transfer(diameter);
86  Pstream::gatherList(procD);
87  Pstream::scatterList(procD);
88  diameter =
89  ListListOps::combine<List<scalar>>
90  (
91  procD, accessOp<List<scalar>>()
92  );
93 
95  procU[Pstream::myProcNo()].transfer(U);
96  Pstream::gatherList(procU);
97  Pstream::scatterList(procU);
98  U =
99  ListListOps::combine<List<vector>>
100  (
101  procU, accessOp<List<vector>>()
102  );
103  }
104 
105  nParticles = time.size();
106 
107  // Reset SOI according to user selection
108  scalar minTime = min(time);
109  forAll(time, i)
110  {
111  time[i] -= minTime;
112  }
113 
114  // Sort and renumber to ensure lists in ascending time
115  labelList sortedIndices(Foam::sortedOrder(time));
116 
117  time_ = UIndirectList<scalar>(time, sortedIndices);
118  position_ = UIndirectList<point>(position, sortedIndices);
119  diameter_ = UIndirectList<scalar>(diameter, sortedIndices);
120  U_ = UIndirectList<vector>(U, sortedIndices);
121 
122  // Pre-calculate injected particle volumes
123  List<scalar> volume(nParticles);
124  scalar sumVolume = 0;
125  forAll(volume, i)
126  {
127  scalar vol = pow3(diameter_[i])*mathematical::pi/16.0;
128  volume[i] = vol;
129  sumVolume += vol;
130  }
131  volume_.transfer(volume);
132 
133  // Set the volume of particles to inject
134  this->volumeTotal_ = sumVolume;
135 
136  // Provide some feedback
137  Info<< " Read " << nParticles << " particles" << endl;
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
143 template<class CloudType>
145 (
146  const dictionary& dict,
147  CloudType& owner,
148  const word& modelName
149 )
150 :
151  InjectionModel<CloudType>(dict, owner, modelName, typeName),
152  cloudName_(this->coeffDict().lookup("cloud")),
153  injectorCells_(),
154  injectorTetFaces_(),
155  injectorTetPts_(),
156  time_(this->template getModelProperty<scalarList>("time")),
157  position_(this->template getModelProperty<vectorList>("position")),
158  positionOffset_(this->coeffDict().lookup("positionOffset")),
159  diameter_(this->template getModelProperty<scalarList>("diameter")),
160  U_(this->template getModelProperty<vectorList>("U")),
161  volume_(this->template getModelProperty<scalarList>("volume")),
162  ignoreOutOfBounds_
163  (
164  this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
165  ),
166  currentParticlei_
167  (
168  this->template getModelProperty<label>
169  (
170  "currentParticlei",
171  -1
172  )
173  )
174 {
175  if (this->parcelBasis_ != InjectionModel<CloudType>::pbFixed)
176  {
178  << "Injector model: " << this->modelName()
179  << " Parcel basis must be set to fixed"
180  << exit(FatalError);
181  }
182 
183  if (!time_.size())
184  {
185  // Clean start
186  initialise();
187  }
188 
189  injectorCells_.setSize(position_.size());
190  injectorTetFaces_.setSize(position_.size());
191  injectorTetPts_.setSize(position_.size());
192 
193  updateMesh();
194 
195  this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
196 }
197 
198 
199 template<class CloudType>
201 (
203 )
204 :
206  cloudName_(im.cloudName_),
207  injectorCells_(im.injectorCells_),
208  injectorTetFaces_(im.injectorTetFaces_),
209  injectorTetPts_(im.injectorTetPts_),
210  time_(im.time_),
211  position_(im.position_),
212  positionOffset_(im.positionOffset_),
213  diameter_(im.diameter_),
214  U_(im.U_),
215  volume_(im.volume_),
216  ignoreOutOfBounds_(im.ignoreOutOfBounds_),
217  currentParticlei_(im.currentParticlei_)
218 {}
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
223 template<class CloudType>
225 {
226  label nRejected = 0;
227 
228  bitSet keep(position_.size(), true);
229 
230  forAll(position_, particlei)
231  {
232  if
233  (
234  !this->findCellAtPosition
235  (
236  injectorCells_[particlei],
237  injectorTetFaces_[particlei],
238  injectorTetPts_[particlei],
239  position_[particlei],
240  !ignoreOutOfBounds_
241  )
242  )
243  {
244  keep.unset(particlei);
245  nRejected++;
246  }
247  }
248 
249  if (nRejected > 0)
250  {
251  inplaceSubset(keep, time_);
252  inplaceSubset(keep, position_);
253  inplaceSubset(keep, diameter_);
254  inplaceSubset(keep, U_);
255  inplaceSubset(keep, volume_);
256  inplaceSubset(keep, injectorCells_);
257  inplaceSubset(keep, injectorTetFaces_);
258  inplaceSubset(keep, injectorTetPts_);
259 
260  Info<< " " << nRejected
261  << " particles ignored, out of bounds" << endl;
262  }
263 }
264 
265 
266 template<class CloudType>
268 {
269  return max(time_);
270 }
271 
272 
273 template<class CloudType>
275 (
276  const scalar time0,
277  const scalar time1
278 )
279 {
280  label nParticles = 0;
281  forAll(time_, particlei)
282  {
283  if ((time_[particlei] >= time0) && (time_[particlei] < time1))
284  {
285  ++nParticles;
286  }
287  }
288 
289  return nParticles;
290 }
291 
292 
293 template<class CloudType>
295 (
296  const scalar time0,
297  const scalar time1
298 )
299 {
300  scalar sumVolume = 0;
301  forAll(time_, particlei)
302  {
303  if ((time_[particlei] >= time0) && (time_[particlei] < time1))
304  {
305  sumVolume += volume_[particlei];
306  }
307  }
308 
309  return sumVolume;
310 }
311 
312 
313 template<class CloudType>
315 (
316  const label parceli,
317  const label nParcels,
318  const scalar time,
319  vector& position,
320  label& cellOwner,
321  label& tetFacei,
322  label& tetPti
323 )
324 {
325  // Note: optimisation - consume particle from lists to reduce storage
326  // as injection proceeds
327 
328  currentParticlei_++;
329 
330  position = position_[currentParticlei_];
331  cellOwner = injectorCells_[currentParticlei_];
332  tetFacei = injectorTetFaces_[currentParticlei_];
333  tetPti = injectorTetPts_[currentParticlei_];
334 }
335 
336 
337 template<class CloudType>
339 (
340  const label parceli,
341  const label,
342  const scalar,
343  typename CloudType::parcelType& parcel
344 )
345 {
346  // Set particle velocity
347  parcel.U() = U_[currentParticlei_];
348 
349  // Set particle diameter
350  parcel.d() = diameter_[currentParticlei_];
351 }
352 
353 
354 template<class CloudType>
356 {
357  return false;
358 }
359 
360 
361 template<class CloudType>
363 (
364  const label
365 )
366 {
367  return true;
368 }
369 
370 
371 template<class CloudType>
373 {
375 
376  if (this->writeTime())
377  {
378  this->setModelProperty("currentParticlei", currentParticlei_);
379  this->setModelProperty("time", time_);
380  this->setModelProperty("position", position_);
381  this->setModelProperty("diameter", diameter_);
382  this->setModelProperty("U", U_);
383  this->setModelProperty("volume", volume_);
384  }
385 }
386 
387 
388 // ************************************************************************* //
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:295
mathematicalConstants.H
Foam::InjectionModel::info
virtual void info(Ostream &os)
Write injection info to stream.
Definition: InjectionModel.C:661
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::accessOp
Definition: UList.H:614
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:355
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:315
Foam::injectedParticleCloud
Definition: injectedParticleCloud.H:53
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
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:350
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
SortableList.H
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:339
Foam::InjectedParticleInjection::validInjection
virtual bool validInjection(const label parceli)
Return flag to identify whether or not injection of parcelI is.
Definition: InjectedParticleInjection.C:363
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
InjectedParticleInjection.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:459
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:121
Foam::InjectedParticleInjection::timeEnd
scalar timeEnd() const
Return the end-of-injection time.
Definition: InjectedParticleInjection.C:267
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:39
Foam::InjectedParticleInjection::info
void info(Ostream &os)
Write injection info to stream.
Definition: InjectedParticleInjection.C:372
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:590
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:145
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
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:464
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
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:275
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:224
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:446
Foam::InjectedParticleInjection::injectorTetFaces_
labelList injectorTetFaces_
List of tetFace label per injector.
Definition: InjectedParticleInjection.H:90