streamLineBase.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 OpenFOAM Foundation
9  Copyright (C) 2015-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "streamLineBase.H"
30 #include "fvMesh.H"
31 #include "ReadFields.H"
32 #include "sampledSet.H"
33 #include "globalIndex.H"
34 #include "mapDistribute.H"
35 #include "interpolationCellPoint.H"
36 #include "wallPolyPatch.H"
37 #include "meshSearchMeshObject.H"
38 #include "mapPolyMesh.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace functionObjects
45 {
46  defineTypeNameAndDebug(streamLineBase, 0);
47 }
48 }
49 
50 
51 const Foam::Enum
52 <
54 >
56 ({
57  { trackDirType::FORWARD, "forward" },
58  { trackDirType::BACKWARD, "backward" },
59  { trackDirType::BIDIRECTIONAL, "bidirectional" }
60 });
61 
62 
63 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64 
65 const Foam::word&
67 {
68  if (!sampledSetPtr_)
69  {
71  }
72 
73  return sampledSetAxis_;
74 }
75 
76 
77 const Foam::sampledSet&
79 {
80  if (!sampledSetPtr_)
81  {
82  sampledSetPtr_ = sampledSet::New
83  (
84  "seedSampleSet",
85  mesh_,
87  dict_.subDict("seedSampleSet")
88  );
89 
90  sampledSetAxis_ = sampledSetPtr_->axis();
91  }
92 
93  return *sampledSetPtr_;
94 }
95 
96 
99 {
100  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
101 
102  label nFaces = 0;
103 
104  for (const polyPatch& pp : patches)
105  {
106  //if (!polyPatch::constraintType(pp.type()))
107  if (isA<wallPolyPatch>(pp))
108  {
109  nFaces += pp.size();
110  }
111  }
112 
113  labelList addressing(nFaces);
114 
115  nFaces = 0;
116 
117  for (const polyPatch& pp : patches)
118  {
119  //if (!polyPatch::constraintType(pp.type()))
120  if (isA<wallPolyPatch>(pp))
121  {
122  forAll(pp, i)
123  {
124  addressing[nFaces++] = pp.start()+i;
125  }
126  }
127  }
128 
130  (
132  (
133  mesh_.faces(),
134  addressing
135  ),
136  mesh_.points()
137  );
138 }
139 
140 
142 (
143  const label nSeeds,
144  label& UIndex,
145  PtrList<volScalarField>& vsFlds,
146  PtrList<interpolation<scalar>>& vsInterp,
147  PtrList<volVectorField>& vvFlds,
148  PtrList<interpolation<vector>>& vvInterp
149 )
150 {
151  // Read fields
152  label nScalar = 0;
153  label nVector = 0;
154 
155  for (const word& fieldName : fields_)
156  {
157  if (foundObject<volScalarField>(fieldName))
158  {
159  ++nScalar;
160  }
161  else if (foundObject<volVectorField>(fieldName))
162  {
163  ++nVector;
164  }
165  else
166  {
168  << "Cannot find field " << fieldName << nl
169  << "Valid scalar fields are:"
170  << flatOutput(mesh_.names(volScalarField::typeName)) << nl
171  << "Valid vector fields are:"
172  << flatOutput(mesh_.names(volVectorField::typeName))
173  << exit(FatalError);
174  }
175  }
176  vsInterp.setSize(nScalar);
177  nScalar = 0;
178  vvInterp.setSize(nVector);
179  nVector = 0;
180 
181  for (const word& fieldName : fields_)
182  {
183  if (foundObject<volScalarField>(fieldName))
184  {
185  const volScalarField& f = lookupObject<volScalarField>(fieldName);
186  vsInterp.set
187  (
188  nScalar++,
190  (
191  interpolationScheme_,
192  f
193  )
194  );
195  }
196  else if (foundObject<volVectorField>(fieldName))
197  {
198  const volVectorField& f = lookupObject<volVectorField>(fieldName);
199 
200  if (f.name() == UName_)
201  {
202  UIndex = nVector;
203  }
204 
205  vvInterp.set
206  (
207  nVector++,
209  (
210  interpolationScheme_,
211  f
212  )
213  );
214  }
215  }
216 
217  // Store the names
218  scalarNames_.setSize(vsInterp.size());
219  forAll(vsInterp, i)
220  {
221  scalarNames_[i] = vsInterp[i].psi().name();
222  }
223  vectorNames_.setSize(vvInterp.size());
224  forAll(vvInterp, i)
225  {
226  vectorNames_[i] = vvInterp[i].psi().name();
227  }
228 
229  // Check that we know the index of U in the interpolators.
230 
231  if (UIndex == -1)
232  {
234  << "Cannot find field to move particles with : " << UName_ << nl
235  << "This field has to be present in the sampled fields " << fields_
236  << " and in the objectRegistry."
237  << exit(FatalError);
238  }
239 
240  // Sampled data
241  // ~~~~~~~~~~~~
242 
243  // Size to maximum expected sizes.
244  allTracks_.clear();
245  allTracks_.setCapacity(nSeeds);
246  allScalars_.setSize(vsInterp.size());
247  forAll(allScalars_, i)
248  {
249  allScalars_[i].clear();
250  allScalars_[i].setCapacity(nSeeds);
251  }
252  allVectors_.setSize(vvInterp.size());
253  forAll(allVectors_, i)
254  {
255  allVectors_[i].clear();
256  allVectors_[i].setCapacity(nSeeds);
257  }
258 }
259 
260 
262 (
263  const label tracki,
264 
265  const scalar w,
266  const label lefti,
267  const label righti,
268 
269  DynamicList<point>& newTrack,
270  DynamicList<scalarList>& newScalars,
271  DynamicList<vectorList>& newVectors
272 ) const
273 {
274  const label sz = newTrack.size();
275 
276  const List<point>& track = allTracks_[tracki];
277 
278  newTrack.append((1.0-w)*track[lefti] + w*track[righti]);
279 
280  // Scalars
281  {
282  newScalars.append(scalarList(allScalars_.size()));
283  scalarList& newVals = newScalars[sz];
284 
285  forAll(allScalars_, scalari)
286  {
287  const scalarList& trackVals = allScalars_[scalari][tracki];
288  newVals[scalari] = (1.0-w)*trackVals[lefti] + w*trackVals[righti];
289  }
290  }
291 
292  // Vectors
293  {
294  newVectors.append(vectorList(allVectors_.size()));
295  vectorList& newVals = newVectors[sz];
296 
297  forAll(allVectors_, vectori)
298  {
299  const vectorList& trackVals = allVectors_[vectori][tracki];
300  newVals[vectori] = (1.0-w)*trackVals[lefti] + w*trackVals[righti];
301  }
302  }
303 }
304 
305 
306 // Can split a track into multiple tracks
308 (
309  const treeBoundBox& bb,
310  const label tracki,
311  PtrList<DynamicList<point>>& newTracks,
312  PtrList<DynamicList<scalarList>>& newScalars,
313  PtrList<DynamicList<vectorList>>& newVectors
314 ) const
315 {
316  const List<point>& track = allTracks_[tracki];
317 
318  if (track.size())
319  {
320  for
321  (
322  label segmenti = 1;
323  segmenti < track.size();
324  segmenti++
325  )
326  {
327  const point& startPt = track[segmenti-1];
328  const point& endPt = track[segmenti];
329 
330  const vector d(endPt-startPt);
331  const scalar magD = mag(d);
332  if (magD > ROOTVSMALL)
333  {
334  if (bb.contains(startPt))
335  {
336  // Store 1.0*track[segmenti-1]+0*track[segmenti]
337  storePoint
338  (
339  tracki,
340 
341  0.0,
342  segmenti-1,
343  segmenti,
344 
345  newTracks.last(),
346  newScalars.last(),
347  newVectors.last()
348  );
349 
350  if (!bb.contains(endPt))
351  {
352  point clipPt;
353  if (bb.intersects(endPt, startPt, clipPt))
354  {
355  // End of track. Store point and interpolated
356  // values
357  storePoint
358  (
359  tracki,
360 
361  mag(clipPt-startPt)/magD,
362  segmenti-1,
363  segmenti,
364 
365  newTracks.last(),
366  newScalars.last(),
367  newVectors.last()
368  );
369 
370  newTracks.last().shrink();
371  newScalars.last().shrink();
372  newVectors.last().shrink();
373  }
374  }
375  }
376  else
377  {
378  // startPt outside box. New track. Get starting point
379 
380  point clipPt;
381  if (bb.intersects(startPt, endPt, clipPt))
382  {
383  // New track
384  newTracks.append
385  (
386  new DynamicList<point>(track.size()/10)
387  );
388  newScalars.append
389  (
390  new DynamicList<scalarList>(track.size()/10)
391  );
392  newVectors.append
393  (
394  new DynamicList<vectorList>(track.size()/10)
395  );
396 
397  // Store point and interpolated values
398  storePoint
399  (
400  tracki,
401 
402  mag(clipPt-startPt)/magD,
403  segmenti-1,
404  segmenti,
405 
406  newTracks.last(),
407  newScalars.last(),
408  newVectors.last()
409  );
410 
411  if (!bb.contains(endPt))
412  {
413  bb.intersects
414  (
415  endPt,
416  point(clipPt),
417  clipPt
418  );
419 
420  // Store point and interpolated values
421  storePoint
422  (
423  tracki,
424 
425  mag(clipPt-startPt)/magD,
426  segmenti-1,
427  segmenti,
428 
429  newTracks.last(),
430  newScalars.last(),
431  newVectors.last()
432  );
433 
434  newTracks.last().shrink();
435  newScalars.last().shrink();
436  newVectors.last().shrink();
437  }
438  }
439  }
440  }
441  }
442 
443  // Last point
444  if (bb.contains(track.last()))
445  {
446  storePoint
447  (
448  tracki,
449 
450  1.0,
451  track.size()-2,
452  track.size()-1,
453 
454  newTracks.last(),
455  newScalars.last(),
456  newVectors.last()
457  );
458  }
459  }
460 }
461 
462 
464 {
465  // Storage for new tracks. Per track, per sample the coordinate (newTracks)
466  // or values for all the sampled fields (newScalars, newVectors)
467  PtrList<DynamicList<point>> newTracks;
470 
471  forAll(allTracks_, tracki)
472  {
473  const List<point>& track = allTracks_[tracki];
474 
475  if (track.size())
476  {
477  // New track. Assume it consists of the whole track
478  newTracks.append(new DynamicList<point>(track.size()));
479  newScalars.append(new DynamicList<scalarList>(track.size()));
480  newVectors.append(new DynamicList<vectorList>(track.size()));
481 
482  // Trim, split and append to newTracks
483  trimToBox(bb, tracki, newTracks, newScalars, newVectors);
484  }
485  }
486 
487  // Transfer newTracks to allTracks_
488  allTracks_.setSize(newTracks.size());
489  forAll(allTracks_, tracki)
490  {
491  allTracks_[tracki].transfer(newTracks[tracki]);
492  }
493  // Replace track scalars
494  forAll(allScalars_, scalari)
495  {
496  DynamicList<scalarList>& fieldVals = allScalars_[scalari];
497  fieldVals.setSize(newTracks.size());
498 
499  forAll(fieldVals, tracki)
500  {
501  scalarList& trackVals = allScalars_[scalari][tracki];
502  trackVals.setSize(newScalars[tracki].size());
503  forAll(trackVals, samplei)
504  {
505  trackVals[samplei] = newScalars[tracki][samplei][scalari];
506  }
507  }
508  }
509  // Replace track vectors
510  forAll(allVectors_, vectori)
511  {
512  DynamicList<vectorList>& fieldVals = allVectors_[vectori];
513  fieldVals.setSize(newTracks.size());
514  forAll(fieldVals, tracki)
515  {
516  vectorList& trackVals = allVectors_[vectori][tracki];
517  trackVals.setSize(newVectors[tracki].size());
518  forAll(trackVals, samplei)
519  {
520  trackVals[samplei] = newVectors[tracki][samplei][vectori];
521  }
522  }
523  }
524 }
525 
526 
528 {
529  if (Pstream::parRun())
530  {
531  // Append slave tracks to master ones
532  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
533 
534  globalIndex globalTrackIDs(allTracks_.size());
535 
536  // Construct a distribution map to pull all to the master.
537  labelListList sendMap(Pstream::nProcs());
538  labelListList recvMap(Pstream::nProcs());
539 
540  if (Pstream::master())
541  {
542  // Master: receive all. My own first, then consecutive
543  // processors.
544  label tracki = 0;
545 
546  forAll(recvMap, proci)
547  {
548  labelList& fromProc = recvMap[proci];
549  fromProc.setSize(globalTrackIDs.localSize(proci));
550  forAll(fromProc, i)
551  {
552  fromProc[i] = tracki++;
553  }
554  }
555  }
556 
557  labelList& toMaster = sendMap[0];
558  toMaster.setSize(globalTrackIDs.localSize());
559  forAll(toMaster, i)
560  {
561  toMaster[i] = i;
562  }
563 
564  const mapDistribute distMap
565  (
566  globalTrackIDs.size(),
567  std::move(sendMap),
568  std::move(recvMap)
569  );
570 
571 
572  // Distribute the track positions. Note: use scheduled comms
573  // to prevent buffering.
574  allTracks_.shrink();
576  (
578  distMap.schedule(),
579  distMap.constructSize(),
580  distMap.subMap(),
581  false,
582  distMap.constructMap(),
583  false,
584  allTracks_,
585  flipOp()
586  );
587  allTracks_.setCapacity(allTracks_.size());
588 
589  // Distribute the scalars
590  forAll(allScalars_, scalari)
591  {
592  allScalars_[scalari].shrink();
594  (
596  distMap.schedule(),
597  distMap.constructSize(),
598  distMap.subMap(),
599  false,
600  distMap.constructMap(),
601  false,
602  allScalars_[scalari],
603  flipOp()
604  );
605  allScalars_[scalari].setCapacity(allScalars_[scalari].size());
606  }
607  // Distribute the vectors
608  forAll(allVectors_, vectori)
609  {
610  allVectors_[vectori].shrink();
612  (
614  distMap.schedule(),
615  distMap.constructSize(),
616  distMap.subMap(),
617  false,
618  distMap.constructMap(),
619  false,
620  allVectors_[vectori],
621  flipOp()
622  );
623  allVectors_[vectori].setCapacity(allVectors_[vectori].size());
624  }
625  }
626 
627 
628  // Note: filenames scattered below since used in global call
629  fileName scalarVtkFile;
630  fileName vectorVtkFile;
631 
632  if (Pstream::master())
633  {
634  if (!bounds_.empty())
635  {
636  // Clip to bounding box
637  trimToBox(treeBoundBox(bounds_));
638  }
639 
640 
641  label nTracks = 0;
642  label n = 0;
643  forAll(allTracks_, tracki)
644  {
645  if (allTracks_[tracki].size())
646  {
647  nTracks++;
648  n += allTracks_[tracki].size();
649  }
650  }
651 
652  Log << " Tracks:" << nTracks << nl
653  << " Total samples:" << n
654  << endl;
655 
656 
657  // Massage into form suitable for writers
658  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
659 
660  // Make output directory
661 
662  fileName vtkPath
663  (
664  time_.globalPath()/functionObject::outputPrefix/"sets"/name()
665  );
666  if (mesh_.name() != polyMesh::defaultRegion)
667  {
668  vtkPath = vtkPath/mesh_.name();
669  }
670  vtkPath = vtkPath/mesh_.time().timeName();
671 
672  mkDir(vtkPath);
673 
674  // Convert track positions (and compact out empty tracks)
675 
676  PtrList<coordSet> tracks(nTracks);
677  nTracks = 0;
678  labelList oldToNewTrack(allTracks_.size(), -1);
679 
680  forAll(allTracks_, tracki)
681  {
682  if (allTracks_[tracki].size())
683  {
684  List<point>& points = allTracks_[tracki];
685  scalarList dist(points.size());
686  dist[0] = 0;
687  for (label pointi = 1; pointi < points.size(); ++pointi)
688  {
689  dist[pointi] =
690  dist[pointi-1] + mag(points[pointi] - points[pointi-1]);
691  }
692 
693  tracks.set
694  (
695  nTracks,
696  new coordSet
697  (
698  "track" + Foam::name(nTracks),
699  sampledSetAxis(), // "xyz"
700  std::move(allTracks_[tracki]),
701  std::move(dist)
702  )
703  );
704  oldToNewTrack[tracki] = nTracks;
705  ++nTracks;
706  }
707  }
708 
709  // Convert scalar values
710 
711  if (!allScalars_.empty() && !tracks.empty())
712  {
713  List<List<scalarField>> scalarValues(allScalars_.size());
714 
715  forAll(allScalars_, scalari)
716  {
717  DynamicList<scalarList>& allTrackVals = allScalars_[scalari];
718  scalarValues[scalari].setSize(nTracks);
719 
720  forAll(allTrackVals, tracki)
721  {
722  scalarList& vals = allTrackVals[tracki];
723  if (vals.size())
724  {
725  const label newTracki = oldToNewTrack[tracki];
726  scalarValues[scalari][newTracki].transfer(vals);
727  }
728  }
729  }
730 
731  scalarVtkFile = fileName
732  (
733  vtkPath
734  / scalarFormatterPtr_().getFileName
735  (
736  tracks[0],
737  scalarNames_
738  )
739  );
740 
741  Log << " Writing data to " << scalarVtkFile.path() << endl;
742 
743  scalarFormatterPtr_().write
744  (
745  true, // writeTracks
746  List<scalarField>(), // times
747  tracks,
748  scalarNames_,
749  scalarValues,
750  OFstream(scalarVtkFile)()
751  );
752  }
753 
754  // Convert vector values
755 
756  if (!allVectors_.empty() && !tracks.empty())
757  {
758  List<List<vectorField>> vectorValues(allVectors_.size());
759 
760  forAll(allVectors_, vectori)
761  {
762  DynamicList<vectorList>& allTrackVals = allVectors_[vectori];
763  vectorValues[vectori].setSize(nTracks);
764 
765  forAll(allTrackVals, tracki)
766  {
767  vectorList& vals = allTrackVals[tracki];
768  if (vals.size())
769  {
770  const label newTracki = oldToNewTrack[tracki];
771  vectorValues[vectori][newTracki].transfer(vals);
772  }
773  }
774  }
775 
776  vectorVtkFile = fileName
777  (
778  vtkPath
779  / vectorFormatterPtr_().getFileName(tracks[0], vectorNames_)
780  );
781 
782  //Info<< " Writing vector data to " << vectorVtkFile << endl;
783 
784  vectorFormatterPtr_().write
785  (
786  true, // writeTracks
787  List<scalarField>(), // times
788  tracks,
789  vectorNames_,
790  vectorValues,
791  OFstream(vectorVtkFile)()
792  );
793  }
794  }
795 
796 
797  // File names are generated on the master but setProperty needs to
798  // be across all procs
799  Pstream::scatter(scalarVtkFile);
800  for (const word& fieldName : scalarNames_)
801  {
803  propsDict.add
804  (
805  "file",
806  time_.relativePath(scalarVtkFile, true)
807  );
808  setProperty(fieldName, propsDict);
809  }
810 
811  Pstream::scatter(vectorVtkFile);
812  for (const word& fieldName : vectorNames_)
813  {
815  propsDict.add
816  (
817  "file",
818  time_.relativePath(vectorVtkFile, true)
819  );
820  setProperty(fieldName, propsDict);
821  }
822 
823  return true;
824 }
825 
826 
828 (
829  const word& newUName,
830  const wordList& newFieldNames
831 )
832 {
833  UName_ = newUName;
834  fields_ = newFieldNames;
835 }
836 
837 
838 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
839 
841 (
842  const word& name,
843  const Time& runTime,
844  const dictionary& dict
845 )
846 :
848  dict_(dict),
849  fields_()
850 {}
851 
852 
854 (
855  const word& name,
856  const Time& runTime,
857  const dictionary& dict,
858  const wordList& fieldNames
859 )
860 :
862  dict_(dict),
863  fields_(fieldNames)
864 {}
865 
866 
867 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
868 
870 {}
871 
872 
873 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
874 
876 {
877  if (&dict_ != &dict)
878  {
879  // Update local copy of dictionary:
880  dict_ = dict;
881  }
882 
884 
885  Info<< type() << " " << name() << ":" << nl;
886 
887  UName_ = dict.getOrDefault<word>("U", "U");
888 
889  if (fields_.empty())
890  {
891  dict.readEntry("fields", fields_);
892 
893  if (!fields_.found(UName_))
894  {
896  << "Velocity field for tracking " << UName_
897  << " should be present in the list of fields " << fields_
898  << exit(FatalIOError);
899  }
900  }
901 
902  Info<< " Employing velocity field " << UName_ << endl;
903 
904  bool trackForward;
905  if (dict.readIfPresent("trackForward", trackForward))
906  {
907  trackDir_ =
908  (
909  trackForward
910  ? trackDirType::FORWARD
911  : trackDirType::BACKWARD
912  );
913 
914  if (dict.found("direction"))
915  {
917  << "Cannot specify both 'trackForward' and 'direction'" << nl
918  << exit(FatalIOError);
919  }
920  }
921  else
922  {
923  trackDir_ = trackDirTypeNames.get("direction", dict);
924  }
925  dict.readEntry("lifeTime", lifeTime_);
926  if (lifeTime_ < 1)
927  {
929  << "Illegal value " << lifeTime_ << " for lifeTime"
930  << exit(FatalIOError);
931  }
932 
933 
934  trackLength_ = VGREAT;
935  if (dict.readIfPresent("trackLength", trackLength_))
936  {
937  Info<< type() << " : fixed track length specified : "
938  << trackLength_ << nl << endl;
939  }
940 
941 
942  bounds_ = boundBox::invertedBox;
943  if (dict.readIfPresent("bounds", bounds_) && !bounds_.empty())
944  {
945  Info<< " clipping all segments to " << bounds_ << nl << endl;
946  }
947 
948 
949  interpolationScheme_ = dict.getOrDefault
950  (
951  "interpolationScheme",
953  );
954 
955  //Info<< " using interpolation " << interpolationScheme_ << endl;
956 
957  cloudName_ = dict.getOrDefault<word>("cloud", type());
958 
959  sampledSetPtr_.clear();
960  sampledSetAxis_.clear();
961 
962  scalarFormatterPtr_ = writer<scalar>::New(dict.get<word>("setFormat"));
963  vectorFormatterPtr_ = writer<vector>::New(dict.get<word>("setFormat"));
964 
965  return true;
966 }
967 
968 
970 {
971  return true;
972 }
973 
974 
976 {
977  Log << type() << " " << name() << " write:" << nl;
978 
979  // Do all injection and tracking
980  track();
981 
982  writeToFile();
983 
984  return true;
985 }
986 
987 
989 {
990  if (&mpm.mesh() == &mesh_)
991  {
992  read(dict_);
993  }
994 }
995 
996 
998 {
999  // Moving mesh affects the search tree
1000  read(dict_);
1001 }
1002 
1003 
1004 // ************************************************************************* //
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::sampledSet::New
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition: sampledSet.C:519
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::functionObjects::streamLineBase::trackDirType
trackDirType
Enumeration defining the track direction.
Definition: streamLineBase.H:74
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:318
Foam::fileName::path
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
Foam::functionObjects::streamLineBase::~streamLineBase
virtual ~streamLineBase()
Destructor.
Definition: streamLineBase.C:869
meshSearchMeshObject.H
Foam::DynamicList< point >
mapPolyMesh.H
globalIndex.H
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::functionObjects::streamLineBase::write
virtual bool write()
Track and write.
Definition: streamLineBase.C:975
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::functionObjects::streamLineBase::sampledSetPoints
const sampledSet & sampledSetPoints() const
Demand driven construction of the sampledSet.
Definition: streamLineBase.C:78
Foam::MeshObject< polyMesh, GeometricMeshObject, meshSearchMeshObject >::New
static const meshSearchMeshObject & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
wallPolyPatch.H
Foam::functionObjects::streamLineBase::storePoint
void storePoint(const label tracki, const scalar w, const label lefti, const label righti, DynamicList< point > &newTrack, DynamicList< List< scalar >> &newScalars, DynamicList< List< vector >> &newVectors) const
Generate point and values by interpolating from existing values.
Definition: streamLineBase.C:262
Foam::functionObjects::streamLineBase::wallPatch
autoPtr< indirectPrimitivePatch > wallPatch() const
Construct patch out of all wall patch faces.
Definition: streamLineBase.C:98
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (can be nullptr),.
Definition: PtrList.H:138
Foam::functionObjects::streamLineBase::movePoints
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: streamLineBase.C:997
Foam::DynamicList::shrink
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::functionObjects::streamLineBase::read
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamLineBase.C:875
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
interpolationCellPoint.H
Foam::functionObjects::streamLineBase::writeToFile
virtual bool writeToFile()
Write tracks to file.
Definition: streamLineBase.C:527
Foam::functionObjects::streamLineBase::sampledSetPtr_
autoPtr< sampledSet > sampledSetPtr_
Seed set engine.
Definition: streamLineBase.H:88
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::functionObjects::streamLineBase::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: streamLineBase.C:988
streamLineBase.H
sampledSet.H
Foam::interpolationCellPoint
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Definition: interpolationCellPoint.H:50
Foam::functionObjects::streamLineBase::execute
virtual bool execute()
Execute the averaging.
Definition: streamLineBase.C:969
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::treeBoundBox::intersects
bool intersects(const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
Intersects segment; set point to intersection position and face,.
Definition: treeBoundBox.C:162
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::mapDistributeBase::distribute
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &, const negateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute data. Note:schedule only used for.
Definition: mapDistributeBaseTemplates.C:122
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
propsDict
IOdictionary propsDict(IOobject("particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED))
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::functionObjects::streamLineBase::initInterpolations
void initInterpolations(const label nSeeds, label &UIndex, PtrList< volScalarField > &vsFlds, PtrList< interpolation< scalar >> &vsInterp, PtrList< volVectorField > &vvFlds, PtrList< interpolation< vector >> &vvInterp)
Initialise fields, interpolators and track storage.
Definition: streamLineBase.C:142
Foam::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:376
Foam::PtrList::append
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:113
Foam::DynamicList::setSize
void setSize(const label n)
Same as resize()
Definition: DynamicList.H:224
Foam::interpolation< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::functionObjects::streamLineBase::trackDirTypeNames
static const Enum< trackDirType > trackDirTypeNames
Names for the trackDir.
Definition: streamLineBase.H:82
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::vectorList
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:64
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
ReadFields.H
Field reading functions for post-processing utilities.
Foam::treeBoundBox::contains
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:320
Foam::nl
constexpr char nl
Definition: Ostream.H:404
mapDistribute.H
Foam::UPstream::commsTypes::scheduled
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::functionObjects::streamLineBase::trimToBox
void trimToBox(const treeBoundBox &bb, const label tracki, PtrList< DynamicList< point >> &newTracks, PtrList< DynamicList< scalarList >> &newScalars, PtrList< DynamicList< vectorList >> &newVectors) const
Trim and possibly split a track.
Definition: streamLineBase.C:308
Foam::functionObjects::streamLineBase::sampledSetAxis_
word sampledSetAxis_
Axis of the sampled points to output.
Definition: streamLineBase.H:91
Foam::DynamicList::transfer
void transfer(List< T > &list)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:464
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::functionObjects::streamLineBase::streamLineBase
streamLineBase(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
Definition: streamLineBase.C:841
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::mapPolyMesh::mesh
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::flipOp
Functor to negate primitives. Dummy for most other types.
Definition: flipOp.H:53
Foam::functionObjects::streamLineBase::resetFieldNames
virtual void resetFieldNames(const word &newUName, const wordList &newFieldNames)
Reset the field names.
Definition: streamLineBase.C:828
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:38
Foam::functionObjects::streamLineBase::sampledSetAxis
const word & sampledSetAxis() const
The axis of the sampledSet. Creates sampledSet if required.
Definition: streamLineBase.C:66