Loading...
Searching...
No Matches
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-2022 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{
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::allGatherList(procTime);
66 time =
68 (
69 procTime, accessOp<List<scalar>>()
70 );
71
72 List<List<point>> procPosition(Pstream::nProcs());
73 procPosition[Pstream::myProcNo()].transfer(position);
74 Pstream::allGatherList(procPosition);
75 position =
77 (
78 procPosition, accessOp<List<point>>()
79 );
80
82 procD[Pstream::myProcNo()].transfer(diameter);
84 diameter =
86 (
87 procD, accessOp<List<scalar>>()
88 );
89
91 procU[Pstream::myProcNo()].transfer(U);
93 U =
95 (
96 procU, accessOp<List<vector>>()
97 );
98 }
99
100 nParticles = time.size();
101
102 // Reset SOI according to user selection
103 scalar minTime = min(time);
104 forAll(time, i)
105 {
106 time[i] -= minTime;
107 }
108
109 // Sort and renumber to ensure lists in ascending time
110 labelList sortedIndices(Foam::sortedOrder(time));
111
112 time_ = UIndirectList<scalar>(time, sortedIndices);
113 position_ = UIndirectList<point>(position, sortedIndices);
114 diameter_ = UIndirectList<scalar>(diameter, sortedIndices);
115 U_ = UIndirectList<vector>(U, sortedIndices);
116
117 // Pre-calculate injected particle volumes
118 List<scalar> volume(nParticles);
119 scalar sumVolume = 0;
120 forAll(volume, i)
121 {
122 scalar vol = pow3(diameter_[i])*mathematical::pi/6.0;
123 volume[i] = vol;
124 sumVolume += vol;
125 }
126 volume_.transfer(volume);
127
128 // Set the volume of particles to inject
129 this->volumeTotal_ = sumVolume;
130
131 // Provide some feedback
132 Info<< " Read " << nParticles << " particles" << endl;
133}
134
135
136// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137
138template<class CloudType>
140(
141 const dictionary& dict,
142 CloudType& owner,
143 const word& modelName
144)
145:
146 InjectionModel<CloudType>(dict, owner, modelName, typeName),
147 cloudName_(this->coeffDict().lookup("cloud")),
148 injectorCells_(),
149 injectorTetFaces_(),
150 injectorTetPts_(),
151 time_(this->template getModelProperty<scalarList>("time")),
152 position_(this->template getModelProperty<vectorList>("position")),
153 positionOffset_(this->coeffDict().lookup("positionOffset")),
154 diameter_(this->template getModelProperty<scalarList>("diameter")),
155 U_(this->template getModelProperty<vectorList>("U")),
156 volume_(this->template getModelProperty<scalarList>("volume")),
157 ignoreOutOfBounds_
158 (
159 this->coeffDict().getOrDefault("ignoreOutOfBounds", false)
160 ),
161 currentParticlei_
162 (
163 this->template getModelProperty<label>
164 (
165 "currentParticlei",
166 -1
167 )
168 )
169{
171 {
173 << "Injector model: " << this->modelName()
174 << " Parcel basis must be set to fixed"
175 << exit(FatalError);
176 }
177
178 if (!time_.size())
179 {
180 // Clean start
181 initialise();
182 }
183
187
189
190 this->massTotal_ = this->volumeTotal_*this->owner().constProps().rho0();
191}
192
193
194template<class CloudType>
196(
198)
199:
201 cloudName_(im.cloudName_),
202 injectorCells_(im.injectorCells_),
203 injectorTetFaces_(im.injectorTetFaces_),
204 injectorTetPts_(im.injectorTetPts_),
205 time_(im.time_),
206 position_(im.position_),
207 positionOffset_(im.positionOffset_),
208 diameter_(im.diameter_),
209 U_(im.U_),
210 volume_(im.volume_),
211 ignoreOutOfBounds_(im.ignoreOutOfBounds_),
213{}
214
215
216// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217
218template<class CloudType>
220{
221 label nRejected = 0;
222
223 bitSet keep(position_.size(), true);
224
225 forAll(position_, particlei)
226 {
227 if
228 (
229 !this->findCellAtPosition
230 (
231 injectorCells_[particlei],
232 injectorTetFaces_[particlei],
233 injectorTetPts_[particlei],
234 position_[particlei],
235 !ignoreOutOfBounds_
236 )
237 )
238 {
239 keep.unset(particlei);
240 nRejected++;
241 }
242 }
243
244 if (nRejected > 0)
245 {
246 inplaceSubset(keep, time_);
247 inplaceSubset(keep, position_);
248 inplaceSubset(keep, diameter_);
249 inplaceSubset(keep, U_);
250 inplaceSubset(keep, volume_);
251 inplaceSubset(keep, injectorCells_);
252 inplaceSubset(keep, injectorTetFaces_);
253 inplaceSubset(keep, injectorTetPts_);
254
255 Info<< " " << nRejected
256 << " particles ignored, out of bounds" << endl;
257 }
258}
259
260
261template<class CloudType>
263{
264 return max(time_);
265}
266
267
268template<class CloudType>
270(
271 const scalar time0,
272 const scalar time1
273)
274{
275 label nParticles = 0;
276 forAll(time_, particlei)
277 {
278 if ((time_[particlei] >= time0) && (time_[particlei] < time1))
279 {
280 ++nParticles;
281 }
283
284 return nParticles;
285}
286
287
288template<class CloudType>
290(
291 const scalar time0,
292 const scalar time1
293)
294{
295 scalar sumVolume = 0;
296 forAll(time_, particlei)
297 {
298 if ((time_[particlei] >= time0) && (time_[particlei] < time1))
299 {
300 sumVolume += volume_[particlei];
301 }
303
304 return sumVolume;
305}
306
307
308template<class CloudType>
310(
311 const label parceli,
312 const label nParcels,
313 const scalar time,
314 vector& position,
315 label& cellOwner,
316 label& tetFacei,
317 label& tetPti
318)
319{
320 // Note: optimisation - consume particle from lists to reduce storage
321 // as injection proceeds
322
323 currentParticlei_++;
324
325 position = position_[currentParticlei_];
329}
330
331
332template<class CloudType>
334(
335 const label parceli,
336 const label,
337 const scalar,
338 typename CloudType::parcelType& parcel
339)
340{
341 // Set particle velocity
342 parcel.U() = U_[currentParticlei_];
344 // Set particle diameter
345 parcel.d() = diameter_[currentParticlei_];
346}
347
348
349template<class CloudType>
351{
352 return false;
353}
354
355
356template<class CloudType>
358(
359 const label
361{
362 return true;
363}
364
365
366template<class CloudType>
368{
370
371 if (this->writeTime())
372 {
373 this->setModelProperty("currentParticlei", currentParticlei_);
374 this->setModelProperty("time", time_);
375 this->setModelProperty("position", position_);
376 this->setModelProperty("diameter", diameter_);
377 this->setModelProperty("U", U_);
378 this->setModelProperty("volume", volume_);
379 }
380}
381
382
383// ************************************************************************* //
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
const List< typename ParcelType::constantProperties > & constProps() const
Return all of the constant properties.
Definition DSMCCloudI.H:91
const vector & U() const
Return const access to velocity.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
Replays an set of particle data based on an injectedParticleCloud, using the assumption of one partic...
void initialise()
Initialise injectors.
const word cloudName_
Name of cloud used to seed the new particles.
virtual scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
vectorList U_
List of velocity per particle [m/s].
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.
vector positionOffset_
Position offset to apply to input positions.
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.
label currentParticlei_
Index of current particle.
vectorList position_
List of position per particle [m].
Switch ignoreOutOfBounds_
Flag to suppress errors if particle injection site is out-of-bounds.
labelList injectorTetFaces_
List of tetFace label per injector.
scalarList time_
List of injection time per particle [s].
InjectedParticleInjection(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
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.
scalarList volume_
List of volume per particle [m3].
scalarList diameter_
List of diameter per particle [m].
scalar timeEnd() const
Return the end-of-injection time.
Templated injection model class.
parcelBasis parcelBasis_
Parcel basis enumeration.
virtual void info()
Write injection info.
scalar volumeTotal_
Total volume of particles introduced by this injector [m^3] Note: scaled to ensure massTotal is achie...
InjectionModel(CloudType &owner)
Construct null from owner.
scalar massTotal_
Total mass to inject [kg].
virtual bool findCellAtPosition(label &celli, label &tetFacei, label &tetPti, vector &position, bool errorOnNotFound=true)
Find the cell that contains the supplied position.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
bitSet & unset(const bitSet &other)
Unset (subtract) the bits specified in the other bitset, which is a set difference corresponds to the...
Definition bitSetI.H:540
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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:60
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition ListListOps.C:62
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:40
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
DSMCCloud< dsmcParcel > CloudType
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< vector > vectorList
List of vector.
Definition vectorList.H:32
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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:26
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Object access operator or list access operator (default is pass-through).
Definition UList.H:1120