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-2023 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 "OFstream.H"
33 #include "sampledSet.H"
34 #include "globalIndex.H"
35 #include "mapDistribute.H"
36 #include "interpolationCellPoint.H"
37 #include "wallPolyPatch.H"
38 #include "meshSearchMeshObject.H"
39 #include "mapPolyMesh.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace functionObjects
46 {
47  defineTypeNameAndDebug(streamLineBase, 0);
48 }
49 }
50 
51 
52 const Foam::Enum
53 <
55 >
57 ({
58  { trackDirType::FORWARD, "forward" },
59  { trackDirType::BACKWARD, "backward" },
60  { trackDirType::BIDIRECTIONAL, "bidirectional" }
61 });
62 
63 
64 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
65 
66 const Foam::word&
68 {
69  if (!sampledSetPtr_)
70  {
72  }
73 
74  return sampledSetAxis_;
75 }
76 
77 
78 const Foam::sampledSet&
80 {
81  if (!sampledSetPtr_)
82  {
83  sampledSetPtr_ = sampledSet::New
84  (
85  "seedSampleSet",
86  mesh_,
88  dict_.subDict("seedSampleSet")
89  );
90 
91  sampledSetAxis_ = sampledSetPtr_->axis();
92  }
93 
94  return *sampledSetPtr_;
95 }
96 
97 
100 {
101  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
102 
103  label nFaces = 0;
104 
105  for (const polyPatch& pp : patches)
106  {
107  //if (!polyPatch::constraintType(pp.type()))
108  if (isA<wallPolyPatch>(pp))
109  {
110  nFaces += pp.size();
111  }
112  }
113 
114  labelList addressing(nFaces);
115 
116  nFaces = 0;
117 
118  for (const polyPatch& pp : patches)
119  {
120  //if (!polyPatch::constraintType(pp.type()))
121  if (isA<wallPolyPatch>(pp))
122  {
123  forAll(pp, i)
124  {
125  addressing[nFaces++] = pp.start()+i;
126  }
127  }
128  }
129 
131  (
132  IndirectList<face>
133  (
134  mesh_.faces(),
135  addressing
136  ),
137  mesh_.points()
138  );
139 }
140 
141 
144 (
145  const label nSeeds,
146  PtrList<interpolation<scalar>>& vsInterp,
147  PtrList<interpolation<vector>>& vvInterp
148 )
149 {
151 
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 scalar/vector field " << fieldName << nl
169  << "Valid scalar fields: "
170  << flatOutput(mesh_.sortedNames<volScalarField>()) << nl
171  << "Valid vector fields: "
172  << flatOutput(mesh_.sortedNames<volVectorField>()) << nl
173  << exit(FatalError);
174  }
175  }
176  vsInterp.resize(nScalar);
177  vvInterp.resize(nVector);
178  nScalar = 0;
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,
189  interpolation<scalar>::New(interpolationScheme_, f)
190  );
191  ++nScalar;
192  }
193  else if (foundObject<volVectorField>(fieldName))
194  {
195  const volVectorField& f = lookupObject<volVectorField>(fieldName);
196 
197  vvInterp.set
198  (
199  nVector,
200  interpolation<vector>::New(interpolationScheme_, f)
201  );
202 
203  if (f.name() == UName_)
204  {
205  // Velocity is part of sampled velocity fields
206  UInterp.cref(vvInterp[nVector]);
207  }
208 
209  ++nVector;
210  }
211  }
212 
213  if (!UInterp)
214  {
215  // Velocity was not in sampled velocity fields
216  UInterp.reset
217  (
219  (
220  interpolationScheme_,
221  // Fatal if missing
222  lookupObject<volVectorField>(UName_)
223  )
224  );
225  }
226 
227  // Store the names
228  scalarNames_.resize(vsInterp.size());
229  forAll(vsInterp, i)
230  {
231  scalarNames_[i] = vsInterp[i].psi().name();
232  }
233  vectorNames_.resize(vvInterp.size());
234  forAll(vvInterp, i)
235  {
236  vectorNames_[i] = vvInterp[i].psi().name();
237  }
238 
239  // Sampled data
240  // ~~~~~~~~~~~~
241 
242  // Size to maximum expected sizes.
243  allTracks_.clear();
244  allTracks_.setCapacity(nSeeds);
245  allScalars_.resize(vsInterp.size());
246  forAll(allScalars_, i)
247  {
248  allScalars_[i].clear();
249  allScalars_[i].setCapacity(nSeeds);
250  }
251  allVectors_.resize(vvInterp.size());
252  forAll(allVectors_, i)
253  {
254  allVectors_[i].clear();
255  allVectors_[i].setCapacity(nSeeds);
256  }
257 
258  return UInterp;
259 }
260 
261 
263 (
264  const label tracki,
265 
266  const scalar w,
267  const label lefti,
268  const label righti,
269 
270  DynamicList<point>& newTrack,
271  DynamicList<scalarList>& newScalars,
272  DynamicList<vectorList>& newVectors
273 ) const
274 {
275  const List<point>& track = allTracks_[tracki];
276 
277  newTrack.push_back(lerp(track[lefti], track[righti], w));
278 
279  // Scalars
280  {
281  scalarList& newVals = newScalars.emplace_back(allScalars_.size());
282 
283  forAll(allScalars_, i)
284  {
285  const scalarList& trackVals = allScalars_[i][tracki];
286  newVals[i] = lerp(trackVals[lefti], trackVals[righti], w);
287  }
288  }
289 
290  // Vectors
291  {
292  vectorList& newVals = newVectors.emplace_back(allVectors_.size());
293 
294  forAll(allVectors_, i)
295  {
296  const vectorList& trackVals = allVectors_[i][tracki];
297  newVals[i] = lerp(trackVals[lefti], trackVals[righti], w);
298  }
299  }
300 }
301 
302 
303 // Can split a track into multiple tracks
305 (
306  const treeBoundBox& bb,
307  const label tracki,
308  PtrList<DynamicList<point>>& newTracks,
309  PtrList<DynamicList<scalarList>>& newScalars,
310  PtrList<DynamicList<vectorList>>& newVectors
311 ) const
312 {
313  const List<point>& track = allTracks_[tracki];
314 
315  if (track.size())
316  {
317  for
318  (
319  label segmenti = 1;
320  segmenti < track.size();
321  segmenti++
322  )
323  {
324  const point& startPt = track[segmenti-1];
325  const point& endPt = track[segmenti];
326 
327  const vector d(endPt-startPt);
328  const scalar magD = mag(d);
329  if (magD > ROOTVSMALL)
330  {
331  if (bb.contains(startPt))
332  {
333  // Store 1.0*track[segmenti-1]+0*track[segmenti]
334  storePoint
335  (
336  tracki,
337 
338  0.0,
339  segmenti-1,
340  segmenti,
341 
342  newTracks.back(),
343  newScalars.back(),
344  newVectors.back()
345  );
346 
347  if (!bb.contains(endPt))
348  {
349  point clipPt;
350  if (bb.intersects(endPt, startPt, clipPt))
351  {
352  // End of track. Store point and interpolated
353  // values
354  storePoint
355  (
356  tracki,
357 
358  mag(clipPt-startPt)/magD,
359  segmenti-1,
360  segmenti,
361 
362  newTracks.back(),
363  newScalars.back(),
364  newVectors.back()
365  );
366 
367  newTracks.back().shrink();
368  newScalars.back().shrink();
369  newVectors.back().shrink();
370  }
371  }
372  }
373  else
374  {
375  // startPt outside box. New track. Get starting point
376 
377  point clipPt;
378  if (bb.intersects(startPt, endPt, clipPt))
379  {
380  // New track
381  const label defltCapacity(track.size()/10);
382 
383  // Store point and interpolated values
384  storePoint
385  (
386  tracki,
387 
388  mag(clipPt-startPt)/magD,
389  segmenti-1,
390  segmenti,
391 
392  newTracks.emplace_back(defltCapacity),
393  newScalars.emplace_back(defltCapacity),
394  newVectors.emplace_back(defltCapacity)
395  );
396 
397  if (!bb.contains(endPt))
398  {
399  bb.intersects
400  (
401  endPt,
402  point(clipPt),
403  clipPt
404  );
405 
406  // Store point and interpolated values
407  storePoint
408  (
409  tracki,
410 
411  mag(clipPt-startPt)/magD,
412  segmenti-1,
413  segmenti,
414 
415  newTracks.back(),
416  newScalars.back(),
417  newVectors.back()
418  );
419 
420  newTracks.back().shrink();
421  newScalars.back().shrink();
422  newVectors.back().shrink();
423  }
424  }
425  }
426  }
427  }
428 
429  // Last point
430  if (bb.contains(track.back()))
431  {
432  storePoint
433  (
434  tracki,
435 
436  1.0,
437  track.size()-2,
438  track.size()-1,
439 
440  newTracks.back(),
441  newScalars.back(),
442  newVectors.back()
443  );
444  }
445  }
446 }
447 
448 
449 void Foam::functionObjects::streamLineBase::trimToBox(const treeBoundBox& bb)
450 {
451  // Storage for new tracks. Per track, per sample the coordinate (newTracks)
452  // or values for all the sampled fields (newScalars, newVectors)
453  PtrList<DynamicList<point>> newTracks;
454  PtrList<DynamicList<scalarList>> newScalars;
455  PtrList<DynamicList<vectorList>> newVectors;
456 
457  forAll(allTracks_, tracki)
458  {
459  const List<point>& track = allTracks_[tracki];
460 
461  if (track.size())
462  {
463  // New track. Assume it consists of the whole track
464  newTracks.emplace_back(track.size());
465  newScalars.emplace_back(track.size());
466  newVectors.emplace_back(track.size());
467 
468  // Trim, split and append to newTracks
469  trimToBox(bb, tracki, newTracks, newScalars, newVectors);
470  }
471  }
472 
473  // Transfer newTracks to allTracks_
474  allTracks_.setSize(newTracks.size());
475  forAll(allTracks_, tracki)
476  {
477  allTracks_[tracki].transfer(newTracks[tracki]);
478  }
479  // Replace track scalars
480  forAll(allScalars_, scalari)
481  {
482  DynamicList<scalarList>& fieldVals = allScalars_[scalari];
483  fieldVals.setSize(newTracks.size());
484 
485  forAll(fieldVals, tracki)
486  {
487  scalarList& trackVals = allScalars_[scalari][tracki];
488  trackVals.setSize(newScalars[tracki].size());
489  forAll(trackVals, samplei)
490  {
491  trackVals[samplei] = newScalars[tracki][samplei][scalari];
492  }
493  }
494  }
495  // Replace track vectors
496  forAll(allVectors_, vectori)
497  {
498  DynamicList<vectorList>& fieldVals = allVectors_[vectori];
499  fieldVals.setSize(newTracks.size());
500  forAll(fieldVals, tracki)
501  {
502  vectorList& trackVals = allVectors_[vectori][tracki];
503  trackVals.setSize(newVectors[tracki].size());
504  forAll(trackVals, samplei)
505  {
506  trackVals[samplei] = newVectors[tracki][samplei][vectori];
507  }
508  }
509  }
510 }
511 
512 
514 {
515  if (Pstream::parRun())
516  {
517  // Append slave tracks to master ones
518  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
519 
520  globalIndex globalTrackIDs(allTracks_.size());
521 
522  // Construct a distribution map to pull all to the master.
523  labelListList sendMap(Pstream::nProcs());
524  labelListList recvMap(Pstream::nProcs());
525 
526  if (Pstream::master())
527  {
528  // Master: receive all. My own first, then consecutive
529  // processors.
530  label tracki = 0;
531 
532  forAll(recvMap, proci)
533  {
534  labelList& fromProc = recvMap[proci];
535  fromProc.setSize(globalTrackIDs.localSize(proci));
536  forAll(fromProc, i)
537  {
538  fromProc[i] = tracki++;
539  }
540  }
541  }
542 
543  labelList& toMaster = sendMap[0];
544  toMaster.setSize(globalTrackIDs.localSize());
545  forAll(toMaster, i)
546  {
547  toMaster[i] = i;
548  }
549 
550  const mapDistribute distMap
551  (
552  globalTrackIDs.totalSize(),
553  std::move(sendMap),
554  std::move(recvMap)
555  );
556 
557 
558  // Distribute the track positions. Note: use scheduled comms
559  // to prevent buffering.
560  allTracks_.shrink();
562  (
564  distMap.schedule(),
565  distMap.constructSize(),
566  distMap.subMap(),
567  false,
568  distMap.constructMap(),
569  false,
570  allTracks_,
571  flipOp()
572  );
573  allTracks_.setCapacity(allTracks_.size());
574 
575  // Distribute the scalars
576  forAll(allScalars_, scalari)
577  {
578  allScalars_[scalari].shrink();
580  (
582  distMap.schedule(),
583  distMap.constructSize(),
584  distMap.subMap(),
585  false,
586  distMap.constructMap(),
587  false,
588  allScalars_[scalari],
589  flipOp()
590  );
591  allScalars_[scalari].setCapacity(allScalars_[scalari].size());
592  }
593  // Distribute the vectors
594  forAll(allVectors_, vectori)
595  {
596  allVectors_[vectori].shrink();
598  (
600  distMap.schedule(),
601  distMap.constructSize(),
602  distMap.subMap(),
603  false,
604  distMap.constructMap(),
605  false,
606  allVectors_[vectori],
607  flipOp()
608  );
609  allVectors_[vectori].setCapacity(allVectors_[vectori].size());
610  }
611  }
612 
613 
614  // Note: filenames scattered below since used in global call
615  HashTable<fileName> outputFileNames;
616 
617  if (Pstream::master())
618  {
619  if (!bounds_.empty())
620  {
621  // Clip to bounding box
622  trimToBox(treeBoundBox(bounds_));
623  }
624 
625 
626  label nTracks = 0;
627  label n = 0;
628  forAll(allTracks_, tracki)
629  {
630  if (allTracks_[tracki].size())
631  {
632  nTracks++;
633  n += allTracks_[tracki].size();
634  }
635  }
636 
637  Log << " Tracks:" << nTracks << nl
638  << " Total samples:" << n
639  << endl;
640 
641 
642  // Massage into form suitable for writers
643  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
644 
645  // Make output directory
646 
647  fileName vtkPath
648  (
649  time_.globalPath()/functionObject::outputPrefix/"sets"
650  / name()/mesh_.regionName()
651  / mesh_.time().timeName()
652  );
653 
654  mkDir(vtkPath);
655 
656  // Convert track positions (and compact out empty tracks)
657 
658  PtrList<coordSet> tracks(nTracks);
659  nTracks = 0;
660  labelList oldToNewTrack(allTracks_.size(), -1);
661 
662  forAll(allTracks_, tracki)
663  {
664  if (allTracks_[tracki].size())
665  {
666  List<point>& points = allTracks_[tracki];
667  scalarList dist(points.size());
668  dist[0] = 0;
669  for (label pointi = 1; pointi < points.size(); ++pointi)
670  {
671  dist[pointi] =
672  dist[pointi-1] + mag(points[pointi] - points[pointi-1]);
673  }
674 
675  tracks.set
676  (
677  nTracks,
678  new coordSet
679  (
680  "track" + Foam::name(nTracks),
681  sampledSetAxis(), // "xyz"
682  std::move(allTracks_[tracki]),
683  std::move(dist)
684  )
685  );
686  oldToNewTrack[tracki] = nTracks;
687  ++nTracks;
688  }
689  }
690 
691 
692  const bool canWrite =
693  (
694  !tracks.empty()
695  && trackWriterPtr_
696  && trackWriterPtr_->enabled()
697  && (!allScalars_.empty() || !allVectors_.empty())
698  );
699 
700  if (canWrite)
701  {
702  auto& writer = trackWriterPtr_();
703 
704  writer.nFields(allScalars_.size() + allVectors_.size());
705 
706  writer.open
707  (
708  tracks,
709  (vtkPath / tracks[0].name())
710  );
711 
712 
713  // Temporary measure
714  if (!allScalars_.empty())
715  {
716  List<List<scalarField>> scalarValues(allScalars_.size());
717 
718  forAll(allScalars_, scalari)
719  {
720  auto& allTrackVals = allScalars_[scalari];
721  scalarValues[scalari].resize(nTracks);
722 
723  forAll(allTrackVals, tracki)
724  {
725  scalarList& vals = allTrackVals[tracki];
726  if (vals.size())
727  {
728  const label newTracki = oldToNewTrack[tracki];
729  scalarValues[scalari][newTracki].transfer(vals);
730  }
731  }
732  }
733 
734  forAll(scalarNames_, i)
735  {
736  fileName outFile =
737  writer.write(scalarNames_[i], scalarValues[i]);
738 
739  outputFileNames.insert
740  (
741  scalarNames_[i],
742  time_.relativePath(outFile, true)
743  );
744  }
745  }
746 
747  if (!allVectors_.empty())
748  {
749  List<List<vectorField>> vectorValues(allVectors_.size());
750 
751  forAll(allVectors_, vectori)
752  {
753  auto& allTrackVals = allVectors_[vectori];
754  vectorValues[vectori].setSize(nTracks);
755 
756  forAll(allTrackVals, tracki)
757  {
758  vectorList& vals = allTrackVals[tracki];
759  if (vals.size())
760  {
761  const label newTracki = oldToNewTrack[tracki];
762  vectorValues[vectori][newTracki].transfer(vals);
763  }
764  }
765  }
766 
767  forAll(vectorNames_, i)
768  {
769  fileName outFile =
770  writer.write(vectorNames_[i], vectorValues[i]);
771 
772  outputFileNames.insert
773  (
774  vectorNames_[i],
775  time_.relativePath(outFile, true)
776  );
777  }
778  }
779 
780  writer.close(true);
781  }
782 
783  // Log << " Writing data to " << scalarVtkFile.path() << endl;
784  }
785 
786 
787  // File names generated on the master but setProperty needed everywher
788  Pstream::broadcast(outputFileNames);
789 
790  forAllConstIters(outputFileNames, iter)
791  {
792  const word& fieldName = iter.key();
793  const fileName& outputName = iter.val();
794 
796  propsDict.add("file", outputName);
797  setProperty(fieldName, propsDict);
798  }
799 
800  return true;
801 }
802 
803 
805 (
806  const word& newUName,
807  const wordList& newFieldNames
808 )
809 {
810  UName_ = newUName;
811  fields_ = newFieldNames;
812 }
813 
814 
815 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
816 
818 (
819  const word& name,
820  const Time& runTime,
821  const dictionary& dict
822 )
823 :
824  functionObjects::fvMeshFunctionObject(name, runTime, dict),
825  dict_(dict),
826  fields_()
827 {}
828 
829 
831 (
832  const word& name,
833  const Time& runTime,
834  const dictionary& dict,
835  const wordList& fieldNames
836 )
837 :
838  functionObjects::fvMeshFunctionObject(name, runTime, dict),
839  dict_(dict),
840  fields_(fieldNames)
841 {}
842 
843 
844 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
847 {}
848 
849 
850 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
851 
853 {
854  if (&dict_ != &dict)
855  {
856  // Update local copy of dictionary:
857  dict_ = dict;
858  }
859 
861 
862  Info<< type() << " " << name() << ":" << nl;
863 
864  UName_ = dict.getOrDefault<word>("U", "U");
865  Info<< " Employing velocity field " << UName_ << endl;
866 
867  if (fields_.empty())
868  {
869  dict.readEntry("fields", fields_);
870  }
871 
872  bool trackForward;
873  if (dict.readIfPresent("trackForward", trackForward))
874  {
875  trackDir_ =
876  (
877  trackForward
878  ? trackDirType::FORWARD
879  : trackDirType::BACKWARD
880  );
881 
882  if (dict.found("direction"))
883  {
885  << "Cannot specify both 'trackForward' and 'direction'" << nl
886  << exit(FatalIOError);
887  }
888  }
889  else
890  {
891  trackDir_ = trackDirTypeNames.getOrDefault
892  (
893  "direction",
894  dict,
895  trackDirType::FORWARD
896  );
897  }
898  dict.readEntry("lifeTime", lifeTime_);
899  if (lifeTime_ < 1)
900  {
902  << "Illegal value " << lifeTime_ << " for lifeTime"
903  << exit(FatalIOError);
904  }
905 
906 
907  trackLength_ = VGREAT;
908  if (dict.readIfPresent("trackLength", trackLength_))
909  {
910  Info<< type() << " : fixed track length specified : "
911  << trackLength_ << nl << endl;
912  }
913 
914 
915  bounds_.reset();
916  if (dict.readIfPresent("bounds", bounds_) && !bounds_.empty())
917  {
918  Info<< " clipping all segments to " << bounds_ << nl << endl;
919  }
920 
921 
922  interpolationScheme_ = dict.getOrDefault
923  (
924  "interpolationScheme",
925  interpolationCellPoint<scalar>::typeName
926  );
927 
928  //Info<< " using interpolation " << interpolationScheme_ << endl;
929 
930  cloudName_ = dict.getOrDefault<word>("cloud", type());
931 
932  sampledSetPtr_.clear();
933  sampledSetAxis_.clear();
934 
935  const word setFormat(dict.get<word>("setFormat"));
936 
937  trackWriterPtr_ = coordSetWriter::New
938  (
939  setFormat,
940  dict.subOrEmptyDict("formatOptions").optionalSubDict(setFormat)
941  );
942 
943  return true;
944 }
945 
948 {
949  return true;
950 }
951 
952 
954 {
955  Log << type() << " " << name() << " write:" << nl;
956 
957  // Do all injection and tracking
958  track();
960  writeToFile();
961 
962  return true;
963 }
964 
965 
967 {
968  if (&mpm.mesh() == &mesh_)
969  {
970  read(dict_);
971  }
972 }
973 
974 
976 {
977  // Moving mesh affects the search tree
978  read(dict_);
979 }
980 
981 
982 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
streamLineBase(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
dictionary dict
const sampledSet & sampledSetPoints() const
Demand driven construction of the sampledSet.
defineTypeNameAndDebug(ObukhovLength, 0)
word sampledSetAxis_
Axis of the sampled points to output.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual bool writeToFile()
Write tracks to file.
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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:477
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
static const meshSearchMeshObject & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
trackDirType
Enumeration defining the track direction.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
refPtr< interpolation< vector > > initInterpolations(const label nSeeds, PtrList< interpolation< scalar >> &vsInterp, PtrList< interpolation< vector >> &vvInterp)
Initialise interpolators and track storage.
Field reading functions for post-processing utilities.
List< vector > vectorList
List of vector.
Definition: vectorList.H:32
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: refPtrI.H:216
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
IOdictionary propsDict(dictIO)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static autoPtr< coordSetWriter > New(const word &writeFormat)
Return a reference to the selected writer.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
word outputName("finiteArea-edges.obj")
void setSize(const label n)
Alias for resize()
Definition: List.H:316
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
virtual bool write()
Track and write.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
"scheduled" : (MPI_Send, MPI_Recv)
const pointField & points
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:79
virtual bool read(const dictionary &)
Read the field average data.
A class for handling words, derived from Foam::string.
Definition: word.H:63
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.
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
Definition: DynamicListI.H:538
static const Enum< trackDirType > trackDirTypeNames
Names for the trackDir.
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.
Vector< scalar > vector
Definition: vector.H:57
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual bool execute()
Execute the averaging.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
autoPtr< sampledSet > sampledSetPtr_
Seed set engine.
labelList f(nPoints)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
void push_back(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:555
virtual void resetFieldNames(const word &newUName, const wordList &newFieldNames)
Reset the field names.
vector point
Point is a vector.
Definition: point.H:37
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:437
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
#define Log
Definition: PDRblock.C:28
static word outputPrefix
Directory prefix.
T & back()
Access last element of the list, position [size()-1].
Definition: UListI.H:251
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
label n
const word & sampledSetAxis() const
The axis of the sampledSet. Creates sampledSet if required.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))
List< label > labelList
A List of labels.
Definition: List.H:62
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
autoPtr< indirectPrimitivePatch > wallPatch() const
Construct patch out of all wall patch faces.
static void distribute(const UPstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &field, const T &nullValue, const CombineOp &cop, const NegateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute combine data with specified combine operation and negate operator (for flips)...
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
void reset(T *p=nullptr) noexcept
Delete managed pointer and set to new given pointer.
Definition: refPtrI.H:314
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...