ParticleCollector.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "ParticleCollector.H"
27 #include "Pstream.H"
28 #include "surfaceWriter.H"
29 #include "randomGenerator.H"
30 #include "triangle.H"
31 #include "cloud.H"
32 #include "axesRotation.H"
33 #include "OSspecific.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class CloudType>
39 (
40  const faceList& faces,
41  const Field<point>& points,
42  const Field<scalar>& area
43 )
44 {
45  // Create the output file if not already created
46  if (log_)
47  {
48  if (debug)
49  {
50  Info<< "Creating output file" << endl;
51  }
52 
53  if (Pstream::master())
54  {
55  // Create directory if does not exist
56  mkDir(this->writeTimeDir());
57 
58  // Open new file at start up
59  outputFilePtr_.reset
60  (
61  new OFstream(this->writeTimeDir()/(type() + ".dat"))
62  );
63 
64  outputFilePtr_()
65  << "# Source : " << type() << nl
66  << "# Bins : " << faces.size() << nl
67  << "# Total area : " << sum(area) << nl;
68 
69  outputFilePtr_()
70  << "# Geometry :" << nl
71  << '#'
72  << tab << "Bin"
73  << tab << "(Centre_x Centre_y Centre_z)"
74  << tab << "Area"
75  << nl;
76 
77  forAll(faces, i)
78  {
79  outputFilePtr_()
80  << '#'
81  << tab << i
82  << tab << faces[i].centre(points)
83  << tab << area[i]
84  << nl;
85  }
86 
87  outputFilePtr_()
88  << '#' << nl
89  << "# Output format:" << nl;
90 
91  forAll(faces, i)
92  {
93  word id = Foam::name(i);
94  word binId = "bin_" + id;
95 
96  outputFilePtr_()
97  << '#'
98  << tab << "Time"
99  << tab << binId
100  << tab << "mass[" << id << "]"
101  << tab << "massFlowRate[" << id << "]"
102  << endl;
103  }
104  }
105  }
106 }
107 
108 
109 template<class CloudType>
111 (
112  const List<Field<point>>& polygons
113 )
114 {
115  mode_ = mtPolygon;
116 
117  label nPoints = 0;
118  forAll(polygons, polyI)
119  {
120  label np = polygons[polyI].size();
121  if (np < 3)
122  {
123  FatalIOErrorInFunction(this->coeffDict())
124  << "polygons must consist of at least 3 points"
125  << exit(FatalIOError);
126  }
127 
128  nPoints += np;
129  }
130 
131  label pointOffset = 0;
132  points_.setSize(nPoints);
133  faces_.setSize(polygons.size());
134  area_.setSize(polygons.size());
135  forAll(faces_, facei)
136  {
137  const Field<point>& polyPoints = polygons[facei];
138  face f(identityMap(polyPoints.size()) + pointOffset);
139  UIndirectList<point>(points_, f) = polyPoints;
140  area_[facei] = f.mag(points_);
141  faces_[facei].transfer(f);
142 
143  pointOffset += polyPoints.size();
144  }
145 }
146 
147 
148 template<class CloudType>
150 {
151  mode_ = mtConcentricCircle;
152 
153  vector origin(this->coeffDict().lookup("origin"));
154 
155  this->coeffDict().lookup("radius") >> radius_;
156  nSector_ = this->coeffDict().template lookup<label>("nSector");
157 
158  label nS = nSector_;
159 
160  vector refDir;
161  if (nSector_ > 1)
162  {
163  refDir = this->coeffDict().lookup("refDir");
164  refDir -= normal_[0]*(normal_[0] & refDir);
165  refDir /= mag(refDir);
166  }
167  else
168  {
169  // Set 4 quadrants for single sector cases
170  nS = 4;
171  refDir = normalised(perpendicular(normal_[0]));
172  }
173 
174  scalar dTheta = degToRad(5.0);
175  scalar dThetaSector = degToRad(360.0)/scalar(nS);
176  label intervalPerSector = max(1, ceil(dThetaSector/dTheta));
177  dTheta = dThetaSector/scalar(intervalPerSector);
178 
179  label nPointPerSector = intervalPerSector + 1;
180 
181  label nPointPerRadius = nS*(nPointPerSector - 1);
182  label nPoint = radius_.size()*nPointPerRadius;
183  label nFace = radius_.size()*nS;
184 
185  // Add origin
186  nPoint++;
187 
188  points_.setSize(nPoint);
189  faces_.setSize(nFace);
190  area_.setSize(nFace);
191 
192  coordSys_ =
193  coordinateSystems::cylindrical("coordSys", origin, normal_[0], refDir);
194 
195  List<label> ptIDs(identityMap(nPointPerRadius));
196 
197  points_[0] = origin;
198 
199  // Points
200  forAll(radius_, radI)
201  {
202  label pointOffset = radI*nPointPerRadius + 1;
203 
204  for (label i = 0; i < nPointPerRadius; i++)
205  {
206  label pI = i + pointOffset;
207  point pCyl(radius_[radI], i*dTheta, 0.0);
208  points_[pI] = coordSys_.globalPosition(pCyl);
209  }
210  }
211 
212  // Faces
213  DynamicList<label> facePts(2*nPointPerSector);
214  forAll(radius_, radI)
215  {
216  if (radI == 0)
217  {
218  for (label secI = 0; secI < nS; secI++)
219  {
220  facePts.clear();
221 
222  // Append origin point
223  facePts.append(0);
224 
225  for (label ptI = 0; ptI < nPointPerSector; ptI++)
226  {
227  label i = ptI + secI*(nPointPerSector - 1);
228  label id = ptIDs.fcIndex(i - 1) + 1;
229  facePts.append(id);
230  }
231 
232  label facei = secI + radI*nS;
233 
234  faces_[facei] = face(facePts);
235  area_[facei] = faces_[facei].mag(points_);
236  }
237  }
238  else
239  {
240  for (label secI = 0; secI < nS; secI++)
241  {
242  facePts.clear();
243 
244  label offset = (radI - 1)*nPointPerRadius + 1;
245 
246  for (label ptI = 0; ptI < nPointPerSector; ptI++)
247  {
248  label i = ptI + secI*(nPointPerSector - 1);
249  label id = offset + ptIDs.fcIndex(i - 1);
250  facePts.append(id);
251  }
252  for (label ptI = nPointPerSector-1; ptI >= 0; ptI--)
253  {
254  label i = ptI + secI*(nPointPerSector - 1);
255  label id = offset + nPointPerRadius + ptIDs.fcIndex(i - 1);
256  facePts.append(id);
257  }
258 
259  label facei = secI + radI*nS;
260 
261  faces_[facei] = face(facePts);
262  area_[facei] = faces_[facei].mag(points_);
263  }
264  }
265  }
266 }
267 
268 
269 template<class CloudType>
271 (
272  const point& p1,
273  const point& p2
274 ) const
275 {
276  forAll(faces_, facei)
277  {
278  const label facePoint0 = faces_[facei][0];
279 
280  const point& pf = points_[facePoint0];
281 
282  const scalar d1 = normal_[facei] & (p1 - pf);
283  const scalar d2 = normal_[facei] & (p2 - pf);
284 
285  if (sign(d1) == sign(d2))
286  {
287  // Did not cross polygon plane
288  continue;
289  }
290 
291  // Intersection point
292  const point pIntersect = p1 + (d1/(d1 - d2))*(p2 - p1);
293 
294  // Identify if point is within the bounds of the face. Create triangles
295  // between the intersection point and each edge of the face. If all the
296  // triangle normals point in the same direction as the face normal, then
297  // the particle is within the face. Note that testing for pointHits on
298  // the face's decomposed triangles does not work due to ambiguity along
299  // the diagonals.
300  const face& f = faces_[facei];
301  const vector a = f.area(points_);
302  bool inside = true;
303  for (label i = 0; i < f.size(); ++ i)
304  {
305  const label j = f.fcIndex(i);
306  const triPointRef t(pIntersect, points_[f[i]], points_[f[j]]);
307  if ((a & t.area()) < 0)
308  {
309  inside = false;
310  break;
311  }
312  }
313 
314  // Add to the list of hits
315  if (inside)
316  {
317  hitFaceIndices_.append(facei);
318  }
319  }
320 }
321 
322 
323 template<class CloudType>
325 (
326  const point& p1,
327  const point& p2
328 ) const
329 {
330  label secI = -1;
331 
332  const scalar d1 = normal_[0] & (p1 - coordSys_.origin());
333  const scalar d2 = normal_[0] & (p2 - coordSys_.origin());
334 
335  if (sign(d1) == sign(d2))
336  {
337  // Did not cross plane
338  return;
339  }
340 
341  // Intersection point in cylindrical co-ordinate system
342  const point pCyl = coordSys_.localPosition(p1 + (d1/(d1 - d2))*(p2 - p1));
343 
344  scalar r = pCyl[0];
345 
346  if (r < radius_.last())
347  {
348  label radI = 0;
349  while (r > radius_[radI])
350  {
351  radI++;
352  }
353 
354  if (nSector_ == 1)
355  {
356  secI = 4*radI;
357  }
358  else
359  {
360  scalar theta = pCyl[1] + constant::mathematical::pi;
361 
362  secI =
363  nSector_*radI
364  + floor
365  (
366  scalar(nSector_)*theta/constant::mathematical::twoPi
367  );
368  }
369  }
370 
371  if (secI != -1)
372  {
373  hitFaceIndices_.append(secI);
374  }
375 }
376 
377 
378 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
379 
380 template<class CloudType>
382 {
383  const fvMesh& mesh = this->owner().mesh();
384  const Time& time = mesh.time();
385  scalar timeNew = time.value();
386  scalar timeElapsed = timeNew - timeOld_;
387 
388  totalTime_ += timeElapsed;
389 
390  const scalar alpha = (totalTime_ - timeElapsed)/totalTime_;
391  const scalar beta = timeElapsed/totalTime_;
392 
393  forAll(faces_, facei)
394  {
395  massFlowRate_[facei] =
396  alpha*massFlowRate_[facei] + beta*mass_[facei]/timeElapsed;
397  massTotal_[facei] += mass_[facei];
398  }
399 
400  const label proci = Pstream::myProcNo();
401 
402  Info<< type() << " output:" << nl;
403 
404  Field<scalar> faceMassTotal(mass_.size(), 0.0);
405  this->getModelProperty("massTotal", faceMassTotal);
406 
407  Field<scalar> faceMassFlowRate(massFlowRate_.size(), 0.0);
408  this->getModelProperty("massFlowRate", faceMassFlowRate);
409 
410 
411  scalar sumTotalMass = 0.0;
412  scalar sumAverageMFR = 0.0;
413  forAll(faces_, facei)
414  {
415  scalarList allProcMass(Pstream::nProcs());
416  allProcMass[proci] = massTotal_[facei];
417  Pstream::gatherList(allProcMass);
418  faceMassTotal[facei] += sum(allProcMass);
419 
420  scalarList allProcMassFlowRate(Pstream::nProcs());
421  allProcMassFlowRate[proci] = massFlowRate_[facei];
422  Pstream::gatherList(allProcMassFlowRate);
423  faceMassFlowRate[facei] += sum(allProcMassFlowRate);
424 
425  sumTotalMass += faceMassTotal[facei];
426  sumAverageMFR += faceMassFlowRate[facei];
427 
428  if (outputFilePtr_.valid())
429  {
430  outputFilePtr_()
431  << time.name()
432  << tab << facei
433  << tab << faceMassTotal[facei]
434  << tab << faceMassFlowRate[facei]
435  << endl;
436  }
437  }
438 
439  Info<< " sum(total mass) = " << sumTotalMass << nl
440  << " sum(average mass flow rate) = " << sumAverageMFR << nl
441  << endl;
442 
443 
444  if (surfaceFormat_ != "none")
445  {
446  if (Pstream::master())
447  {
449  (
450  surfaceWriter::New(surfaceFormat_, this->coeffDict())
451  );
452 
453  writer->write
454  (
455  this->writeTimeDir(),
456  "collector",
457  points_,
458  faces_,
459  false,
460  "massTotal",
461  faceMassTotal,
462  "massFlowRate",
463  faceMassFlowRate
464  );
465  }
466  }
467 
468 
469  if (resetOnWrite_)
470  {
471  Field<scalar> dummy(faceMassTotal.size(), 0.0);
472  this->setModelProperty("massTotal", dummy);
473  this->setModelProperty("massFlowRate", dummy);
474 
475  timeOld_ = timeNew;
476  totalTime_ = 0.0;
477  }
478  else
479  {
480  this->setModelProperty("massTotal", faceMassTotal);
481  this->setModelProperty("massFlowRate", faceMassFlowRate);
482  }
483 
484  forAll(faces_, facei)
485  {
486  mass_[facei] = 0.0;
487  massTotal_[facei] = 0.0;
488  massFlowRate_[facei] = 0.0;
489  }
490 }
491 
492 
493 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
494 
495 template<class CloudType>
497 (
498  const dictionary& dict,
499  CloudType& owner,
500  const word& modelName
501 )
502 :
503  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
504  mode_(mtUnknown),
505  parcelType_(this->coeffDict().lookupOrDefault("parcelType", -1)),
506  removeCollected_(this->coeffDict().lookup("removeCollected")),
507  points_(),
508  faces_(),
509  nSector_(0),
510  radius_(),
511  coordSys_("coordSys", vector::zero, axesRotation(sphericalTensor::I)),
512  normal_(),
513  negateParcelsOppositeNormal_
514  (
515  readBool(this->coeffDict().lookup("negateParcelsOppositeNormal"))
516  ),
517  surfaceFormat_(this->coeffDict().lookup("surfaceFormat")),
518  resetOnWrite_(this->coeffDict().lookup("resetOnWrite")),
519  totalTime_(0.0),
520  mass_(),
521  massTotal_(),
522  massFlowRate_(),
523  log_(this->coeffDict().lookup("log")),
524  outputFilePtr_(),
525  timeOld_(owner.mesh().time().value()),
526  hitFaceIndices_()
527 {
528  normal_ /= mag(normal_);
529 
530  word mode(this->coeffDict().lookup("mode"));
531  if (mode == "polygon")
532  {
533  List<Field<point>> polygons(this->coeffDict().lookup("polygons"));
534 
535  initPolygons(polygons);
536 
537  vector n0(this->coeffDict().lookup("normal"));
538  normal_ = vectorField(faces_.size(), n0);
539  }
540  else if (mode == "polygonWithNormal")
541  {
542  List<Tuple2<Field<point>, vector>> polygonAndNormal
543  (
544  this->coeffDict().lookup("polygons")
545  );
546 
547  List<Field<point>> polygons(polygonAndNormal.size());
548  normal_.setSize(polygonAndNormal.size());
549 
550  forAll(polygons, polyI)
551  {
552  polygons[polyI] = polygonAndNormal[polyI].first();
553  normal_[polyI] = polygonAndNormal[polyI].second();
554  normal_[polyI] /= mag(normal_[polyI]) + rootVSmall;
555  }
556 
557  initPolygons(polygons);
558  }
559  else if (mode == "concentricCircle")
560  {
561  vector n0(this->coeffDict().lookup("normal"));
562  normal_ = vectorField(1, n0);
563 
564  initConcentricCircles();
565  }
566  else
567  {
569  << "Unknown mode " << mode << ". Available options are "
570  << "polygon, polygonWithNormal and concentricCircle"
571  << exit(FatalIOError);
572  }
573 
574  mass_.setSize(faces_.size(), 0.0);
575  massTotal_.setSize(faces_.size(), 0.0);
576  massFlowRate_.setSize(faces_.size(), 0.0);
577 
578  makeLogFile(faces_, points_, area_);
579 }
580 
581 
582 template<class CloudType>
584 (
586 )
587 :
589  mode_(pc.mode_),
590  parcelType_(pc.parcelType_),
591  removeCollected_(pc.removeCollected_),
592  points_(pc.points_),
593  faces_(pc.faces_),
594  nSector_(pc.nSector_),
595  radius_(pc.radius_),
596  coordSys_(pc.coordSys_),
597  normal_(pc.normal_),
598  negateParcelsOppositeNormal_(pc.negateParcelsOppositeNormal_),
599  surfaceFormat_(pc.surfaceFormat_),
600  resetOnWrite_(pc.resetOnWrite_),
601  totalTime_(pc.totalTime_),
602  mass_(pc.mass_),
603  massTotal_(pc.massTotal_),
604  massFlowRate_(pc.massFlowRate_),
605  log_(pc.log_),
606  outputFilePtr_(),
607  timeOld_(0.0),
608  hitFaceIndices_()
609 {}
610 
611 
612 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
613 
614 template<class CloudType>
616 {}
617 
618 
619 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
620 
621 template<class CloudType>
623 (
624  parcelType& p,
625  const scalar dt,
626  const point& position0,
627  bool& keepParticle
628 )
629 {
630  if ((parcelType_ != -1) && (parcelType_ != p.typeId()))
631  {
632  return;
633  }
634 
635  hitFaceIndices_.clear();
636 
637  switch (mode_)
638  {
639  case mtPolygon:
640  {
641  collectParcelPolygon
642  (
643  position0,
644  p.position(this->owner().mesh())
645  );
646  break;
647  }
648  case mtConcentricCircle:
649  {
650  collectParcelConcentricCircles
651  (
652  position0,
653  p.position(this->owner().mesh())
654  );
655  break;
656  }
657  default:
658  {}
659  }
660 
661  forAll(hitFaceIndices_, i)
662  {
663  label facei = hitFaceIndices_[i];
664  scalar m = p.nParticle()*p.mass();
665 
666  if (negateParcelsOppositeNormal_)
667  {
668  scalar Unormal = 0;
669  vector Uhat = p.U();
670  switch (mode_)
671  {
672  case mtPolygon:
673  {
674  Unormal = Uhat & normal_[facei];
675  break;
676  }
677  case mtConcentricCircle:
678  {
679  Unormal = Uhat & normal_[0];
680  break;
681  }
682  default:
683  {}
684  }
685 
686  Uhat /= mag(Uhat) + rootVSmall;
687 
688  if (Unormal < 0)
689  {
690  m = -m;
691  }
692  }
693 
694  // Add mass contribution
695  mass_[facei] += m;
696 
697  if (nSector_ == 1)
698  {
699  mass_[facei + 1] += m;
700  mass_[facei + 2] += m;
701  mass_[facei + 3] += m;
702  }
703 
704  if (removeCollected_)
705  {
706  keepParticle = false;
707  }
708  }
709 }
710 
711 
712 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated cloud function object base class.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Function object to collect the parcel mass- and mass flow rate over a set of polygons....
virtual void postMove(parcelType &p, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
virtual ~ParticleCollector()
Destructor.
ParticleCollector(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
void write()
Write post-processing info.
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
T & first()
Return the first element of the list.
Definition: UListI.H:114
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A coordinate rotation specified using global axis.
Definition: axesRotation.H:67
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
A class for handling words, derived from string.
Definition: word.H:62
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const pointField & points
label nPoints
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const scalar twoPi(2 *pi)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool readBool(Istream &)
Definition: boolIO.C:66
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:461
scalar degToRad(const scalar deg)
Convert degrees to radians.
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static const char tab
Definition: Ostream.H:265
messageStream Info
static const Identity< scalar > I
Definition: Identity.H:93
vector point
Point is a vector.
Definition: point.H:41
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:510
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
dimensionSet perpendicular(const dimensionSet &)
Definition: dimensionSet.C:516
void offset(label &lst, const label o)
List< face > faceList
Definition: faceListFwd.H:41
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
dictionary dict
volScalarField & p
mkDir(pdfPath)