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-2021 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 "unitConversion.H"
30 #include "Random.H"
31 #include "triangle.H"
32 #include "cloud.H"
33 #include "axesRotation.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(identity(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 = 5.0;
175  scalar dThetaSector = 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_ = coordinateSystems::cylindrical
193  (
194  "coordSys",
195  origin,
196  normal_[0],
197  refDir,
198  false
199  );
200 
201  List<label> ptIDs(identity(nPointPerRadius));
202 
203  points_[0] = origin;
204 
205  // Points
206  forAll(radius_, radI)
207  {
208  label pointOffset = radI*nPointPerRadius + 1;
209 
210  for (label i = 0; i < nPointPerRadius; i++)
211  {
212  label pI = i + pointOffset;
213  point pCyl(radius_[radI], degToRad(i*dTheta), 0.0);
214  points_[pI] = coordSys_.globalPosition(pCyl);
215  }
216  }
217 
218  // Faces
219  DynamicList<label> facePts(2*nPointPerSector);
220  forAll(radius_, radI)
221  {
222  if (radI == 0)
223  {
224  for (label secI = 0; secI < nS; secI++)
225  {
226  facePts.clear();
227 
228  // Append origin point
229  facePts.append(0);
230 
231  for (label ptI = 0; ptI < nPointPerSector; ptI++)
232  {
233  label i = ptI + secI*(nPointPerSector - 1);
234  label id = ptIDs.fcIndex(i - 1) + 1;
235  facePts.append(id);
236  }
237 
238  label facei = secI + radI*nS;
239 
240  faces_[facei] = face(facePts);
241  area_[facei] = faces_[facei].mag(points_);
242  }
243  }
244  else
245  {
246  for (label secI = 0; secI < nS; secI++)
247  {
248  facePts.clear();
249 
250  label offset = (radI - 1)*nPointPerRadius + 1;
251 
252  for (label ptI = 0; ptI < nPointPerSector; ptI++)
253  {
254  label i = ptI + secI*(nPointPerSector - 1);
255  label id = offset + ptIDs.fcIndex(i - 1);
256  facePts.append(id);
257  }
258  for (label ptI = nPointPerSector-1; ptI >= 0; ptI--)
259  {
260  label i = ptI + secI*(nPointPerSector - 1);
261  label id = offset + nPointPerRadius + ptIDs.fcIndex(i - 1);
262  facePts.append(id);
263  }
264 
265  label facei = secI + radI*nS;
266 
267  faces_[facei] = face(facePts);
268  area_[facei] = faces_[facei].mag(points_);
269  }
270  }
271  }
272 }
273 
274 
275 template<class CloudType>
277 (
278  const point& p1,
279  const point& p2
280 ) const
281 {
282  forAll(faces_, facei)
283  {
284  const label facePoint0 = faces_[facei][0];
285 
286  const point& pf = points_[facePoint0];
287 
288  const scalar d1 = normal_[facei] & (p1 - pf);
289  const scalar d2 = normal_[facei] & (p2 - pf);
290 
291  if (sign(d1) == sign(d2))
292  {
293  // Did not cross polygon plane
294  continue;
295  }
296 
297  // Intersection point
298  const point pIntersect = p1 + (d1/(d1 - d2))*(p2 - p1);
299 
300  // Identify if point is within the bounds of the face. Create triangles
301  // between the intersection point and each edge of the face. If all the
302  // triangle normals point in the same direction as the face normal, then
303  // the particle is within the face. Note that testing for pointHits on
304  // the face's decomposed triangles does not work due to ambiguity along
305  // the diagonals.
306  const face& f = faces_[facei];
307  const vector a = f.area(points_);
308  bool inside = true;
309  for (label i = 0; i < f.size(); ++ i)
310  {
311  const label j = f.fcIndex(i);
312  const triPointRef t(pIntersect, points_[f[i]], points_[f[j]]);
313  if ((a & t.area()) < 0)
314  {
315  inside = false;
316  break;
317  }
318  }
319 
320  // Add to the list of hits
321  if (inside)
322  {
323  hitFaceIDs_.append(facei);
324  }
325  }
326 }
327 
328 
329 template<class CloudType>
331 (
332  const point& p1,
333  const point& p2
334 ) const
335 {
336  label secI = -1;
337 
338  const scalar d1 = normal_[0] & (p1 - coordSys_.origin());
339  const scalar d2 = normal_[0] & (p2 - coordSys_.origin());
340 
341  if (sign(d1) == sign(d2))
342  {
343  // Did not cross plane
344  return;
345  }
346 
347  // Intersection point in cylindrical co-ordinate system
348  const point pCyl = coordSys_.localPosition(p1 + (d1/(d1 - d2))*(p2 - p1));
349 
350  scalar r = pCyl[0];
351 
352  if (r < radius_.last())
353  {
354  label radI = 0;
355  while (r > radius_[radI])
356  {
357  radI++;
358  }
359 
360  if (nSector_ == 1)
361  {
362  secI = 4*radI;
363  }
364  else
365  {
366  scalar theta = pCyl[1] + constant::mathematical::pi;
367 
368  secI =
369  nSector_*radI
370  + floor
371  (
372  scalar(nSector_)*theta/constant::mathematical::twoPi
373  );
374  }
375  }
376 
377  if (secI != -1)
378  {
379  hitFaceIDs_.append(secI);
380  }
381 }
382 
383 
384 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
385 
386 template<class CloudType>
388 {
389  const fvMesh& mesh = this->owner().mesh();
390  const Time& time = mesh.time();
391  scalar timeNew = time.value();
392  scalar timeElapsed = timeNew - timeOld_;
393 
394  totalTime_ += timeElapsed;
395 
396  const scalar alpha = (totalTime_ - timeElapsed)/totalTime_;
397  const scalar beta = timeElapsed/totalTime_;
398 
399  forAll(faces_, facei)
400  {
401  massFlowRate_[facei] =
402  alpha*massFlowRate_[facei] + beta*mass_[facei]/timeElapsed;
403  massTotal_[facei] += mass_[facei];
404  }
405 
406  const label proci = Pstream::myProcNo();
407 
408  Info<< type() << " output:" << nl;
409 
410  Field<scalar> faceMassTotal(mass_.size(), 0.0);
411  this->getModelProperty("massTotal", faceMassTotal);
412 
413  Field<scalar> faceMassFlowRate(massFlowRate_.size(), 0.0);
414  this->getModelProperty("massFlowRate", faceMassFlowRate);
415 
416 
417  scalar sumTotalMass = 0.0;
418  scalar sumAverageMFR = 0.0;
419  forAll(faces_, facei)
420  {
421  scalarList allProcMass(Pstream::nProcs());
422  allProcMass[proci] = massTotal_[facei];
423  Pstream::gatherList(allProcMass);
424  faceMassTotal[facei] += sum(allProcMass);
425 
426  scalarList allProcMassFlowRate(Pstream::nProcs());
427  allProcMassFlowRate[proci] = massFlowRate_[facei];
428  Pstream::gatherList(allProcMassFlowRate);
429  faceMassFlowRate[facei] += sum(allProcMassFlowRate);
430 
431  sumTotalMass += faceMassTotal[facei];
432  sumAverageMFR += faceMassFlowRate[facei];
433 
434  if (outputFilePtr_.valid())
435  {
436  outputFilePtr_()
437  << time.timeName()
438  << tab << facei
439  << tab << faceMassTotal[facei]
440  << tab << faceMassFlowRate[facei]
441  << endl;
442  }
443  }
444 
445  Info<< " sum(total mass) = " << sumTotalMass << nl
446  << " sum(average mass flow rate) = " << sumAverageMFR << nl
447  << endl;
448 
449 
450  if (surfaceFormat_ != "none")
451  {
452  if (Pstream::master())
453  {
455  (
456  surfaceWriter::New(surfaceFormat_, time.writeFormat())
457  );
458 
459  writer->write
460  (
461  this->writeTimeDir(),
462  "collector",
463  points_,
464  faces_,
465  "massTotal",
466  faceMassTotal,
467  false
468  );
469 
470  writer->write
471  (
472  this->writeTimeDir(),
473  "collector",
474  points_,
475  faces_,
476  "massFlowRate",
477  faceMassFlowRate,
478  false
479  );
480  }
481  }
482 
483 
484  if (resetOnWrite_)
485  {
486  Field<scalar> dummy(faceMassTotal.size(), 0.0);
487  this->setModelProperty("massTotal", dummy);
488  this->setModelProperty("massFlowRate", dummy);
489 
490  timeOld_ = timeNew;
491  totalTime_ = 0.0;
492  }
493  else
494  {
495  this->setModelProperty("massTotal", faceMassTotal);
496  this->setModelProperty("massFlowRate", faceMassFlowRate);
497  }
498 
499  forAll(faces_, facei)
500  {
501  mass_[facei] = 0.0;
502  massTotal_[facei] = 0.0;
503  massFlowRate_[facei] = 0.0;
504  }
505 }
506 
507 
508 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
509 
510 template<class CloudType>
512 (
513  const dictionary& dict,
514  CloudType& owner,
515  const word& modelName
516 )
517 :
518  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
519  mode_(mtUnknown),
520  parcelType_(this->coeffDict().lookupOrDefault("parcelType", -1)),
521  removeCollected_(this->coeffDict().lookup("removeCollected")),
522  points_(),
523  faces_(),
524  nSector_(0),
525  radius_(),
526  coordSys_("coordSys", vector::zero, axesRotation(sphericalTensor::I)),
527  normal_(),
528  negateParcelsOppositeNormal_
529  (
530  readBool(this->coeffDict().lookup("negateParcelsOppositeNormal"))
531  ),
532  surfaceFormat_(this->coeffDict().lookup("surfaceFormat")),
533  resetOnWrite_(this->coeffDict().lookup("resetOnWrite")),
534  totalTime_(0.0),
535  mass_(),
536  massTotal_(),
537  massFlowRate_(),
538  log_(this->coeffDict().lookup("log")),
539  outputFilePtr_(),
540  timeOld_(owner.mesh().time().value()),
541  hitFaceIDs_()
542 {
543  normal_ /= mag(normal_);
544 
545  word mode(this->coeffDict().lookup("mode"));
546  if (mode == "polygon")
547  {
548  List<Field<point>> polygons(this->coeffDict().lookup("polygons"));
549 
550  initPolygons(polygons);
551 
552  vector n0(this->coeffDict().lookup("normal"));
553  normal_ = vectorField(faces_.size(), n0);
554  }
555  else if (mode == "polygonWithNormal")
556  {
557  List<Tuple2<Field<point>, vector>> polygonAndNormal
558  (
559  this->coeffDict().lookup("polygons")
560  );
561 
562  List<Field<point>> polygons(polygonAndNormal.size());
563  normal_.setSize(polygonAndNormal.size());
564 
565  forAll(polygons, polyI)
566  {
567  polygons[polyI] = polygonAndNormal[polyI].first();
568  normal_[polyI] = polygonAndNormal[polyI].second();
569  normal_[polyI] /= mag(normal_[polyI]) + rootVSmall;
570  }
571 
572  initPolygons(polygons);
573  }
574  else if (mode == "concentricCircle")
575  {
576  vector n0(this->coeffDict().lookup("normal"));
577  normal_ = vectorField(1, n0);
578 
579  initConcentricCircles();
580  }
581  else
582  {
583  FatalIOErrorInFunction(this->coeffDict())
584  << "Unknown mode " << mode << ". Available options are "
585  << "polygon, polygonWithNormal and concentricCircle"
586  << exit(FatalIOError);
587  }
588 
589  mass_.setSize(faces_.size(), 0.0);
590  massTotal_.setSize(faces_.size(), 0.0);
591  massFlowRate_.setSize(faces_.size(), 0.0);
592 
593  makeLogFile(faces_, points_, area_);
594 }
595 
596 
597 template<class CloudType>
599 (
601 )
602 :
604  mode_(pc.mode_),
605  parcelType_(pc.parcelType_),
606  removeCollected_(pc.removeCollected_),
607  points_(pc.points_),
608  faces_(pc.faces_),
609  nSector_(pc.nSector_),
610  radius_(pc.radius_),
611  coordSys_(pc.coordSys_),
612  normal_(pc.normal_),
613  negateParcelsOppositeNormal_(pc.negateParcelsOppositeNormal_),
614  surfaceFormat_(pc.surfaceFormat_),
615  resetOnWrite_(pc.resetOnWrite_),
616  totalTime_(pc.totalTime_),
617  mass_(pc.mass_),
618  massTotal_(pc.massTotal_),
619  massFlowRate_(pc.massFlowRate_),
620  log_(pc.log_),
621  outputFilePtr_(),
622  timeOld_(0.0),
623  hitFaceIDs_()
624 {}
625 
626 
627 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
628 
629 template<class CloudType>
631 {}
632 
633 
634 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
635 
636 template<class CloudType>
638 (
639  parcelType& p,
640  const scalar dt,
641  const point& position0,
642  bool& keepParticle
643 )
644 {
645  if ((parcelType_ != -1) && (parcelType_ != p.typeId()))
646  {
647  return;
648  }
649 
650  hitFaceIDs_.clear();
651 
652  switch (mode_)
653  {
654  case mtPolygon:
655  {
656  collectParcelPolygon(position0, p.position());
657  break;
658  }
659  case mtConcentricCircle:
660  {
661  collectParcelConcentricCircles(position0, p.position());
662  break;
663  }
664  default:
665  {}
666  }
667 
668  forAll(hitFaceIDs_, i)
669  {
670  label facei = hitFaceIDs_[i];
671  scalar m = p.nParticle()*p.mass();
672 
673  if (negateParcelsOppositeNormal_)
674  {
675  scalar Unormal = 0;
676  vector Uhat = p.U();
677  switch (mode_)
678  {
679  case mtPolygon:
680  {
681  Unormal = Uhat & normal_[facei];
682  break;
683  }
684  case mtConcentricCircle:
685  {
686  Unormal = Uhat & normal_[0];
687  break;
688  }
689  default:
690  {}
691  }
692 
693  Uhat /= mag(Uhat) + rootVSmall;
694 
695  if (Unormal < 0)
696  {
697  m = -m;
698  }
699  }
700 
701  // Add mass contribution
702  mass_[facei] += m;
703 
704  if (nSector_ == 1)
705  {
706  mass_[facei + 1] += m;
707  mass_[facei + 2] += m;
708  mass_[facei + 3] += m;
709  }
710 
711  if (removeCollected_)
712  {
713  keepParticle = false;
714  }
715  }
716 }
717 
718 
719 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:461
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:259
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
ParticleCollector(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Unit conversion functions.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
volVectorField vectorField(fieldObject, mesh)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
bool readBool(Istream &)
Definition: boolIO.C:60
void write()
Write post-processing info.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
T & first()
Return the first element of the list.
Definition: UListI.H:114
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
static const Identity< scalar > I
Definition: Identity.H:93
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
const Type & value() const
Return const reference to value.
virtual void postMove(parcelType &p, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
const scalar twoPi(2 *pi)
static const char nl
Definition: Ostream.H:260
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:285
virtual ~ParticleCollector()
Destructor.
labelList f(nPoints)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Vector< Cmpt > perpendicular(const Vector< Cmpt > &v)
Definition: VectorI.H:166
void setSize(const label)
Reset size of List.
Definition: List.C:281
vector point
Point is a vector.
Definition: point.H:41
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
Function object to collect the parcel mass- and mass flow rate over a set of polygons. The polygons can either be specified by sets of user- supplied points, or in a concentric circles arrangement. If a parcel is &#39;collected&#39;, it can be flagged to be removed from the domain using the removeCollected entry.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A coordinate rotation specified using global axis.
Definition: axesRotation.H:64
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Templated cloud function object base class.
IOerror FatalIOError