particleTracks.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 Application
25  particleTracks
26 
27 Description
28  Generates a VTK file of particle tracks for cases that were computed using
29  a tracked-parcel-type cloud.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "Cloud.H"
35 #include "IOdictionary.H"
36 #include "fvMesh.H"
37 #include "Time.H"
38 #include "timeSelector.H"
39 #include "OFstream.H"
40 #include "passiveParticleCloud.H"
41 #include "writer.H"
42 
43 using namespace Foam;
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
50  #include "addRegionOption.H"
51 
52  #include "setRootCase.H"
53 
54  #include "createTime.H"
55  instantList timeDirs = timeSelector::select0(runTime, args);
56  #include "createNamedMesh.H"
57  #include "createFields.H"
58 
59  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61  fileName vtkPath(runTime.path()/"VTK");
62  mkDir(vtkPath);
63 
64  Info<< "Scanning times to determine track data for cloud " << cloudName
65  << nl << endl;
66 
67  labelList maxIds(Pstream::nProcs(), -1);
68  forAll(timeDirs, timeI)
69  {
70  runTime.setTime(timeDirs[timeI], timeI);
71  Info<< "Time = " << runTime.timeName() << endl;
72 
73  Info<< " Reading particle positions" << endl;
75  Info<< " Read " << returnReduce(myCloud.size(), sumOp<label>())
76  << " particles" << endl;
77 
79  {
80  label origId = iter().origId();
81  label origProc = iter().origProc();
82 
83  if (origProc >= maxIds.size())
84  {
85  // Expand size
86  maxIds.setSize(origProc+1, -1);
87  }
88 
89  maxIds[origProc] = max(maxIds[origProc], origId);
90  }
91  }
92 
93  label maxNProcs = returnReduce(maxIds.size(), maxOp<label>());
94 
95  Info<< "Detected particles originating from " << maxNProcs
96  << " processors." << nl << endl;
97 
98  maxIds.setSize(maxNProcs, -1);
99 
102 
103  labelList numIds = maxIds + 1;
104 
105  Info<< nl << "Particle statistics:" << endl;
106  forAll(maxIds, proci)
107  {
108  Info<< " Found " << numIds[proci] << " particles originating"
109  << " from processor " << proci << endl;
110  }
111  Info<< nl << endl;
112 
113 
114  // calc starting ids for particles on each processor
115  List<label> startIds(numIds.size(), 0);
116  for (label i = 0; i < numIds.size()-1; i++)
117  {
118  startIds[i+1] += startIds[i] + numIds[i];
119  }
120  label nParticle = startIds.last() + numIds[startIds.size()-1];
121 
122 
123 
124  // number of tracks to generate
125  label nTracks = nParticle/sampleFrequency;
126 
127  // storage for all particle tracks
128  List<DynamicList<vector>> allTracks(nTracks);
129 
130  Info<< "\nGenerating " << nTracks << " particle tracks for cloud "
131  << cloudName << nl << endl;
132 
133  forAll(timeDirs, timeI)
134  {
135  runTime.setTime(timeDirs[timeI], timeI);
136  Info<< "Time = " << runTime.timeName() << endl;
137 
138  List<pointField> allPositions(Pstream::nProcs());
139  List<labelField> allOrigIds(Pstream::nProcs());
140  List<labelField> allOrigProcs(Pstream::nProcs());
141 
142  // Read particles. Will be size 0 if no particles.
143  Info<< " Reading particle positions" << endl;
145 
146  // collect the track data on all processors that have positions
147  allPositions[Pstream::myProcNo()].setSize
148  (
149  myCloud.size(),
151  );
152  allOrigIds[Pstream::myProcNo()].setSize(myCloud.size(), 0);
153  allOrigProcs[Pstream::myProcNo()].setSize(myCloud.size(), 0);
154 
155  label i = 0;
156  forAllConstIter(passiveParticleCloud, myCloud, iter)
157  {
158  allPositions[Pstream::myProcNo()][i] = iter().position();
159  allOrigIds[Pstream::myProcNo()][i] = iter().origId();
160  allOrigProcs[Pstream::myProcNo()][i] = iter().origProc();
161  i++;
162  }
163 
164  // collect the track data on the master processor
165  Pstream::gatherList(allPositions);
166  Pstream::gatherList(allOrigIds);
167  Pstream::gatherList(allOrigProcs);
168 
169  Info<< " Constructing tracks" << nl << endl;
170  if (Pstream::master())
171  {
172  forAll(allPositions, proci)
173  {
174  forAll(allPositions[proci], i)
175  {
176  label globalId =
177  startIds[allOrigProcs[proci][i]]
178  + allOrigIds[proci][i];
179 
180  if (globalId % sampleFrequency == 0)
181  {
182  label trackId = globalId/sampleFrequency;
183  if (allTracks[trackId].size() < maxPositions)
184  {
185  allTracks[trackId].append
186  (
187  allPositions[proci][i]
188  );
189  }
190  }
191  }
192  }
193  }
194  }
195 
196 
197  if (Pstream::master())
198  {
199  PtrList<coordSet> tracks(allTracks.size());
200  forAll(allTracks, trackI)
201  {
202  tracks.set
203  (
204  trackI,
205  new coordSet
206  (
207  "track" + Foam::name(trackI),
208  "distance"
209  )
210  );
211  tracks[trackI].transfer(allTracks[trackI]);
212  }
213 
214  autoPtr<writer<scalar>> scalarFormatterPtr = writer<scalar>::New
215  (
216  setFormat
217  );
218 
219  //OFstream vtkTracks(vtkPath/"particleTracks.vtk");
220  fileName vtkFile
221  (
222  scalarFormatterPtr().getFileName
223  (
224  tracks[0],
225  wordList(0)
226  )
227  );
228 
229  OFstream vtkTracks
230  (
231  vtkPath
232  / "particleTracks." + vtkFile.ext()
233  );
234 
235  Info<< "\nWriting particle tracks in " << setFormat
236  << " format to " << vtkTracks.name()
237  << nl << endl;
238 
239 
240  scalarFormatterPtr().write
241  (
242  true, // writeTracks
243  tracks,
244  wordList(0),
246  vtkTracks
247  );
248  }
249 
250  return 0;
251 }
252 
253 
254 // ************************************************************************* //
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Output to file stream.
Definition: OFstream.H:81
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:417
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
const word cloudName(propsDict.lookup("cloudName"))
word setFormat(propsDict.lookupOrDefault< word >("setFormat","vtk"))
label sampleFrequency(readLabel(propsDict.lookup("sampleFrequency")))
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
A Cloud of passive particles.
dynamicFvMesh & mesh
Holds list of sampling positions.
Definition: coordSet.H:49
label maxPositions(readLabel(propsDict.lookup("maxPositions")))
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:253
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:295
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< word > wordList
A List of words.
Definition: fileName.H:54
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:399
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
messageStream Info
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:53
T & last()
Return the last element of the list.
Definition: UListI.H:128
Foam::argList args(argc, argv)
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Namespace for OpenFOAM.