CellZoneInjection.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "CellZoneInjection.H"
32#include "globalIndex.H"
33#include "Pstream.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class CloudType>
39(
40 const labelList& cellZoneCells
41)
42{
43 const fvMesh& mesh = this->owner().mesh();
44 const scalarField& V = mesh.V();
45 const label nCells = cellZoneCells.size();
46 Random& rnd = this->owner().rndGen();
47
48 DynamicList<vector> positions(nCells); // initial size only
49 DynamicList<label> injectorCells(nCells); // initial size only
50 DynamicList<label> injectorTetFaces(nCells); // initial size only
51 DynamicList<label> injectorTetPts(nCells); // initial size only
52
53 scalar newParticlesTotal = 0.0;
54 label addParticlesTotal = 0;
55
56 forAll(cellZoneCells, i)
57 {
58 const label celli = cellZoneCells[i];
59
60 // Calc number of particles to add
61 const scalar newParticles = V[celli]*numberDensity_;
62 newParticlesTotal += newParticles;
63 label addParticles = floor(newParticles);
64 addParticlesTotal += addParticles;
65
66 const scalar diff = newParticlesTotal - addParticlesTotal;
67 if (diff > 1)
68 {
69 label corr = floor(diff);
70 addParticles += corr;
71 addParticlesTotal += corr;
72 }
73
74 // Construct cell tet indices
75 const List<tetIndices> cellTetIs =
76 polyMeshTetDecomposition::cellTetIndices(mesh, celli);
77
78 // Construct cell tet volume fractions
79 scalarList cTetVFrac(cellTetIs.size(), Zero);
80 for (label tetI = 1; tetI < cellTetIs.size() - 1; tetI++)
81 {
82 cTetVFrac[tetI] =
83 cTetVFrac[tetI-1] + cellTetIs[tetI].tet(mesh).mag()/V[celli];
84 }
85 cTetVFrac.last() = 1.0;
86
87 // Set new particle position and cellId
88 for (label pI = 0; pI < addParticles; pI++)
89 {
90 const scalar volFrac = rnd.sample01<scalar>();
91 label tetI = 0;
92 forAll(cTetVFrac, vfI)
93 {
94 if (cTetVFrac[vfI] > volFrac)
95 {
96 tetI = vfI;
97 break;
98 }
99 }
100 positions.append(cellTetIs[tetI].tet(mesh).randomPoint(rnd));
101
102 injectorCells.append(celli);
103 injectorTetFaces.append(cellTetIs[tetI].face());
104 injectorTetPts.append(cellTetIs[tetI].tetPt());
105 }
106 }
107
108 // Parallel operation manipulations
109 globalIndex globalPositions(positions.size());
110 List<vector> allPositions(globalPositions.totalSize(), point::max);
111 List<label> allInjectorCells(globalPositions.totalSize(), -1);
112 List<label> allInjectorTetFaces(globalPositions.totalSize(), -1);
113 List<label> allInjectorTetPts(globalPositions.totalSize(), -1);
114
115 // Gather all positions on to all processors
116 SubList<vector>
117 (
118 allPositions,
119 globalPositions.localSize(Pstream::myProcNo()),
120 globalPositions.localStart(Pstream::myProcNo())
121 ) = positions;
122
123 Pstream::listCombineAllGather(allPositions, minEqOp<point>());
124
125 // Gather local cell tet and tet-point Ids, but leave non-local ids set -1
126 SubList<label>
127 (
128 allInjectorCells,
129 globalPositions.localSize(Pstream::myProcNo()),
130 globalPositions.localStart(Pstream::myProcNo())
131 ) = injectorCells;
132 SubList<label>
133 (
134 allInjectorTetFaces,
135 globalPositions.localSize(Pstream::myProcNo()),
136 globalPositions.localStart(Pstream::myProcNo())
137 ) = injectorTetFaces;
138 SubList<label>
139 (
140 allInjectorTetPts,
141 globalPositions.localSize(Pstream::myProcNo()),
142 globalPositions.localStart(Pstream::myProcNo())
143 ) = injectorTetPts;
144
145 // Transfer data
146 positions_.transfer(allPositions);
147 injectorCells_.transfer(allInjectorCells);
148 injectorTetFaces_.transfer(allInjectorTetFaces);
149 injectorTetPts_.transfer(allInjectorTetPts);
150
151
152 if (debug)
153 {
154 OFstream points("points.obj");
155 forAll(positions_, i)
156 {
157 meshTools::writeOBJ(points, positions_[i]);
158 }
159 }
160}
161
162
163// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
164
165template<class CloudType>
167(
168 const dictionary& dict,
169 CloudType& owner,
170 const word& modelName
171)
172:
173 InjectionModel<CloudType>(dict, owner, modelName, typeName),
174 cellZoneName_(this->coeffDict().lookup("cellZone")),
175 numberDensity_(this->coeffDict().getScalar("numberDensity")),
176 positions_(),
177 injectorCells_(),
178 injectorTetFaces_(),
179 injectorTetPts_(),
180 diameters_(),
181 U0_(this->coeffDict().lookup("U0")),
182 sizeDistribution_
183 (
185 (
186 this->coeffDict().subDict("sizeDistribution"), owner.rndGen()
187 )
188 )
189{
190 updateMesh();
191}
192
193
194template<class CloudType>
196(
198)
199:
201 cellZoneName_(im.cellZoneName_),
202 numberDensity_(im.numberDensity_),
203 positions_(im.positions_),
204 injectorCells_(im.injectorCells_),
205 injectorTetFaces_(im.injectorTetFaces_),
206 injectorTetPts_(im.injectorTetPts_),
207 diameters_(im.diameters_),
208 U0_(im.U0_),
209 sizeDistribution_(im.sizeDistribution_.clone())
210{}
211
212
213// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
214
215template<class CloudType>
217{}
218
219
220// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221
222template<class CloudType>
224{
225 // Set/cache the injector cells
226 const fvMesh& mesh = this->owner().mesh();
227 const label zoneI = mesh.cellZones().findZoneID(cellZoneName_);
228
229 if (zoneI < 0)
230 {
232 << "Unknown cell zone name: " << cellZoneName_
233 << ". Valid cell zones are: " << mesh.cellZones().names()
234 << nl << exit(FatalError);
235 }
236
237 const labelList& cellZoneCells = mesh.cellZones()[zoneI];
238 const label nCells = cellZoneCells.size();
239 const scalar nCellsTotal = returnReduce(nCells, sumOp<label>());
240 const scalar VCells = sum(scalarField(mesh.V(), cellZoneCells));
241 const scalar VCellsTotal = returnReduce(VCells, sumOp<scalar>());
242 Info<< " cell zone size = " << nCellsTotal << endl;
243 Info<< " cell zone volume = " << VCellsTotal << endl;
244
245 if ((nCellsTotal == 0) || (VCellsTotal*numberDensity_ < 1))
246 {
248 << "Number of particles to be added to cellZone " << cellZoneName_
249 << " is zero" << endl;
250 }
251 else
252 {
253 setPositions(cellZoneCells);
254
255 Info<< " number density = " << numberDensity_ << nl
256 << " number of particles = " << positions_.size() << endl;
257
258 // Construct parcel diameters
259 diameters_.setSize(positions_.size());
260 forAll(diameters_, i)
261 {
262 diameters_[i] = sizeDistribution_->sample();
263 }
264 }
265
266 // Determine volume of particles to inject
267 this->volumeTotal_ = sum(pow3(diameters_))*constant::mathematical::pi/6.0;
268}
269
270
271template<class CloudType>
273{
274 // Not used
275 return this->SOI_;
276}
277
278
279template<class CloudType>
281(
282 const scalar time0,
283 const scalar time1
284)
285{
286 if ((0.0 >= time0) && (0.0 < time1))
287 {
288 return positions_.size();
289 }
290
291 return 0;
292}
293
294
295template<class CloudType>
297(
298 const scalar time0,
299 const scalar time1
300)
301{
302 // All parcels introduced at SOI
303 if ((0.0 >= time0) && (0.0 < time1))
304 {
305 return this->volumeTotal_;
306 }
307
308 return 0.0;
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 position = positions_[parcelI];
325 cellOwner = injectorCells_[parcelI];
326 tetFacei = injectorTetFaces_[parcelI];
327 tetPti = injectorTetPts_[parcelI];
328}
329
330
331template<class CloudType>
333(
334 const label parcelI,
335 const label,
336 const scalar,
337 typename CloudType::parcelType& parcel
338)
339{
340 // set particle velocity
341 parcel.U() = U0_;
342
343 // set particle diameter
344 parcel.d() = diameters_[parcelI];
345}
346
347
348template<class CloudType>
350{
351 return false;
352}
353
354
355template<class CloudType>
357{
358 return true;
359}
360
361
362// ************************************************************************* //
Injection positions specified by a particle number density within a cell set.
scalar volumeToInject(const scalar time0, const scalar time1)
Volume of parcels to introduce relative to SOI.
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 bool validInjection(const label parcelI)
Return flag to identify whether or not injection of parcelI is.
virtual void setProperties(const label parcelI, const label nParcels, const scalar time, typename CloudType::parcelType &parcel)
Set the parcel properties.
virtual ~CellZoneInjection()
Destructor.
virtual void updateMesh()
Set injector locations when mesh is updated.
virtual bool fullyDescribed() const
Flag to identify whether model fully describes the parcel.
scalar timeEnd() const
Return the end-of-injection time.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Templated injection model class.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:525
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:304
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar pi(M_PI)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23