streamLine.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Pstream.H"
27 #include "functionObjectList.H"
28 #include "streamLine.H"
29 #include "fvMesh.H"
31 #include "ReadFields.H"
32 #include "meshSearch.H"
33 #include "sampledSet.H"
34 #include "globalIndex.H"
35 #include "mapDistribute.H"
36 #include "interpolationCellPoint.H"
37 #include "PatchTools.H"
38 #include "mapPolyMesh.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace functionObjects
46 {
47  defineTypeNameAndDebug(streamLine, 0);
48  addToRunTimeSelectionTable(functionObject, streamLine, dictionary);
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
56 Foam::functionObjects::streamLine::wallPatch() const
57 {
58  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
59 
60  const polyBoundaryMesh& patches = mesh.boundaryMesh();
61 
62  label nFaces = 0;
63 
64  forAll(patches, patchi)
65  {
66  if (isA<wallPolyPatch>(patches[patchi]))
67  {
68  nFaces += patches[patchi].size();
69  }
70  }
71 
72  labelList addressing(nFaces);
73 
74  nFaces = 0;
75 
76  forAll(patches, patchi)
77  {
78  if (isA<wallPolyPatch>(patches[patchi]))
79  {
80  const polyPatch& pp = patches[patchi];
81 
82  forAll(pp, i)
83  {
84  addressing[nFaces++] = pp.start()+i;
85  }
86  }
87  }
88 
89  return autoPtr<indirectPrimitivePatch>
90  (
92  (
93  IndirectList<face>
94  (
95  mesh.faces(),
96  addressing
97  ),
98  mesh.points()
99  )
100  );
101 }
102 
103 
104 void Foam::functionObjects::streamLine::track()
105 {
106  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
107 
108  IDLList<streamLineParticle> initialParticles;
109  streamLineParticleCloud particles
110  (
111  mesh,
112  cloudName_,
113  initialParticles
114  );
115 
116  const sampledSet& seedPoints = sampledSetPtr_();
117 
118  forAll(seedPoints, i)
119  {
120  particles.addParticle
121  (
122  new streamLineParticle
123  (
124  mesh,
125  seedPoints[i],
126  seedPoints.cells()[i],
127  lifeTime_
128  )
129  );
130  }
131 
132  label nSeeds = returnReduce(particles.size(), sumOp<label>());
133 
134  Info << " seeded " << nSeeds << " particles" << endl;
135 
136  // Read or lookup fields
137  PtrList<volScalarField> vsFlds;
138  PtrList<interpolation<scalar>> vsInterp;
139  PtrList<volVectorField> vvFlds;
140  PtrList<interpolation<vector>> vvInterp;
141 
142  label UIndex = -1;
143 
144  label nScalar = 0;
145  label nVector = 0;
146 
147  forAll(fields_, i)
148  {
149  if (mesh.foundObject<volScalarField>(fields_[i]))
150  {
151  nScalar++;
152  }
153  else if (mesh.foundObject<volVectorField>(fields_[i]))
154  {
155  nVector++;
156  }
157  else
158  {
160  << "Cannot find field " << fields_[i] << nl
161  << "Valid scalar fields are:"
162  << mesh.names(volScalarField::typeName) << nl
163  << "Valid vector fields are:"
164  << mesh.names(volVectorField::typeName)
165  << exit(FatalError);
166  }
167  }
168  vsInterp.setSize(nScalar);
169  nScalar = 0;
170  vvInterp.setSize(nVector);
171  nVector = 0;
172 
173  forAll(fields_, i)
174  {
175  if (mesh.foundObject<volScalarField>(fields_[i]))
176  {
177  const volScalarField& f = mesh.lookupObject<volScalarField>
178  (
179  fields_[i]
180  );
181  vsInterp.set
182  (
183  nScalar++,
185  (
186  interpolationScheme_,
187  f
188  )
189  );
190  }
191  else if (mesh.foundObject<volVectorField>(fields_[i]))
192  {
193  const volVectorField& f = mesh.lookupObject<volVectorField>
194  (
195  fields_[i]
196  );
197 
198  if (f.name() == UName_)
199  {
200  UIndex = nVector;
201  }
202 
203  vvInterp.set
204  (
205  nVector++,
207  (
208  interpolationScheme_,
209  f
210  )
211  );
212  }
213  }
214 
215  // Store the names
216  scalarNames_.setSize(vsInterp.size());
217  forAll(vsInterp, i)
218  {
219  scalarNames_[i] = vsInterp[i].psi().name();
220  }
221  vectorNames_.setSize(vvInterp.size());
222  forAll(vvInterp, i)
223  {
224  vectorNames_[i] = vvInterp[i].psi().name();
225  }
226 
227  // Check that we know the index of U in the interpolators.
228 
229  if (UIndex == -1)
230  {
232  << "Cannot find field to move particles with : " << UName_ << nl
233  << "This field has to be present in the sampled fields " << fields_
234  << " and in the objectRegistry."
235  << exit(FatalError);
236  }
237 
238  // Sampled data
239  // ~~~~~~~~~~~~
240 
241  // Size to maximum expected sizes.
242  allTracks_.clear();
243  allTracks_.setCapacity(nSeeds);
244  allScalars_.setSize(vsInterp.size());
245  forAll(allScalars_, i)
246  {
247  allScalars_[i].clear();
248  allScalars_[i].setCapacity(nSeeds);
249  }
250  allVectors_.setSize(vvInterp.size());
251  forAll(allVectors_, i)
252  {
253  allVectors_[i].clear();
254  allVectors_[i].setCapacity(nSeeds);
255  }
256 
257 
258  // Additional particle info
259  streamLineParticle::trackingData td
260  (
261  particles,
262  vsInterp,
263  vvInterp,
264  UIndex, // index of U in vvInterp
265  trackForward_, // track in +u direction?
266  nSubCycle_, // automatic track control:step through cells in steps?
267  trackLength_, // fixed track length
268 
269  allTracks_,
270  allScalars_,
271  allVectors_
272  );
273 
274 
275  // Set very large dt. Note: cannot use GREAT since 1/GREAT is SMALL
276  // which is a trigger value for the tracking...
277  const scalar trackTime = Foam::sqrt(GREAT);
278 
279  // Track
280  particles.move(td, trackTime);
281 }
282 
283 
284 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
285 
286 Foam::functionObjects::streamLine::streamLine
287 (
288  const word& name,
289  const Time& runTime,
290  const dictionary& dict
291 )
292 :
293  functionObject(name),
294  obr_
295  (
297  (
299  )
300  ),
301  dict_(dict),
302  nSubCycle_(0)
303 {
304  if (!isA<fvMesh>(obr_))
305  {
307  << "objectRegistry is not an fvMesh" << exit(FatalError);
308  }
309 
310  read(dict_);
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
315 
317 {}
318 
319 
320 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
321 
323 {
324  Info<< type() << " " << name() << ":" << nl;
325 
326  dict.lookup("fields") >> fields_;
327  if (dict.found("U"))
328  {
329  dict.lookup("U") >> UName_;
330  }
331  else
332  {
333  UName_ = "U";
334  if (dict.found("U"))
335  {
336  IOWarningInFunction(dict)
337  << "Using deprecated entry \"U\"."
338  << " Please use \"UName\" instead."
339  << endl;
340  dict.lookup("U") >> UName_;
341  }
342  }
343 
344  if (findIndex(fields_, UName_) == -1)
345  {
347  << "Velocity field for tracking " << UName_
348  << " should be present in the list of fields " << fields_
349  << exit(FatalIOError);
350  }
351 
352 
353  dict.lookup("trackForward") >> trackForward_;
354  dict.lookup("lifeTime") >> lifeTime_;
355  if (lifeTime_ < 1)
356  {
358  << "Illegal value " << lifeTime_ << " for lifeTime"
359  << exit(FatalError);
360  }
361 
362 
363  bool subCycling = dict.found("nSubCycle");
364  bool fixedLength = dict.found("trackLength");
365 
366  if (subCycling && fixedLength)
367  {
369  << "Cannot both specify automatic time stepping (through '"
370  << "nSubCycle' specification) and fixed track length (through '"
371  << "trackLength')"
372  << exit(FatalIOError);
373  }
374 
375 
376  nSubCycle_ = 1;
377  if (dict.readIfPresent("nSubCycle", nSubCycle_))
378  {
379  trackLength_ = VGREAT;
380  if (nSubCycle_ < 1)
381  {
382  nSubCycle_ = 1;
383  }
384  Info<< " automatic track length specified through"
385  << " number of sub cycles : " << nSubCycle_ << nl << endl;
386  }
387  else
388  {
389  dict.lookup("trackLength") >> trackLength_;
390 
391  Info<< " fixed track length specified : "
392  << trackLength_ << nl << endl;
393  }
394 
395 
396  interpolationScheme_ = dict.lookupOrDefault
397  (
398  "interpolationScheme",
400  );
401 
402  cloudName_ = dict.lookupOrDefault<word>("cloudName", "streamLine");
403  dict.lookup("seedSampleSet") >> seedSet_;
404 
405  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
406 
407  meshSearchPtr_.reset(new meshSearch(mesh));
408 
409  const dictionary& coeffsDict = dict.subDict(seedSet_ + "Coeffs");
410  sampledSetPtr_ = sampledSet::New
411  (
412  seedSet_,
413  mesh,
414  meshSearchPtr_(),
415  coeffsDict
416  );
417  coeffsDict.lookup("axis") >> sampledSetAxis_;
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  const fvMesh& mesh = dynamic_cast<const fvMesh&>(obr_);
438 
439 
440  // Do all injection and tracking
441  track();
442 
443 
444  if (Pstream::parRun())
445  {
446  // Append slave tracks to master ones
447  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
448 
449  globalIndex globalTrackIDs(allTracks_.size());
450 
451  // Construct a distribution map to pull all to the master.
452  labelListList sendMap(Pstream::nProcs());
453  labelListList recvMap(Pstream::nProcs());
454 
455  if (Pstream::master())
456  {
457  // Master: receive all. My own first, then consecutive
458  // processors.
459  label trackI = 0;
460 
461  forAll(recvMap, proci)
462  {
463  labelList& fromProc = recvMap[proci];
464  fromProc.setSize(globalTrackIDs.localSize(proci));
465  forAll(fromProc, i)
466  {
467  fromProc[i] = trackI++;
468  }
469  }
470  }
471 
472  labelList& toMaster = sendMap[0];
473  toMaster.setSize(globalTrackIDs.localSize());
474  forAll(toMaster, i)
475  {
476  toMaster[i] = i;
477  }
478 
479  const mapDistribute distMap
480  (
481  globalTrackIDs.size(),
482  sendMap.xfer(),
483  recvMap.xfer()
484  );
485 
486 
487  // Distribute the track positions. Note: use scheduled comms
488  // to prevent buffering.
490  (
492  distMap.schedule(),
493  distMap.constructSize(),
494  distMap.subMap(),
495  false,
496  distMap.constructMap(),
497  false,
498  allTracks_,
499  flipOp()
500  );
501 
502  // Distribute the scalars
503  forAll(allScalars_, scalarI)
504  {
505  allScalars_[scalarI].shrink();
507  (
509  distMap.schedule(),
510  distMap.constructSize(),
511  distMap.subMap(),
512  false,
513  distMap.constructMap(),
514  false,
515  allScalars_[scalarI],
516  flipOp()
517  );
518  allScalars_[scalarI].setCapacity(allScalars_[scalarI].size());
519  }
520  // Distribute the vectors
521  forAll(allVectors_, vectorI)
522  {
523  allVectors_[vectorI].shrink();
525  (
527  distMap.schedule(),
528  distMap.constructSize(),
529  distMap.subMap(),
530  false,
531  distMap.constructMap(),
532  false,
533  allVectors_[vectorI],
534  flipOp()
535  );
536  allVectors_[vectorI].setCapacity(allVectors_[vectorI].size());
537  }
538  }
539 
540 
541  label n = 0;
542  forAll(allTracks_, trackI)
543  {
544  n += allTracks_[trackI].size();
545  }
546 
547  Info<< " Tracks:" << allTracks_.size() << nl
548  << " Total samples:" << n
549  << endl;
550 
551 
552  // Massage into form suitable for writers
553  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
554 
555  if (Pstream::master() && allTracks_.size())
556  {
557  // Make output directory
558 
559  fileName vtkPath
560  (
562  ? runTime.path()/".."/"postProcessing"/"sets"/name()
563  : runTime.path()/"postProcessing"/"sets"/name()
564  );
565  if (mesh.name() != fvMesh::defaultRegion)
566  {
567  vtkPath = vtkPath/mesh.name();
568  }
569  vtkPath = vtkPath/mesh.time().timeName();
570 
571  mkDir(vtkPath);
572 
573  // Convert track positions
574 
575  PtrList<coordSet> tracks(allTracks_.size());
576  forAll(allTracks_, trackI)
577  {
578  tracks.set
579  (
580  trackI,
581  new coordSet
582  (
583  "track" + Foam::name(trackI),
584  sampledSetAxis_ //"xyz"
585  )
586  );
587  tracks[trackI].transfer(allTracks_[trackI]);
588  }
589 
590  // Convert scalar values
591 
592  if (allScalars_.size() > 0)
593  {
594  List<List<scalarField>> scalarValues(allScalars_.size());
595 
596  forAll(allScalars_, scalarI)
597  {
598  DynamicList<scalarList>& allTrackVals =
599  allScalars_[scalarI];
600  scalarValues[scalarI].setSize(allTrackVals.size());
601 
602  forAll(allTrackVals, trackI)
603  {
604  scalarList& trackVals = allTrackVals[trackI];
605  scalarValues[scalarI][trackI].transfer(trackVals);
606  }
607  }
608 
609  fileName vtkFile
610  (
611  vtkPath
612  / scalarFormatterPtr_().getFileName
613  (
614  tracks[0],
615  scalarNames_
616  )
617  );
618 
619  Info<< " Writing data to " << vtkFile.path() << endl;
620 
621  scalarFormatterPtr_().write
622  (
623  true, // writeTracks
624  tracks,
625  scalarNames_,
626  scalarValues,
627  OFstream(vtkFile)()
628  );
629  }
630 
631  // Convert vector values
632 
633  if (allVectors_.size() > 0)
634  {
635  List<List<vectorField>> vectorValues(allVectors_.size());
636 
637  forAll(allVectors_, vectorI)
638  {
639  DynamicList<vectorList>& allTrackVals =
640  allVectors_[vectorI];
641  vectorValues[vectorI].setSize(allTrackVals.size());
642 
643  forAll(allTrackVals, trackI)
644  {
645  vectorList& trackVals = allTrackVals[trackI];
646  vectorValues[vectorI][trackI].transfer(trackVals);
647  }
648  }
649 
650  fileName vtkFile
651  (
652  vtkPath
653  / vectorFormatterPtr_().getFileName
654  (
655  tracks[0],
656  vectorNames_
657  )
658  );
659 
660  vectorFormatterPtr_().write
661  (
662  true, // writeTracks
663  tracks,
664  vectorNames_,
665  vectorValues,
666  OFstream(vtkFile)()
667  );
668  }
669  }
670 
671  return true;
672 }
673 
674 
676 {
677  const fvMesh& mesh_ = dynamic_cast<const fvMesh&>(obr_);
678 
679  if (&mpm.mesh() == &mesh_)
680  {
681  read(dict_);
682  }
683 }
684 
685 
687 {
688  const fvMesh& mesh_ = dynamic_cast<const fvMesh&>(obr_);
689 
690  if (&mesh == &mesh_)
691  {
692  // Moving mesh affects the search tree
693  read(dict_);
694  }
695 }
696 
697 
698 // ************************************************************************* //
const Time & time() const
Return time.
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A class for handling file names.
Definition: fileName.H:69
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:137
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:442
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:90
static const char *const typeName
Definition: Field.H:94
Output to file stream.
Definition: OFstream.H:81
static autoPtr< interpolation< scalar > > New(const word &interpolationType, const GeometricField< scalar, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
virtual ~streamLine()
Destructor.
Definition: streamLine.C:316
Abstract base-class for Time/database function objects.
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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:163
virtual bool read(const dictionary &)
Read the field average data.
Definition: streamLine.C:322
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
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
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 Pstream::scheduled.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
List< label > labelList
A List of labels.
Definition: labelList.H:56
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
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:262
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence 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:295
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: streamLine.C:675
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setSize(const label)
Reset size of List.
Definition: List.C:295
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:393
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:399
label patchi
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Class containing processor-to-processor mapping information.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
addToRunTimeSelectionTable(functionObject, blendingFactor, dictionary)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
label n
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: streamLine.C:686
virtual bool execute()
Do nothing.
Definition: streamLine.C:426
Registry of regIOobjects.
fileName path() const
Return path.
Definition: Time.H:269
void transfer(List< T > &)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:259
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:249
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
virtual bool write()
Calculate and write the steamlines.
Definition: streamLine.C:432
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:357
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
IOerror FatalIOError