streamlines.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) 2011-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 "Pstream.H"
27 #include "functionObjectList.H"
28 #include "streamlines.H"
29 #include "streamlinesCloud.H"
30 #include "ReadFields.H"
31 #include "meshSearch.H"
32 #include "sampledSet.H"
33 #include "globalIndex.H"
34 #include "mapDistribute.H"
35 #include "interpolationCellPoint.H"
36 #include "PatchTools.H"
37 #include "mapPolyMesh.H"
38 #include "writeFile.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  template<>
46  const char*
48  {"forward", "backward", "both"};
49 
50  namespace functionObjects
51  {
52  defineTypeNameAndDebug(streamlines, 0);
54 
57  }
58 }
59 
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
64 Foam::functionObjects::streamlines::wallPatch() const
65 {
66  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
67 
68  label nFaces = 0;
69 
70  forAll(patches, patchi)
71  {
72  if (isA<wallPolyPatch>(patches[patchi]))
73  {
74  nFaces += patches[patchi].size();
75  }
76  }
77 
78  labelList addressing(nFaces);
79 
80  nFaces = 0;
81 
82  forAll(patches, patchi)
83  {
84  if (isA<wallPolyPatch>(patches[patchi]))
85  {
86  const polyPatch& pp = patches[patchi];
87 
88  forAll(pp, i)
89  {
90  addressing[nFaces++] = pp.start()+i;
91  }
92  }
93  }
94 
96  (
98  (
100  (
101  mesh_.faces(),
102  addressing
103  ),
104  mesh_.points()
105  )
106  );
107 }
108 
109 
110 void Foam::functionObjects::streamlines::track()
111 {
112  IDLList<streamlinesParticle> initialParticles;
113  streamlinesCloud particles
114  (
115  mesh_,
116  cloudName_,
117  initialParticles
118  );
119 
120  const sampledSet& seedPoints = sampledSetPtr_();
121 
122  forAll(seedPoints, i)
123  {
124  particles.addParticle
125  (
127  (
128  mesh_,
129  seedPoints[i],
130  seedPoints.cells()[i],
131  lifeTime_
132  )
133  );
134  }
135 
136  label nSeeds = returnReduce(particles.size(), sumOp<label>());
137 
138  Info << " seeded " << nSeeds << " particles" << endl;
139 
140  // Read or lookup fields
145 
146  label UIndex = -1;
147 
148  label nScalar = 0;
149  label nVector = 0;
150 
151  forAll(fields_, i)
152  {
153  if (mesh_.foundObject<volScalarField>(fields_[i]))
154  {
155  nScalar++;
156  }
157  else if (mesh_.foundObject<volVectorField>(fields_[i]))
158  {
159  nVector++;
160  }
161  else
162  {
164  << "Cannot find field " << fields_[i] << nl
165  << "Valid scalar fields are:"
166  << mesh_.names(volScalarField::typeName) << nl
167  << "Valid vector fields are:"
168  << mesh_.names(volVectorField::typeName)
169  << exit(FatalError);
170  }
171  }
172  vsInterp.setSize(nScalar);
173  nScalar = 0;
174  vvInterp.setSize(nVector);
175  nVector = 0;
176 
177  forAll(fields_, i)
178  {
179  if (mesh_.foundObject<volScalarField>(fields_[i]))
180  {
181  const volScalarField& f = mesh_.lookupObject<volScalarField>
182  (
183  fields_[i]
184  );
185  vsInterp.set
186  (
187  nScalar++,
189  (
190  interpolationScheme_,
191  f
192  )
193  );
194  }
195  else if (mesh_.foundObject<volVectorField>(fields_[i]))
196  {
197  const volVectorField& f = mesh_.lookupObject<volVectorField>
198  (
199  fields_[i]
200  );
201 
202  if (f.name() == UName_)
203  {
204  UIndex = nVector;
205  }
206 
207  vvInterp.set
208  (
209  nVector++,
211  (
212  interpolationScheme_,
213  f
214  )
215  );
216  }
217  }
218 
219  // Store the names
220  scalarNames_.setSize(vsInterp.size());
221  forAll(vsInterp, i)
222  {
223  scalarNames_[i] = vsInterp[i].psi().name();
224  }
225  vectorNames_.setSize(vvInterp.size());
226  forAll(vvInterp, i)
227  {
228  vectorNames_[i] = vvInterp[i].psi().name();
229  }
230 
231  // Check that we know the index of U in the interpolators.
232 
233  if (UIndex == -1)
234  {
236  << "Cannot find field to move particles with : " << UName_ << nl
237  << "This field has to be present in the sampled fields " << fields_
238  << " and in the objectRegistry."
239  << exit(FatalError);
240  }
241 
242  // Sampled data
243  // ~~~~~~~~~~~~
244 
245  // Size to maximum expected sizes.
246  allTracks_.clear();
247  allTracks_.setCapacity(nSeeds);
248  allAges_.clear();
249  allAges_.setCapacity(nSeeds);
250  allScalars_.setSize(vsInterp.size());
251  forAll(allScalars_, i)
252  {
253  allScalars_[i].clear();
254  allScalars_[i].setCapacity(nSeeds);
255  }
256  allVectors_.setSize(vvInterp.size());
257  forAll(allVectors_, i)
258  {
259  allVectors_[i].clear();
260  allVectors_[i].setCapacity(nSeeds);
261  }
262 
263 
264  // Additional particle info
266  (
267  particles,
268  vsInterp,
269  vvInterp,
270  UIndex, // index of U in vvInterp
271 
272  trackDirection_ == trackDirection::forward,
273 
274  nSubCycle_, // automatic track control:step through cells in steps?
275  trackLength_, // fixed track length
276 
277  allTracks_,
278  allAges_,
279  allScalars_,
280  allVectors_
281  );
282 
283  // Set very large dt. Note: cannot use great since 1/great is small
284  // which is a trigger value for the tracking...
285  const scalar trackTime = Foam::sqrt(great);
286 
287  // Track
288  if (trackDirection_ == trackDirection::both)
289  {
290  initialParticles = particles;
291  }
292 
293  particles.move(particles, td, trackTime);
294 
295  if (trackDirection_ == trackDirection::both)
296  {
297  particles.IDLList<streamlinesParticle>::operator=(initialParticles);
298  td.trackForward_ = !td.trackForward_;
299  particles.move(particles, td, trackTime);
300  }
301 }
302 
303 
304 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
305 
307 (
308  const word& name,
309  const Time& runTime,
310  const dictionary& dict
311 )
312 :
313  fvMeshFunctionObject(name, runTime, dict),
314  dict_(dict),
315  nSubCycle_(0)
316 {
317  read(dict_);
318 }
319 
320 
321 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
322 
324 {}
325 
326 
327 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
328 
330 {
331  if (dict != dict_)
332  {
333  dict_ = dict;
334  }
335 
336  Info<< type() << " " << name() << ":" << nl;
337 
338  dict.lookup("fields") >> fields_;
339  UName_ = dict.lookupOrDefault("U", word("U"));
340 
341  if (findIndex(fields_, UName_) == -1)
342  {
344  << "Velocity field for tracking " << UName_
345  << " should be present in the list of fields " << fields_
346  << exit(FatalIOError);
347  }
348 
349  writeAge_ = dict.lookupOrDefault<Switch>("writeAge", true);
350 
351  // The trackForward entry is maintained here for backwards compatibility
352  if (!dict.found("direction") && dict.found("trackForward"))
353  {
354  trackDirection_ =
355  dict.lookup<bool>("trackForward")
356  ? trackDirection::forward
357  : trackDirection::backward;
358  }
359  else
360  {
361  trackDirection_ = trackDirectionNames_[word(dict.lookup("direction"))];
362  }
363 
364  dict.lookup("lifeTime") >> lifeTime_;
365  if (lifeTime_ < 1)
366  {
368  << "Illegal value " << lifeTime_ << " for lifeTime"
369  << exit(FatalError);
370  }
371 
372 
373  bool subCycling = dict.found("nSubCycle");
374  bool fixedLength = dict.found("trackLength");
375 
376  if (subCycling && fixedLength)
377  {
379  << "Cannot both specify automatic time stepping (through '"
380  << "nSubCycle' specification) and fixed track length (through '"
381  << "trackLength')"
382  << exit(FatalIOError);
383  }
384 
385 
386  nSubCycle_ = 1;
387  if (dict.readIfPresent("nSubCycle", nSubCycle_))
388  {
389  trackLength_ = vGreat;
390  if (nSubCycle_ < 1)
391  {
392  nSubCycle_ = 1;
393  }
394  Info<< " automatic track length specified through"
395  << " number of sub cycles : " << nSubCycle_ << nl << endl;
396  }
397  else
398  {
399  dict.lookup("trackLength") >> trackLength_;
400 
401  Info<< " fixed track length specified : "
402  << trackLength_ << nl << endl;
403  }
404 
405 
406  interpolationScheme_ = dict.lookupOrDefault
407  (
408  "interpolationScheme",
410  );
411 
412  cloudName_ = dict.lookupOrDefault<word>("cloudName", "streamlines");
413 
414  meshSearchPtr_.reset(new meshSearch(mesh_));
415 
416  sampledSetPtr_ = sampledSet::New
417  (
418  "seedSampleSet",
419  mesh_,
420  meshSearchPtr_(),
421  dict.subDict("seedSampleSet")
422  );
423  sampledSetAxis_ = sampledSetPtr_->axis();
424 
425  scalarFormatterPtr_ = setWriter<scalar>::New(dict.lookup("setFormat"));
426  vectorFormatterPtr_ = setWriter<vector>::New(dict.lookup("setFormat"));
427 
428  return true;
429 }
430 
431 
433 {
434  return true;
435 }
436 
437 
439 {
440  Info<< type() << " " << name() << " write:" << nl;
441 
442  const Time& runTime = obr_.time();
443 
444  // Do all injection and tracking
445  track();
446 
447 
448  if (Pstream::parRun())
449  {
450  // Append slave tracks to master ones
451  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
452 
453  globalIndex globalTrackIDs(allTracks_.size());
454 
455  // Construct a distribution map to pull all to the master.
456  labelListList sendMap(Pstream::nProcs());
457  labelListList recvMap(Pstream::nProcs());
458 
459  if (Pstream::master())
460  {
461  // Master: receive all. My own first, then consecutive
462  // processors.
463  label trackI = 0;
464 
465  forAll(recvMap, proci)
466  {
467  labelList& fromProc = recvMap[proci];
468  fromProc.setSize(globalTrackIDs.localSize(proci));
469  forAll(fromProc, i)
470  {
471  fromProc[i] = trackI++;
472  }
473  }
474  }
475 
476  labelList& toMaster = sendMap[0];
477  toMaster.setSize(globalTrackIDs.localSize());
478  forAll(toMaster, i)
479  {
480  toMaster[i] = i;
481  }
482 
483  const mapDistribute distMap
484  (
485  globalTrackIDs.size(),
486  move(sendMap),
487  move(recvMap)
488  );
489 
490 
491  // Distribute the track positions. Note: use scheduled comms
492  // to prevent buffering.
494  (
496  distMap.schedule(),
497  distMap.constructSize(),
498  distMap.subMap(),
499  false,
500  distMap.constructMap(),
501  false,
502  allTracks_,
503  flipOp()
504  );
505 
506  // Distribute the ages
508  (
510  distMap.schedule(),
511  distMap.constructSize(),
512  distMap.subMap(),
513  false,
514  distMap.constructMap(),
515  false,
516  allAges_,
517  flipOp()
518  );
519 
520  // Distribute the scalars
521  forAll(allScalars_, scalari)
522  {
523  allScalars_[scalari].shrink();
525  (
527  distMap.schedule(),
528  distMap.constructSize(),
529  distMap.subMap(),
530  false,
531  distMap.constructMap(),
532  false,
533  allScalars_[scalari],
534  flipOp()
535  );
536  allScalars_[scalari].setCapacity(allScalars_[scalari].size());
537  }
538  // Distribute the vectors
539  forAll(allVectors_, vectori)
540  {
541  allVectors_[vectori].shrink();
543  (
545  distMap.schedule(),
546  distMap.constructSize(),
547  distMap.subMap(),
548  false,
549  distMap.constructMap(),
550  false,
551  allVectors_[vectori],
552  flipOp()
553  );
554  allVectors_[vectori].setCapacity(allVectors_[vectori].size());
555  }
556  }
557 
558 
559  label n = 0;
560  forAll(allTracks_, trackI)
561  {
562  n += allTracks_[trackI].size();
563  }
564 
565  Info<< " Tracks:" << allTracks_.size() << nl
566  << " Total samples:" << n
567  << endl;
568 
569 
570  // Massage into form suitable for writers
571  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
572 
573  if (Pstream::master() && allTracks_.size())
574  {
575  // Make output directory
576 
577  fileName vtkPath
578  (
579  runTime.globalPath()/writeFile::outputPrefix/"sets"/name()
580  );
581  if (mesh_.name() != fvMesh::defaultRegion)
582  {
583  vtkPath = vtkPath/mesh_.name();
584  }
585  vtkPath = vtkPath/mesh_.time().timeName();
586  mkDir(vtkPath);
587 
588  // Convert track positions
589 
590  PtrList<coordSet> tracks(allTracks_.size());
591  forAll(allTracks_, trackI)
592  {
593  tracks.set
594  (
595  trackI,
596  new coordSet
597  (
598  "track" + Foam::name(trackI),
599  sampledSetAxis_ //"xyz"
600  )
601  );
602  tracks[trackI].transfer(allTracks_[trackI]);
603  }
604 
605  // Convert scalar values
606 
607  if (allScalars_.size() > 0 || writeAge_)
608  {
609  List<List<scalarField>> ageAndScalarValues
610  (
611  allScalars_.size() + writeAge_
612  );
613  wordList ageAndScalarNames(allScalars_.size() + writeAge_);
614 
615  if (writeAge_)
616  {
617  DynamicList<scalarList>& allTrackAges = allAges_;
618 
619  ageAndScalarValues[0].setSize(allTrackAges.size());
620  forAll(allTrackAges, trackI)
621  {
622  ageAndScalarValues[0][trackI].transfer
623  (
624  allTrackAges[trackI]
625  );
626  }
627 
628  ageAndScalarNames[0] = "age";
629  }
630 
631  forAll(allScalars_, scalari)
632  {
633  const label ageAndScalari = scalari + writeAge_;
634 
635  DynamicList<scalarList>& allTrackVals = allScalars_[scalari];
636 
637  ageAndScalarValues[ageAndScalari].setSize(allTrackVals.size());
638  forAll(allTrackVals, trackI)
639  {
640  ageAndScalarValues[ageAndScalari][trackI].transfer
641  (
642  allTrackVals[trackI]
643  );
644  }
645 
646  ageAndScalarNames[ageAndScalari] = scalarNames_[scalari];
647  }
648 
649  fileName vtkFile
650  (
651  vtkPath
652  / scalarFormatterPtr_().getFileName
653  (
654  tracks[0],
655  ageAndScalarNames
656  )
657  );
658 
659  Info<< " Writing data to " << vtkFile.path() << endl;
660 
661  scalarFormatterPtr_().write
662  (
663  true, // writeTracks
664  tracks,
665  ageAndScalarNames,
666  ageAndScalarValues,
667  OFstream(vtkFile)()
668  );
669  }
670 
671  // Convert vector values
672 
673  if (allVectors_.size() > 0)
674  {
675  List<List<vectorField>> vectorValues(allVectors_.size());
676 
677  forAll(allVectors_, vectori)
678  {
679  DynamicList<vectorList>& allTrackVals = allVectors_[vectori];
680 
681  vectorValues[vectori].setSize(allTrackVals.size());
682  forAll(allTrackVals, trackI)
683  {
684  vectorValues[vectori][trackI].transfer
685  (
686  allTrackVals[trackI]
687  );
688  }
689  }
690 
691  fileName vtkFile
692  (
693  vtkPath
694  / vectorFormatterPtr_().getFileName
695  (
696  tracks[0],
697  vectorNames_
698  )
699  );
700 
701  Info<< " Writing data to " << vtkFile.path() << endl;
702 
703  vectorFormatterPtr_().write
704  (
705  true, // writeTracks
706  tracks,
707  vectorNames_,
708  vectorValues,
709  OFstream(vtkFile)()
710  );
711  }
712  }
713 
714  return true;
715 }
716 
717 
719 {
720  if (&mpm.mesh() == &mesh_)
721  {
722  read(dict_);
723  }
724 }
725 
726 
728 {
729  if (&mesh == &mesh_)
730  {
731  // Moving mesh affects the search tree
732  read(dict_);
733  }
734 }
735 
736 
737 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
const word & name() const
Return name.
Definition: IOobject.H:303
A class for handling file names.
Definition: fileName.H:79
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:136
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
static const char *const typeName
Definition: Field.H:105
streamlines(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: streamlines.C:307
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:148
Field reading functions for post-processing utilities.
patches[0]
Abstract base-class for Time/database functionObjects.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
PtrList< volVectorField > vvFlds
Definition: readVolFields.H:6
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: streamlines.C:718
virtual bool execute()
Do nothing.
Definition: streamlines.C:432
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void setSize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:175
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:97
dynamicFvMesh & mesh
Holds list of sampling positions.
Definition: coordSet.H:49
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:61
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
A class for handling words, derived from string.
Definition: word.H:59
A Cloud of streamlines particles.
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())
Distribute data. Note:schedule only used for.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Particle class that samples fields as it passes through. Used in streamlines calculation.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:260
virtual bool write()
Calculate and write the streamlines.
Definition: streamlines.C:438
const Time & time() const
Return time.
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
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
defineTypeNameAndDebug(Qdot, 0)
static autoPtr< setWriter > New(const word &writeFormat)
Return a reference to the selected setWriter.
Definition: setWriter.C:35
fileName globalPath() const
Return the global path.
Definition: TimePaths.H:133
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: streamlines.C:727
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
virtual ~streamlines()
Destructor.
Definition: streamlines.C:323
void setSize(const label)
Reset size of List.
Definition: List.C:281
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:142
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Class containing processor-to-processor mapping information.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:352
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Abstract base class for interpolation.
messageStream Info
label n
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:265
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamlines.C:329
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
static const NamedEnum< trackDirection, 3 > trackDirectionNames_
Track direction enumeration names.
Definition: streamlines.H:219
const labelList & cells() const
Definition: sampledSet.H:209
A List with indirect addressing.
Definition: IndirectList.H:101
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
IOerror FatalIOError