ParticleCollector.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 | Copyright (C) 2012-2016 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 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 template<class CloudType>
38 (
39  const faceList& faces,
40  const Field<point>& points,
41  const Field<scalar>& area
42 )
43 {
44  // Create the output file if not already created
45  if (log_)
46  {
47  if (debug)
48  {
49  Info<< "Creating output file" << endl;
50  }
51 
52  if (Pstream::master())
53  {
54  // Create directory if does not exist
55  mkDir(this->writeTimeDir());
56 
57  // Open new file at start up
58  outputFilePtr_.reset
59  (
60  new OFstream(this->writeTimeDir()/(type() + ".dat"))
61  );
62 
63  outputFilePtr_()
64  << "# Source : " << type() << nl
65  << "# Bins : " << faces.size() << nl
66  << "# Total area : " << sum(area) << nl;
67 
68  outputFilePtr_()
69  << "# Geometry :" << nl
70  << '#'
71  << tab << "Bin"
72  << tab << "(Centre_x Centre_y Centre_z)"
73  << tab << "Area"
74  << nl;
75 
76  forAll(faces, i)
77  {
78  outputFilePtr_()
79  << '#'
80  << tab << i
81  << tab << faces[i].centre(points)
82  << tab << area[i]
83  << nl;
84  }
85 
86  outputFilePtr_()
87  << '#' << nl
88  << "# Output format:" << nl;
89 
90  forAll(faces, i)
91  {
92  word id = Foam::name(i);
93  word binId = "bin_" + id;
94 
95  outputFilePtr_()
96  << '#'
97  << tab << "Time"
98  << tab << binId
99  << tab << "mass[" << id << "]"
100  << tab << "massFlowRate[" << id << "]"
101  << endl;
102  }
103  }
104  }
105 }
106 
107 
108 template<class CloudType>
110 (
111  const List<Field<point>>& polygons
112 )
113 {
114  mode_ = mtPolygon;
115 
116  label nPoints = 0;
117  forAll(polygons, polyI)
118  {
119  label np = polygons[polyI].size();
120  if (np < 3)
121  {
122  FatalIOErrorInFunction(this->coeffDict())
123  << "polygons must consist of at least 3 points"
124  << exit(FatalIOError);
125  }
126 
127  nPoints += np;
128  }
129 
130  label pointOffset = 0;
131  points_.setSize(nPoints);
132  faces_.setSize(polygons.size());
133  faceTris_.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 
142  DynamicList<face> tris;
143  f.triangles(points_, tris);
144  faceTris_[facei].transfer(tris);
145 
146  faces_[facei].transfer(f);
147 
148  pointOffset += polyPoints.size();
149  }
150 }
151 
152 
153 template<class CloudType>
155 {
156  mode_ = mtConcentricCircle;
157 
158  vector origin(this->coeffDict().lookup("origin"));
159 
160  radius_ = this->coeffDict().lookup("radius");
161  nSector_ = readLabel(this->coeffDict().lookup("nSector"));
162 
163  label nS = nSector_;
164 
165  vector refDir;
166  if (nSector_ > 1)
167  {
168  refDir = this->coeffDict().lookup("refDir");
169  refDir -= normal_[0]*(normal_[0] & refDir);
170  refDir /= mag(refDir);
171  }
172  else
173  {
174  // set 4 quadrants for single sector cases
175  nS = 4;
176 
177  vector tangent = Zero;
178  scalar magTangent = 0.0;
179 
180  Random rnd(1234);
181  while (magTangent < SMALL)
182  {
183  vector v = rnd.vector01();
184 
185  tangent = v - (v & normal_[0])*normal_[0];
186  magTangent = mag(tangent);
187  }
188 
189  refDir = tangent/magTangent;
190  }
191 
192  scalar dTheta = 5.0;
193  scalar dThetaSector = 360.0/scalar(nS);
194  label intervalPerSector = max(1, ceil(dThetaSector/dTheta));
195  dTheta = dThetaSector/scalar(intervalPerSector);
196 
197  label nPointPerSector = intervalPerSector + 1;
198 
199  label nPointPerRadius = nS*(nPointPerSector - 1);
200  label nPoint = radius_.size()*nPointPerRadius;
201  label nFace = radius_.size()*nS;
202 
203  // add origin
204  nPoint++;
205 
206  points_.setSize(nPoint);
207  faces_.setSize(nFace);
208  area_.setSize(nFace);
209 
210  coordSys_ = cylindricalCS("coordSys", origin, normal_[0], refDir, false);
211 
212  List<label> ptIDs(identity(nPointPerRadius));
213 
214  points_[0] = origin;
215 
216  // points
217  forAll(radius_, radI)
218  {
219  label pointOffset = radI*nPointPerRadius + 1;
220 
221  for (label i = 0; i < nPointPerRadius; i++)
222  {
223  label pI = i + pointOffset;
224  point pCyl(radius_[radI], degToRad(i*dTheta), 0.0);
225  points_[pI] = coordSys_.globalPosition(pCyl);
226  }
227  }
228 
229  // faces
230  DynamicList<label> facePts(2*nPointPerSector);
231  forAll(radius_, radI)
232  {
233  if (radI == 0)
234  {
235  for (label secI = 0; secI < nS; secI++)
236  {
237  facePts.clear();
238 
239  // append origin point
240  facePts.append(0);
241 
242  for (label ptI = 0; ptI < nPointPerSector; ptI++)
243  {
244  label i = ptI + secI*(nPointPerSector - 1);
245  label id = ptIDs.fcIndex(i - 1) + 1;
246  facePts.append(id);
247  }
248 
249  label facei = secI + radI*nS;
250 
251  faces_[facei] = face(facePts);
252  area_[facei] = faces_[facei].mag(points_);
253  }
254  }
255  else
256  {
257  for (label secI = 0; secI < nS; secI++)
258  {
259  facePts.clear();
260 
261  label offset = (radI - 1)*nPointPerRadius + 1;
262 
263  for (label ptI = 0; ptI < nPointPerSector; ptI++)
264  {
265  label i = ptI + secI*(nPointPerSector - 1);
266  label id = offset + ptIDs.fcIndex(i - 1);
267  facePts.append(id);
268  }
269  for (label ptI = nPointPerSector-1; ptI >= 0; ptI--)
270  {
271  label i = ptI + secI*(nPointPerSector - 1);
272  label id = offset + nPointPerRadius + ptIDs.fcIndex(i - 1);
273  facePts.append(id);
274  }
275 
276  label facei = secI + radI*nS;
277 
278  faces_[facei] = face(facePts);
279  area_[facei] = faces_[facei].mag(points_);
280  }
281  }
282  }
283 }
284 
285 
286 template<class CloudType>
288 (
289  const point& p1,
290  const point& p2
291 ) const
292 {
293  label dummyNearType = -1;
294  label dummyNearLabel = -1;
295 
296  forAll(faces_, facei)
297  {
298  const label facePoint0 = faces_[facei][0];
299 
300  const point& pf = points_[facePoint0];
301 
302  const scalar d1 = normal_[facei] & (p1 - pf);
303  const scalar d2 = normal_[facei] & (p2 - pf);
304 
305  if (sign(d1) == sign(d2))
306  {
307  // did not cross polygon plane
308  continue;
309  }
310 
311  // intersection point
312  const point pIntersect = p1 + (d1/(d1 - d2))*(p2 - p1);
313 
314  const List<face>& tris = faceTris_[facei];
315 
316  // identify if point is within poly bounds
317  forAll(tris, triI)
318  {
319  const face& tri = tris[triI];
320  triPointRef t
321  (
322  points_[tri[0]],
323  points_[tri[1]],
324  points_[tri[2]]
325  );
326 
327  if (t.classify(pIntersect, dummyNearType, dummyNearLabel))
328  {
329  hitFaceIDs_.append(facei);
330  }
331  }
332  }
333 }
334 
335 
336 template<class CloudType>
338 (
339  const point& p1,
340  const point& p2
341 ) const
342 {
343  label secI = -1;
344 
345  const scalar d1 = normal_[0] & (p1 - coordSys_.origin());
346  const scalar d2 = normal_[0] & (p2 - coordSys_.origin());
347 
348  if (sign(d1) == sign(d2))
349  {
350  // did not cross plane
351  return;
352  }
353 
354  // intersection point in cylindrical co-ordinate system
355  const point pCyl = coordSys_.localPosition(p1 + (d1/(d1 - d2))*(p2 - p1));
356 
357  scalar r = pCyl[0];
358 
359  if (r < radius_.last())
360  {
361  label radI = 0;
362  while (r > radius_[radI])
363  {
364  radI++;
365  }
366 
367  if (nSector_ == 1)
368  {
369  secI = 4*radI;
370  }
371  else
372  {
373  scalar theta = pCyl[1] + constant::mathematical::pi;
374 
375  secI =
376  nSector_*radI
377  + floor
378  (
379  scalar(nSector_)*theta/constant::mathematical::twoPi
380  );
381  }
382  }
383 
384  hitFaceIDs_.append(secI);
385 }
386 
387 
388 template<class CloudType>
390 {
391  const fvMesh& mesh = this->owner().mesh();
392  const Time& time = mesh.time();
393  scalar timeNew = time.value();
394  scalar timeElapsed = timeNew - timeOld_;
395 
396  totalTime_ += timeElapsed;
397 
398  const scalar alpha = (totalTime_ - timeElapsed)/totalTime_;
399  const scalar beta = timeElapsed/totalTime_;
400 
401  forAll(faces_, facei)
402  {
403  massFlowRate_[facei] =
404  alpha*massFlowRate_[facei] + beta*mass_[facei]/timeElapsed;
405  massTotal_[facei] += mass_[facei];
406  }
407 
408  const label proci = Pstream::myProcNo();
409 
410  Info<< type() << " output:" << nl;
411 
412  Field<scalar> faceMassTotal(mass_.size(), 0.0);
413  this->getModelProperty("massTotal", faceMassTotal);
414 
415  Field<scalar> faceMassFlowRate(massFlowRate_.size(), 0.0);
416  this->getModelProperty("massFlowRate", faceMassFlowRate);
417 
418 
419  scalar sumTotalMass = 0.0;
420  scalar sumAverageMFR = 0.0;
421  forAll(faces_, facei)
422  {
423  scalarList allProcMass(Pstream::nProcs());
424  allProcMass[proci] = massTotal_[facei];
425  Pstream::gatherList(allProcMass);
426  faceMassTotal[facei] += sum(allProcMass);
427 
428  scalarList allProcMassFlowRate(Pstream::nProcs());
429  allProcMassFlowRate[proci] = massFlowRate_[facei];
430  Pstream::gatherList(allProcMassFlowRate);
431  faceMassFlowRate[facei] += sum(allProcMassFlowRate);
432 
433  sumTotalMass += faceMassTotal[facei];
434  sumAverageMFR += faceMassFlowRate[facei];
435 
436  if (outputFilePtr_.valid())
437  {
438  outputFilePtr_()
439  << time.timeName()
440  << tab << facei
441  << tab << faceMassTotal[facei]
442  << tab << faceMassFlowRate[facei]
443  << endl;
444  }
445  }
446 
447  Info<< " sum(total mass) = " << sumTotalMass << nl
448  << " sum(average mass flow rate) = " << sumAverageMFR << nl
449  << endl;
450 
451 
452  if (surfaceFormat_ != "none")
453  {
454  if (Pstream::master())
455  {
457 
458  writer->write
459  (
460  this->writeTimeDir(),
461  "collector",
462  points_,
463  faces_,
464  "massTotal",
465  faceMassTotal,
466  false
467  );
468 
469  writer->write
470  (
471  this->writeTimeDir(),
472  "collector",
473  points_,
474  faces_,
475  "massFlowRate",
476  faceMassFlowRate,
477  false
478  );
479  }
480  }
481 
482 
483  if (resetOnWrite_)
484  {
485  Field<scalar> dummy(faceMassTotal.size(), 0.0);
486  this->setModelProperty("massTotal", dummy);
487  this->setModelProperty("massFlowRate", dummy);
488 
489  timeOld_ = timeNew;
490  totalTime_ = 0.0;
491  }
492  else
493  {
494  this->setModelProperty("massTotal", faceMassTotal);
495  this->setModelProperty("massFlowRate", faceMassFlowRate);
496  }
497 
498  forAll(faces_, facei)
499  {
500  mass_[facei] = 0.0;
501  massTotal_[facei] = 0.0;
502  massFlowRate_[facei] = 0.0;
503  }
504 }
505 
506 
507 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
508 
509 template<class CloudType>
511 (
512  const dictionary& dict,
513  CloudType& owner,
514  const word& modelName
515 )
516 :
517  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
518  mode_(mtUnknown),
519  parcelType_(this->coeffDict().lookupOrDefault("parcelType", -1)),
520  removeCollected_(this->coeffDict().lookup("removeCollected")),
521  points_(),
522  faces_(),
523  faceTris_(),
524  nSector_(0),
525  radius_(),
526  coordSys_(false),
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  faceTris_(pc.faceTris_),
610  nSector_(pc.nSector_),
611  radius_(pc.radius_),
612  coordSys_(pc.coordSys_),
613  normal_(pc.normal_),
614  negateParcelsOppositeNormal_(pc.negateParcelsOppositeNormal_),
615  surfaceFormat_(pc.surfaceFormat_),
616  resetOnWrite_(pc.resetOnWrite_),
617  totalTime_(pc.totalTime_),
618  mass_(pc.mass_),
619  massTotal_(pc.massTotal_),
620  massFlowRate_(pc.massFlowRate_),
621  log_(pc.log_),
622  outputFilePtr_(),
623  timeOld_(0.0),
624  hitFaceIDs_()
625 {}
626 
627 
628 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
629 
630 template<class CloudType>
632 {}
633 
634 
635 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
636 
637 template<class CloudType>
639 (
640  parcelType& p,
641  const label celli,
642  const scalar dt,
643  const point& position0,
644  bool& keepParticle
645 )
646 {
647  if ((parcelType_ != -1) && (parcelType_ != p.typeId()))
648  {
649  return;
650  }
651 
652  // slightly extend end position to avoid falling within tracking tolerances
653  const point position1 = position0 + 1.0001*(p.position() - position0);
654 
655  hitFaceIDs_.clear();
656 
657  switch (mode_)
658  {
659  case mtPolygon:
660  {
661  collectParcelPolygon(position0, position1);
662  break;
663  }
664  case mtConcentricCircle:
665  {
666  collectParcelConcentricCircles(position0, position1);
667  break;
668  }
669  default:
670  {
671  }
672  }
673 
674 
675  forAll(hitFaceIDs_, i)
676  {
677  label facei = hitFaceIDs_[i];
678  scalar m = p.nParticle()*p.mass();
679 
680  if (negateParcelsOppositeNormal_)
681  {
682  vector Uhat = p.U();
683  Uhat /= mag(Uhat) + ROOTVSMALL;
684  if ((Uhat & normal_[facei]) < 0)
685  {
686  m *= -1.0;
687  }
688  }
689 
690  // add mass contribution
691  mass_[facei] += m;
692 
693  if (nSector_ == 1)
694  {
695  mass_[facei + 1] += m;
696  mass_[facei + 2] += m;
697  mass_[facei + 3] += m;
698  }
699 
700  if (removeCollected_)
701  {
702  keepParticle = false;
703  }
704  }
705 }
706 
707 
708 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:261
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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:76
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
volVectorField vectorField(fieldObject, mesh)
Base class for graphics format writing. Entry points are.
Definition: writer.H:78
const Type & value() const
Return const reference to value.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
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.
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)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
static const zero Zero
Definition: zero.H:91
label readLabel(Istream &is)
Definition: label.H:64
const scalar twoPi(2 *pi)
static const char nl
Definition: Ostream.H:262
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:295
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setSize(const label)
Reset size of List.
Definition: List.C:295
virtual void postMove(parcelType &p, const label celli, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
vector point
Point is a vector.
Definition: point.H:41
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
mode_t mode(const fileName &)
Return the file mode.
Definition: POSIX.C:447
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
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
messageStream Info
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
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:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Templated cloud function object base class.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
IOerror FatalIOError