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