steadyParticleTracks.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  steadyParticleTracks
26 
27 Description
28  Generates a VTK file of particle tracks for cases that were computed using
29  a steady-state cloud
30  NOTE: case must be re-constructed (if running in parallel) before use
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "argList.H"
35 #include "Cloud.H"
36 #include "IOdictionary.H"
37 #include "fvMesh.H"
38 #include "Time.H"
39 #include "timeSelector.H"
40 #include "OFstream.H"
41 #include "passiveParticleCloud.H"
42 
43 #include "SortableList.H"
44 #include "IOobjectList.H"
45 #include "PtrList.H"
46 #include "Field.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 using namespace Foam;
52 
53 namespace Foam
54 {
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 label validateFields
59 (
60  const List<word>& userFields,
61  const IOobjectList& cloudObjs
62 )
63 {
64  List<bool> ok(userFields.size(), false);
65 
66  forAll(userFields, i)
67  {
68  ok[i] = ok[i] || fieldOk<label>(cloudObjs, userFields[i]);
69  ok[i] = ok[i] || fieldOk<scalar>(cloudObjs, userFields[i]);
70  ok[i] = ok[i] || fieldOk<vector>(cloudObjs, userFields[i]);
71  ok[i] = ok[i] || fieldOk<sphericalTensor>(cloudObjs, userFields[i]);
72  ok[i] = ok[i] || fieldOk<symmTensor>(cloudObjs, userFields[i]);
73  ok[i] = ok[i] || fieldOk<tensor>(cloudObjs, userFields[i]);
74  }
75 
76  label nOk = 0;
77  forAll(ok, i)
78  {
79  if (ok[i])
80  {
81  nOk++;
82  }
83  else
84  {
85  Info << "\n*** Warning: user specified field '" << userFields[i]
86  << "' unavailable" << endl;
87  }
88  }
89 
90  return nOk;
91 }
92 
93 
94 template<>
95 void writeVTK(OFstream& os, const label& value)
96 {
97  os << value;
98 }
99 
100 
101 template<>
102 void writeVTK(OFstream& os, const scalar& value)
103 {
104  os << value;
105 }
106 
107 }
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 int main(int argc, char *argv[])
112 {
115  #include "addRegionOption.H"
116  #include "addDictOption.H"
117 
118  #include "setRootCase.H"
119 
120  #include "createTime.H"
121  instantList timeDirs = timeSelector::select0(runTime, args);
122  #include "createNamedMesh.H"
123  #include "createFields.H"
124 
125  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 
127  fileName vtkPath(runTime.path()/"VTK");
128  mkDir(vtkPath);
129 
130  typedef HashTable<label, labelPair, labelPair::Hash<>> trackTableType;
131 
132  forAll(timeDirs, timeI)
133  {
134  runTime.setTime(timeDirs[timeI], timeI);
135  Info<< "Time = " << runTime.timeName() << endl;
136 
137  fileName vtkTimePath(runTime.path()/"VTK"/runTime.timeName());
138  mkDir(vtkTimePath);
139 
140  Info<< " Reading particle positions" << endl;
141 
142  PtrList<passiveParticle> particles(0);
143 
144  // Transfer particles to (more convenient) list
145  {
147  Info<< "\n Read " << returnReduce(ppc.size(), sumOp<label>())
148  << " particles" << endl;
149 
150  particles.setSize(ppc.size());
151 
152  label i = 0;
153  forAllIter(passiveParticleCloud, ppc, iter)
154  {
155  particles.set(i++, ppc.remove(&iter()));
156  }
157 
158  // myCloud should now be empty
159  }
160 
161  List<label> particleToTrack(particles.size());
162  label nTracks = 0;
163 
164  {
165  trackTableType trackTable;
166  forAll(particles, i)
167  {
168  const label origProc = particles[i].origProc();
169  const label origId = particles[i].origId();
170 
171  const trackTableType::const_iterator& iter =
172  trackTable.find(labelPair(origProc, origId));
173 
174  if (iter == trackTable.end())
175  {
176  particleToTrack[i] = nTracks;
177  trackTable.insert(labelPair(origProc, origId), nTracks);
178  nTracks++;
179  }
180  else
181  {
182  particleToTrack[i] = iter();
183  }
184  }
185  }
186 
187 
188  if (nTracks == 0)
189  {
190  Info<< "\n No track data" << endl;
191  }
192  else
193  {
194  Info<< "\n Generating " << nTracks << " tracks" << endl;
195 
196  // Determine length of each track
197  labelList trackLengths(nTracks, 0);
198  forAll(particleToTrack, i)
199  {
200  const label trackI = particleToTrack[i];
201  trackLengths[trackI]++;
202  }
203 
204  // Particle "age" property used to sort the tracks
205  List<SortableList<scalar>> agePerTrack(nTracks);
206  List<List<label>> particleMap(nTracks);
207 
208  forAll(trackLengths, i)
209  {
210  const label length = trackLengths[i];
211  agePerTrack[i].setSize(length);
212  particleMap[i].setSize(length);
213  }
214 
215  // Store the particle age per track
216  IOobjectList cloudObjs
217  (
218  mesh,
219  runTime.timeName(),
221  );
222 
223  // TODO: gather age across all procs
224  {
225  tmp<scalarField> tage =
226  readParticleField<scalar>("age", cloudObjs);
227 
228  const scalarField& age = tage();
229 
230  List<label> trackSamples(nTracks, 0);
231 
232  forAll(particleToTrack, i)
233  {
234  const label trackI = particleToTrack[i];
235  const label sampleI = trackSamples[trackI];
236  agePerTrack[trackI][sampleI] = age[i];
237  particleMap[trackI][sampleI] = i;
238  trackSamples[trackI]++;
239  }
240  tage.clear();
241  }
242 
243 
244  if (Pstream::master())
245  {
246  OFstream os(vtkTimePath/"particleTracks.vtk");
247 
248  Info<< "\n Writing particle tracks to " << os.name() << endl;
249 
250  label nPoints = sum(trackLengths);
251 
252  os << "# vtk DataFile Version 2.0" << nl
253  << "particleTracks" << nl
254  << "ASCII" << nl
255  << "DATASET POLYDATA" << nl
256  << "POINTS " << nPoints << " float" << nl;
257 
258  Info<< "\n Writing points" << endl;
259 
260  {
261  forAll(agePerTrack, i)
262  {
263  agePerTrack[i].sort();
264 
265  const labelList& ids = agePerTrack[i].indices();
266  labelList& particleIds = particleMap[i];
267 
268  {
269  // Update addressing
270  List<label> sortedIds(ids);
271  forAll(sortedIds, j)
272  {
273  sortedIds[j] = particleIds[ids[j]];
274  }
275  particleIds = sortedIds;
276  }
277 
278  forAll(ids, j)
279  {
280  const label localId = particleIds[j];
281  const vector& pos = particles[localId].position();
282  os << pos.x() << ' ' << pos.y() << ' ' << pos.z()
283  << nl;
284  }
285  }
286  }
287 
288 
289  // Write track (line) connectivity to file
290 
291  Info<< "\n Writing track lines" << endl;
292  os << "\nLINES " << nTracks << ' ' << nPoints + nTracks << nl;
293 
294  // Write ids of track points to file
295  {
296  label globalPtI = 0;
297  forAll(particleMap, i)
298  {
299  os << particleMap[i].size() << nl;
300 
301  forAll(particleMap[i], j)
302  {
303  os << ' ' << globalPtI++;
304 
305  if (((j + 1) % 10 == 0) && (j != 0))
306  {
307  os << nl;
308  }
309  }
310 
311  os << nl;
312  }
313  }
314 
315 
316  const label nFields = validateFields(userFields, cloudObjs);
317 
318  os << "POINT_DATA " << nPoints << nl
319  << "FIELD attributes " << nFields << nl;
320 
321  Info<< "\n Processing fields" << nl << endl;
322 
323  processFields<label>(os, particleMap, userFields, cloudObjs);
324  processFields<scalar>(os, particleMap, userFields, cloudObjs);
325  processFields<vector>(os, particleMap, userFields, cloudObjs);
326  processFields<sphericalTensor>
327  (os, particleMap, userFields, cloudObjs);
328  processFields<symmTensor>
329  (os, particleMap, userFields, cloudObjs);
330  processFields<tensor>(os, particleMap, userFields, cloudObjs);
331 
332  }
333  }
334  Info<< endl;
335  }
336 
337  Info<< "\ndone" << endl;
338 
339  return 0;
340 }
341 
342 
343 // ************************************************************************* //
#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
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
Output to file stream.
Definition: OFstream.H:82
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
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
static void noParallel()
Remove the parallel options.
Definition: argList.C:148
const Cmpt & z() const
Definition: VectorI.H:87
A Cloud of passive particles.
const Cmpt & y() const
Definition: VectorI.H:81
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dynamicFvMesh & mesh
List< word > userFields(propsDict.lookup("fields"))
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
An STL-conforming hash table.
Definition: HashTable.H:62
const Cmpt & x() const
Definition: VectorI.H:75
static const char nl
Definition: Ostream.H:262
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:257
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:297
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
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const word cloudName(propsDict.lookup("cloudName"))
void writeVTK(OFstream &os, const Type &value)
A class for managing temporary objects.
Definition: PtrList.H:53
Foam::argList args(argc, argv)
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.