streamLine.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-2019 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 "streamLine.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"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  template<>
45  const char*
47  {"forward", "backward", "both"};
48 
49  namespace functionObjects
50  {
51  defineTypeNameAndDebug(streamLine, 0);
53 
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61 
63 Foam::functionObjects::streamLine::wallPatch() const
64 {
65  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
66 
67  label nFaces = 0;
68 
69  forAll(patches, patchi)
70  {
71  if (isA<wallPolyPatch>(patches[patchi]))
72  {
73  nFaces += patches[patchi].size();
74  }
75  }
76 
77  labelList addressing(nFaces);
78 
79  nFaces = 0;
80 
81  forAll(patches, patchi)
82  {
83  if (isA<wallPolyPatch>(patches[patchi]))
84  {
85  const polyPatch& pp = patches[patchi];
86 
87  forAll(pp, i)
88  {
89  addressing[nFaces++] = pp.start()+i;
90  }
91  }
92  }
93 
95  (
97  (
99  (
100  mesh_.faces(),
101  addressing
102  ),
103  mesh_.points()
104  )
105  );
106 }
107 
108 
109 void Foam::functionObjects::streamLine::track()
110 {
111  IDLList<streamLineParticle> initialParticles;
112  streamLineParticleCloud particles
113  (
114  mesh_,
115  cloudName_,
116  initialParticles
117  );
118 
119  const sampledSet& seedPoints = sampledSetPtr_();
120 
121  forAll(seedPoints, i)
122  {
123  particles.addParticle
124  (
126  (
127  mesh_,
128  seedPoints[i],
129  seedPoints.cells()[i],
130  lifeTime_
131  )
132  );
133  }
134 
135  label nSeeds = returnReduce(particles.size(), sumOp<label>());
136 
137  Info << " seeded " << nSeeds << " particles" << endl;
138 
139  // Read or lookup fields
144 
145  label UIndex = -1;
146 
147  label nScalar = 0;
148  label nVector = 0;
149 
150  forAll(fields_, i)
151  {
152  if (mesh_.foundObject<volScalarField>(fields_[i]))
153  {
154  nScalar++;
155  }
156  else if (mesh_.foundObject<volVectorField>(fields_[i]))
157  {
158  nVector++;
159  }
160  else
161  {
163  << "Cannot find field " << fields_[i] << nl
164  << "Valid scalar fields are:"
165  << mesh_.names(volScalarField::typeName) << nl
166  << "Valid vector fields are:"
167  << mesh_.names(volVectorField::typeName)
168  << exit(FatalError);
169  }
170  }
171  vsInterp.setSize(nScalar);
172  nScalar = 0;
173  vvInterp.setSize(nVector);
174  nVector = 0;
175 
176  forAll(fields_, i)
177  {
178  if (mesh_.foundObject<volScalarField>(fields_[i]))
179  {
180  const volScalarField& f = mesh_.lookupObject<volScalarField>
181  (
182  fields_[i]
183  );
184  vsInterp.set
185  (
186  nScalar++,
188  (
189  interpolationScheme_,
190  f
191  )
192  );
193  }
194  else if (mesh_.foundObject<volVectorField>(fields_[i]))
195  {
196  const volVectorField& f = mesh_.lookupObject<volVectorField>
197  (
198  fields_[i]
199  );
200 
201  if (f.name() == UName_)
202  {
203  UIndex = nVector;
204  }
205 
206  vvInterp.set
207  (
208  nVector++,
210  (
211  interpolationScheme_,
212  f
213  )
214  );
215  }
216  }
217 
218  // Store the names
219  scalarNames_.setSize(vsInterp.size());
220  forAll(vsInterp, i)
221  {
222  scalarNames_[i] = vsInterp[i].psi().name();
223  }
224  vectorNames_.setSize(vvInterp.size());
225  forAll(vvInterp, i)
226  {
227  vectorNames_[i] = vvInterp[i].psi().name();
228  }
229 
230  // Check that we know the index of U in the interpolators.
231 
232  if (UIndex == -1)
233  {
235  << "Cannot find field to move particles with : " << UName_ << nl
236  << "This field has to be present in the sampled fields " << fields_
237  << " and in the objectRegistry."
238  << exit(FatalError);
239  }
240 
241  // Sampled data
242  // ~~~~~~~~~~~~
243 
244  // Size to maximum expected sizes.
245  allTracks_.clear();
246  allTracks_.setCapacity(nSeeds);
247  allScalars_.setSize(vsInterp.size());
248  forAll(allScalars_, i)
249  {
250  allScalars_[i].clear();
251  allScalars_[i].setCapacity(nSeeds);
252  }
253  allVectors_.setSize(vvInterp.size());
254  forAll(allVectors_, i)
255  {
256  allVectors_[i].clear();
257  allVectors_[i].setCapacity(nSeeds);
258  }
259 
260 
261  // Additional particle info
263  (
264  particles,
265  vsInterp,
266  vvInterp,
267  UIndex, // index of U in vvInterp
268 
269  trackDirection_ == trackDirection::forward,
270 
271  nSubCycle_, // automatic track control:step through cells in steps?
272  trackLength_, // fixed track length
273 
274  allTracks_,
275  allScalars_,
276  allVectors_
277  );
278 
279  // Set very large dt. Note: cannot use great since 1/great is small
280  // which is a trigger value for the tracking...
281  const scalar trackTime = Foam::sqrt(great);
282 
283  // Track
284  if (trackDirection_ == trackDirection::both)
285  {
286  initialParticles = particles;
287  }
288 
289  particles.move(particles, td, trackTime);
290 
291  if (trackDirection_ == trackDirection::both)
292  {
293  particles.IDLList<streamLineParticle>::operator=(initialParticles);
294  td.trackForward_ = !td.trackForward_;
295  particles.move(particles, td, trackTime);
296  }
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
301 
303 (
304  const word& name,
305  const Time& runTime,
306  const dictionary& dict
307 )
308 :
309  fvMeshFunctionObject(name, runTime, dict),
310  dict_(dict),
311  nSubCycle_(0)
312 {
313  read(dict_);
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
318 
320 {}
321 
322 
323 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
324 
326 {
327  if (dict != dict_)
328  {
329  dict_ = dict;
330  }
331 
332  Info<< type() << " " << name() << ":" << nl;
333 
334  dict.lookup("fields") >> fields_;
335  UName_ = dict.lookupOrDefault("U", word("U"));
336 
337  if (findIndex(fields_, UName_) == -1)
338  {
340  << "Velocity field for tracking " << UName_
341  << " should be present in the list of fields " << fields_
342  << exit(FatalIOError);
343  }
344 
345  // The trackForward entry is maintained here for backwards compatibility
346  if (!dict.found("direction") && dict.found("trackForward"))
347  {
348  trackDirection_ =
349  dict.lookup<bool>("trackForward")
350  ? trackDirection::forward
351  : trackDirection::backward;
352  }
353  else
354  {
355  trackDirection_ = trackDirectionNames_[word(dict.lookup("direction"))];
356  }
357 
358  dict.lookup("lifeTime") >> lifeTime_;
359  if (lifeTime_ < 1)
360  {
362  << "Illegal value " << lifeTime_ << " for lifeTime"
363  << exit(FatalError);
364  }
365 
366 
367  bool subCycling = dict.found("nSubCycle");
368  bool fixedLength = dict.found("trackLength");
369 
370  if (subCycling && fixedLength)
371  {
373  << "Cannot both specify automatic time stepping (through '"
374  << "nSubCycle' specification) and fixed track length (through '"
375  << "trackLength')"
376  << exit(FatalIOError);
377  }
378 
379 
380  nSubCycle_ = 1;
381  if (dict.readIfPresent("nSubCycle", nSubCycle_))
382  {
383  trackLength_ = vGreat;
384  if (nSubCycle_ < 1)
385  {
386  nSubCycle_ = 1;
387  }
388  Info<< " automatic track length specified through"
389  << " number of sub cycles : " << nSubCycle_ << nl << endl;
390  }
391  else
392  {
393  dict.lookup("trackLength") >> trackLength_;
394 
395  Info<< " fixed track length specified : "
396  << trackLength_ << nl << endl;
397  }
398 
399 
400  interpolationScheme_ = dict.lookupOrDefault
401  (
402  "interpolationScheme",
404  );
405 
406  cloudName_ = dict.lookupOrDefault<word>("cloudName", "streamLine");
407 
408  meshSearchPtr_.reset(new meshSearch(mesh_));
409 
410  sampledSetPtr_ = sampledSet::New
411  (
412  "seedSampleSet",
413  mesh_,
414  meshSearchPtr_(),
415  dict.subDict("seedSampleSet")
416  );
417  sampledSetAxis_ = sampledSetPtr_->axis();
418 
419  scalarFormatterPtr_ = writer<scalar>::New(dict.lookup("setFormat"));
420  vectorFormatterPtr_ = writer<vector>::New(dict.lookup("setFormat"));
421 
422  return true;
423 }
424 
425 
427 {
428  return true;
429 }
430 
431 
433 {
434  Info<< type() << " " << name() << " write:" << nl;
435 
436  const Time& runTime = obr_.time();
437 
438  // Do all injection and tracking
439  track();
440 
441 
442  if (Pstream::parRun())
443  {
444  // Append slave tracks to master ones
445  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
446 
447  globalIndex globalTrackIDs(allTracks_.size());
448 
449  // Construct a distribution map to pull all to the master.
450  labelListList sendMap(Pstream::nProcs());
451  labelListList recvMap(Pstream::nProcs());
452 
453  if (Pstream::master())
454  {
455  // Master: receive all. My own first, then consecutive
456  // processors.
457  label trackI = 0;
458 
459  forAll(recvMap, proci)
460  {
461  labelList& fromProc = recvMap[proci];
462  fromProc.setSize(globalTrackIDs.localSize(proci));
463  forAll(fromProc, i)
464  {
465  fromProc[i] = trackI++;
466  }
467  }
468  }
469 
470  labelList& toMaster = sendMap[0];
471  toMaster.setSize(globalTrackIDs.localSize());
472  forAll(toMaster, i)
473  {
474  toMaster[i] = i;
475  }
476 
477  const mapDistribute distMap
478  (
479  globalTrackIDs.size(),
480  move(sendMap),
481  move(recvMap)
482  );
483 
484 
485  // Distribute the track positions. Note: use scheduled comms
486  // to prevent buffering.
488  (
490  distMap.schedule(),
491  distMap.constructSize(),
492  distMap.subMap(),
493  false,
494  distMap.constructMap(),
495  false,
496  allTracks_,
497  flipOp()
498  );
499 
500  // Distribute the scalars
501  forAll(allScalars_, scalarI)
502  {
503  allScalars_[scalarI].shrink();
505  (
507  distMap.schedule(),
508  distMap.constructSize(),
509  distMap.subMap(),
510  false,
511  distMap.constructMap(),
512  false,
513  allScalars_[scalarI],
514  flipOp()
515  );
516  allScalars_[scalarI].setCapacity(allScalars_[scalarI].size());
517  }
518  // Distribute the vectors
519  forAll(allVectors_, vectorI)
520  {
521  allVectors_[vectorI].shrink();
523  (
525  distMap.schedule(),
526  distMap.constructSize(),
527  distMap.subMap(),
528  false,
529  distMap.constructMap(),
530  false,
531  allVectors_[vectorI],
532  flipOp()
533  );
534  allVectors_[vectorI].setCapacity(allVectors_[vectorI].size());
535  }
536  }
537 
538 
539  label n = 0;
540  forAll(allTracks_, trackI)
541  {
542  n += allTracks_[trackI].size();
543  }
544 
545  Info<< " Tracks:" << allTracks_.size() << nl
546  << " Total samples:" << n
547  << endl;
548 
549 
550  // Massage into form suitable for writers
551  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
552 
553  if (Pstream::master() && allTracks_.size())
554  {
555  // Make output directory
556 
557  fileName vtkPath
558  (
560  ? runTime.path()/".."/"postProcessing"/"sets"/name()
561  : runTime.path()/"postProcessing"/"sets"/name()
562  );
563  if (mesh_.name() != fvMesh::defaultRegion)
564  {
565  vtkPath = vtkPath/mesh_.name();
566  }
567  vtkPath = vtkPath/mesh_.time().timeName();
568  vtkPath.clean();
569  mkDir(vtkPath);
570 
571  // Convert track positions
572 
573  PtrList<coordSet> tracks(allTracks_.size());
574  forAll(allTracks_, trackI)
575  {
576  tracks.set
577  (
578  trackI,
579  new coordSet
580  (
581  "track" + Foam::name(trackI),
582  sampledSetAxis_ //"xyz"
583  )
584  );
585  tracks[trackI].transfer(allTracks_[trackI]);
586  }
587 
588  // Convert scalar values
589 
590  if (allScalars_.size() > 0)
591  {
592  List<List<scalarField>> scalarValues(allScalars_.size());
593 
594  forAll(allScalars_, scalarI)
595  {
596  DynamicList<scalarList>& allTrackVals =
597  allScalars_[scalarI];
598  scalarValues[scalarI].setSize(allTrackVals.size());
599 
600  forAll(allTrackVals, trackI)
601  {
602  scalarList& trackVals = allTrackVals[trackI];
603  scalarValues[scalarI][trackI].transfer(trackVals);
604  }
605  }
606 
607  fileName vtkFile
608  (
609  vtkPath
610  / scalarFormatterPtr_().getFileName
611  (
612  tracks[0],
613  scalarNames_
614  )
615  );
616 
617  Info<< " Writing data to " << vtkFile.path() << endl;
618 
619  scalarFormatterPtr_().write
620  (
621  true, // writeTracks
622  tracks,
623  scalarNames_,
624  scalarValues,
625  OFstream(vtkFile)()
626  );
627  }
628 
629  // Convert vector values
630 
631  if (allVectors_.size() > 0)
632  {
633  List<List<vectorField>> vectorValues(allVectors_.size());
634 
635  forAll(allVectors_, vectorI)
636  {
637  DynamicList<vectorList>& allTrackVals =
638  allVectors_[vectorI];
639  vectorValues[vectorI].setSize(allTrackVals.size());
640 
641  forAll(allTrackVals, trackI)
642  {
643  vectorList& trackVals = allTrackVals[trackI];
644  vectorValues[vectorI][trackI].transfer(trackVals);
645  }
646  }
647 
648  fileName vtkFile
649  (
650  vtkPath
651  / vectorFormatterPtr_().getFileName
652  (
653  tracks[0],
654  vectorNames_
655  )
656  );
657 
658  vectorFormatterPtr_().write
659  (
660  true, // writeTracks
661  tracks,
662  vectorNames_,
663  vectorValues,
664  OFstream(vtkFile)()
665  );
666  }
667  }
668 
669  return true;
670 }
671 
672 
674 {
675  if (&mpm.mesh() == &mesh_)
676  {
677  read(dict_);
678  }
679 }
680 
681 
683 {
684  if (&mesh == &mesh_)
685  {
686  // Moving mesh affects the search tree
687  read(dict_);
688  }
689 }
690 
691 
692 // ************************************************************************* //
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:667
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
fileName path() const
Return path.
Definition: Time.H:265
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:158
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:319
static const char *const typeName
Definition: Field.H:105
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
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:151
Field reading functions for post-processing utilities.
patches[0]
virtual ~streamLine()
Destructor.
Definition: streamLine.C:319
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.
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
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
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamLine.C:325
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
A class for handling words, derived from string.
Definition: word.H:59
Particle class that samples fields as it passes through. Used in streamline calculation.
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:183
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:35
static const char nl
Definition: Ostream.H:260
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
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: streamLine.C:673
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A Cloud of streamLine particles.
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:331
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:303
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.
static const NamedEnum< trackDirection, 3 > trackDirectionNames_
Track direction enumeration names.
Definition: streamLine.H:208
messageStream Info
label n
streamLine(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: streamLine.C:303
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:253
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
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...
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: streamLine.C:682
virtual bool execute()
Do nothing.
Definition: streamLine.C:426
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void transfer(List< T > &)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:285
const labelList & cells() const
Definition: sampledSet.H:209
A List with indirect addressing.
Definition: IndirectList.H:101
virtual bool write()
Calculate and write the streamlines.
Definition: streamLine.C:432
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
IOerror FatalIOError