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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  steadyParticleTracks
28 
29 Group
30  grpPostProcessingUtilities
31 
32 Description
33  Generate a legacy VTK file of particle tracks for cases that were
34  computed using a steady-state cloud
35 
36  Note:
37  - Case must be re-constructed (if running in parallel) before use
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "argList.H"
42 #include "Cloud.H"
43 #include "IOdictionary.H"
44 #include "fvMesh.H"
45 #include "Time.H"
46 #include "timeSelector.H"
47 #include "OFstream.H"
48 #include "passiveParticleCloud.H"
49 #include "labelPairHashes.H"
50 #include "SortableList.H"
51 #include "IOobjectList.H"
52 #include "PtrList.H"
53 #include "Field.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 using namespace Foam;
59 
60 namespace Foam
61 {
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 label validateFields
66 (
67  const List<word>& userFields,
68  const IOobjectList& cloudObjs
69 )
70 {
71  List<bool> ok(userFields.size(), false);
72 
73  forAll(userFields, i)
74  {
75  ok[i] = ok[i] || fieldOk<label>(cloudObjs, userFields[i]);
76  ok[i] = ok[i] || fieldOk<scalar>(cloudObjs, userFields[i]);
77  ok[i] = ok[i] || fieldOk<vector>(cloudObjs, userFields[i]);
78  ok[i] = ok[i] || fieldOk<sphericalTensor>(cloudObjs, userFields[i]);
79  ok[i] = ok[i] || fieldOk<symmTensor>(cloudObjs, userFields[i]);
80  ok[i] = ok[i] || fieldOk<tensor>(cloudObjs, userFields[i]);
81  }
82 
83  label nOk = 0;
84  forAll(ok, i)
85  {
86  if (ok[i])
87  {
88  ++nOk;
89  }
90  else
91  {
92  Info << "\n*** Warning: user specified field '" << userFields[i]
93  << "' unavailable" << endl;
94  }
95  }
96 
97  return nOk;
98 }
99 
100 
101 template<>
102 void writeVTK(OFstream& os, const label& value)
103 {
104  os << value;
105 }
106 
107 
108 template<>
109 void writeVTK(OFstream& os, const scalar& value)
110 {
111  os << value;
112 }
113 
114 }
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 int main(int argc, char *argv[])
119 {
121  (
122  "Generate a legacy VTK file of particle tracks for cases that were"
123  " computed using a steady-state cloud"
124  );
125 
128  #include "addRegionOption.H"
129  #include "addDictOption.H"
130 
131  #include "setRootCase.H"
132 
133  #include "createTime.H"
135  #include "createNamedMesh.H"
136  #include "createFields.H"
137 
138  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140  const fileName vtkPath(runTime.rootPath()/runTime.globalCaseName()/"VTK");
141  mkDir(vtkPath);
142 
143  forAll(timeDirs, timeI)
144  {
145  runTime.setTime(timeDirs[timeI], timeI);
146  Info<< "Time = " << runTime.timeName() << endl;
147 
148  const fileName vtkTimePath(vtkPath/runTime.timeName());
149  mkDir(vtkTimePath);
150 
151  Info<< " Reading particle positions" << endl;
152 
153  PtrList<passiveParticle> particles(0);
154 
155  // Transfer particles to (more convenient) list
156  {
158  Info<< "\n Read " << returnReduce(ppc.size(), sumOp<label>())
159  << " particles" << endl;
160 
161  particles.setSize(ppc.size());
162 
163  label i = 0;
164  forAllIters(ppc, iter)
165  {
166  particles.set(i++, ppc.remove(&iter()));
167  }
168 
169  // myCloud should now be empty
170  }
171 
172  List<label> particleToTrack(particles.size());
173  label nTracks = 0;
174 
175  {
176  labelPairLookup trackTable;
177 
178  forAll(particles, i)
179  {
180  const label origProc = particles[i].origProc();
181  const label origId = particles[i].origId();
182 
183  const labelPair key(origProc, origId);
184 
185  const auto iter = trackTable.cfind(key);
186 
187  if (iter.found())
188  {
189  particleToTrack[i] = *iter;
190  }
191  else
192  {
193  particleToTrack[i] = nTracks;
194  trackTable.insert(key, nTracks);
195  ++nTracks;
196  }
197  }
198  }
199 
200 
201  if (nTracks == 0)
202  {
203  Info<< "\n No track data" << endl;
204  }
205  else
206  {
207  Info<< "\n Generating " << nTracks << " tracks" << endl;
208 
209  // Determine length of each track
210  labelList trackLengths(nTracks, Zero);
211  for (const label tracki : particleToTrack)
212  {
213  ++trackLengths[tracki];
214  }
215 
216  // Particle "age" property used to sort the tracks
217  List<SortableList<scalar>> agePerTrack(nTracks);
218  List<List<label>> particleMap(nTracks);
219 
220  forAll(trackLengths, i)
221  {
222  const label length = trackLengths[i];
223  agePerTrack[i].setSize(length);
224  particleMap[i].setSize(length);
225  }
226 
227  // Store the particle age per track
228  IOobjectList cloudObjs
229  (
230  mesh,
231  runTime.timeName(),
233  );
234 
235  // TODO: gather age across all procs
236  {
237  tmp<scalarField> tage =
238  readParticleField<scalar>("age", cloudObjs);
239 
240  const scalarField& age = tage();
241 
242  labelList trackSamples(nTracks, Zero);
243 
244  forAll(particleToTrack, i)
245  {
246  const label trackI = particleToTrack[i];
247  const label sampleI = trackSamples[trackI];
248  agePerTrack[trackI][sampleI] = age[i];
249  particleMap[trackI][sampleI] = i;
250  trackSamples[trackI]++;
251  }
252  tage.clear();
253  }
254 
255 
256  if (Pstream::master())
257  {
258  OFstream os(vtkTimePath/"particleTracks.vtk");
259 
260  Info<< "\n Writing particle tracks to " << os.name() << endl;
261 
262  label nPoints = sum(trackLengths);
263 
264  os << "# vtk DataFile Version 2.0" << nl
265  << "particleTracks" << nl
266  << "ASCII" << nl
267  << "DATASET POLYDATA" << nl
268  << "POINTS " << nPoints << " float" << nl;
269 
270  Info<< "\n Writing points" << endl;
271 
272  {
273  forAll(agePerTrack, i)
274  {
275  agePerTrack[i].sort();
276 
277  const labelList& ids = agePerTrack[i].indices();
278  labelList& particleIds = particleMap[i];
279 
280  {
281  // Update addressing
282  List<label> sortedIds(ids);
283  forAll(sortedIds, j)
284  {
285  sortedIds[j] = particleIds[ids[j]];
286  }
287  particleIds = sortedIds;
288  }
289 
290  forAll(ids, j)
291  {
292  const label localId = particleIds[j];
293  const vector pos(particles[localId].position());
294  os << pos.x() << ' ' << pos.y() << ' ' << pos.z()
295  << nl;
296  }
297  }
298  }
299 
300 
301  // Write track (line) connectivity to file
302 
303  Info<< "\n Writing track lines" << endl;
304  os << "\nLINES " << nTracks << ' ' << nPoints + nTracks << nl;
305 
306  // Write ids of track points to file
307  {
308  label globalPtI = 0;
309  forAll(particleMap, i)
310  {
311  os << particleMap[i].size() << nl;
312 
313  forAll(particleMap[i], j)
314  {
315  os << ' ' << globalPtI++;
316 
317  if (((j + 1) % 10 == 0) && (j != 0))
318  {
319  os << nl;
320  }
321  }
322 
323  os << nl;
324  }
325  }
326 
327 
328  const label nFields = validateFields(userFields, cloudObjs);
329 
330  os << "POINT_DATA " << nPoints << nl
331  << "FIELD attributes " << nFields << nl;
332 
333  Info<< "\n Processing fields" << nl << endl;
334 
335  processFields<label>(os, particleMap, userFields, cloudObjs);
336  processFields<scalar>(os, particleMap, userFields, cloudObjs);
337  processFields<vector>(os, particleMap, userFields, cloudObjs);
338  processFields<sphericalTensor>
339  (os, particleMap, userFields, cloudObjs);
340  processFields<symmTensor>
341  (os, particleMap, userFields, cloudObjs);
342  processFields<tensor>(os, particleMap, userFields, cloudObjs);
343 
344  }
345  }
346  Info<< endl;
347  }
348 
349  Info<< "End\n" << endl;
350 
351  return 0;
352 }
353 
354 
355 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::cloud::prefix
static const word prefix
Definition: cloud.H:83
Foam::TimePaths::globalCaseName
const fileName & globalCaseName() const
Definition: TimePathsI.H:49
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:280
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::OFstream::name
virtual const fileName & name() const
Definition: OSstream.H:103
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::glTF::key
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
Cloud.H
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::HashTable::insert
bool insert(const Key &key, const T &obj)
Definition: HashTableI.H:173
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Definition: UPstream.H:453
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::sumOp
Definition: ops.H:207
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::writeVTK
void writeVTK(OFstream &os, const Type &value)
OFstream.H
addDictOption.H
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
SortableList.H
labelPairHashes.H
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
steadyParticleTracksTemplates.H
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::Info
messageStream Info
argList.H
Field.H
addRegionOption.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::passiveParticleCloud
A Cloud of passive particles.
Definition: passiveParticleCloud.H:48
forAllIters
#define forAllIters(container, iter)
Definition: stdFoam.H:264
createNamedMesh.H
Required Variables.
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
fvMesh.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:49
userFields
List< word > userFields(propsDict.lookup("fields"))
Foam::HashTable< label, labelPair, Foam::Hash< labelPair > >
IOdictionary.H
Time.H
setRootCase.H
Foam::HashTable::cfind
const_iterator cfind(const Key &key) const
Definition: HashTableI.H:134
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:50
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Definition: timeSelector.C:95
passiveParticleCloud.H
Foam::Vector< scalar >
Foam::Time::rootPath
const fileName & rootPath() const
Definition: TimePathsI.H:43
Foam::List< word >
Foam::Time::setTime
virtual void setTime(const Time &t)
Definition: Time.C:996
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:320
timeSelector.H
createTime.H
PtrList.H
Foam::argList::noParallel
static void noParallel()
Definition: argList.C:503
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Definition: timeSelector.C:228
args
Foam::argList args(argc, argv)
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Definition: POSIX.C:533
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:170