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-------------------------------------------------------------------------------
10License
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
30#include "bitSet.H"
32
33using namespace Foam::constant;
34
35// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36
37template<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 {
63 List<List<scalar>> procTime(Pstream::nProcs());
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
83 List<List<scalar>> procD(Pstream::nProcs());
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
93 List<List<vector>> procU(Pstream::nProcs());
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
142template<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{
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
191
192 updateMesh();
193
194 this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
195}
196
197
198template<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
222template<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
265template<class CloudType>
267{
268 return max(time_);
269}
270
271
272template<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
292template<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
312template<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
336template<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
353template<class CloudType>
355{
356 return false;
357}
358
359
360template<class CloudType>
362(
363 const label
364)
365{
366 return true;
367}
368
369
370template<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// ************************************************************************* //
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition: DSMCCloudI.H:98
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Replays an set of particle data based on an injectedParticleCloud, using the assumption of one partic...
void initialise()
Initialise injectors.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
labelList injectorCells_
List of cell label per injector.
virtual label parcelsToInject(const scalar time0, const scalar time1)
Number of parcels to introduce relative to SOI.
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.
virtual void setProperties(const label parceli, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
labelList injectorTetPts_
List of tetPt label per injector.
vectorList position_
List of position per particle [m].
labelList injectorTetFaces_
List of tetFace label per injector.
scalarList time_
List of injection time per particle [s].
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool validInjection(const label parceli)
Return flag to identify whether or not injection of parcelI is.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Templated injection model class.
parcelBasis parcelBasis_
Parcel basis enumeration.
scalar massTotal_
Total mass to inject [kg].
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
bitSet & unset(const bitSet &other)
Definition: bitSetI.H:628
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Lookup type of boundary radiation properties.
Definition: lookup.H:66
const word & modelName() const
Return const access to the name of the sub-model.
Definition: subModelBase.C:107
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
constexpr scalar pi(M_PI)
Different types of constants.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333