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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  particleTracks
29 
30 Group
31  grpPostProcessingUtilities
32 
33 Description
34  Generate particle tracks for cases that were computed using a
35  tracked-parcel-type cloud.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "argList.H"
40 #include "Cloud.H"
41 #include "IOdictionary.H"
42 #include "fvMesh.H"
43 #include "Time.H"
44 #include "timeSelector.H"
45 #include "OFstream.H"
46 #include "passiveParticleCloud.H"
47 #include "writer.H"
48 #include "ListOps.H"
49 
50 #define createTrack(field, trackValues) \
51 createTrackField \
52 ( \
53  field, \
54  sampleFrequency, \
55  maxPositions, \
56  startIds, \
57  allOrigProcs, \
58  allOrigIds, \
59  trackValues \
60 );
61 
62 #define setFields(fields, fieldNames) \
63 setTrackFields \
64 ( \
65  obr, \
66  fields, \
67  fieldNames, \
68  nTracks, \
69  startIds, \
70  allOrigProcs, \
71  allOrigIds, \
72  maxPositions, \
73  sampleFrequency \
74 );
75 
76 #define writeFields(fields, fieldNames, tracks, times, dirs) \
77 writeTrackFields \
78 ( \
79  fields, \
80  fieldNames, \
81  tracks, \
82  times, \
83  dirs, \
84  setFormat, \
85  formatOptions, \
86  cloudName \
87 );
88 
89 using namespace Foam;
90 
91 template<class Type>
92 void createTrackField
93 (
94  const Field<Type>& values,
95  const label sampleFrequency,
96  const label maxPositions,
97  const labelList& startIds,
98  const List<labelField>& allOrigProcs,
99  const List<labelField>& allOrigIds,
100  List<DynamicList<Type>>& trackValues
101 )
102 {
103  List<Field<Type>> procField(Pstream::nProcs());
104  procField[Pstream::myProcNo()] = values;
105  Pstream::gatherList(procField);
106 
107  if (!Pstream::master())
108  {
109  return;
110  }
111 
112  const label nTracks = trackValues.size();
113 
114  forAll(procField, proci)
115  {
116  forAll(procField[proci], i)
117  {
118  const label globalId =
119  startIds[allOrigProcs[proci][i]] + allOrigIds[proci][i];
120 
121  if (globalId % sampleFrequency == 0)
122  {
123  const label trackId = globalId/sampleFrequency;
124 
125  if
126  (
127  trackId < nTracks
128  && trackValues[trackId].size() < maxPositions
129  )
130  {
131  trackValues[trackId].append(procField[proci][i]);
132  }
133  }
134  }
135  }
136 }
137 
138 
139 template<class Type>
140 void writeTrackFields
141 (
142  List<List<DynamicList<Type>>>& fieldValues,
143  const wordList& fieldNames,
144  const PtrList<coordSet>& tracks,
145  const List<scalarField>& times,
146  const List<vectorField>& dirs,
147  const word& setFormat,
148  const dictionary& formatOptions,
149  const word& cloudName
150 )
151 {
152  if (fieldValues.empty())
153  {
154  return;
155  }
156 
157  auto writerPtr = writer<Type>::New(setFormat, formatOptions);
158 
159  const fileName outFile(writerPtr().getFileName(tracks[0], wordList(0)));
160 
161  const fileName outPath
162  (
164  );
165  mkDir(outPath);
166 
167  OFstream os(outPath/(pTraits<Type>::typeName & "tracks." + outFile.ext()));
168 
169  Info<< "Writing " << pTraits<Type>::typeName << " particle tracks in "
170  << setFormat << " format to " << os.name() << endl;
171 
172 
173  List<List<Field<Type>>> fields(fieldValues.size());
174  forAll(fields, fieldi)
175  {
176  fields[fieldi].setSize(fieldValues[fieldi].size());
177  forAll(fields[fieldi], tracki)
178  {
179  fields[fieldi][tracki].transfer(fieldValues[fieldi][tracki]);
180  }
181  }
182 
183  writerPtr().write(true, times, tracks, fieldNames, fields, os);
184 }
185 
186 
187 template<class Type>
188 Foam::label setTrackFields
189 (
190  const objectRegistry& obr,
193  const label nTracks,
194  const labelList& startIds,
195  const List<labelField>& allOrigProcs,
196  const List<labelField>& allOrigIds,
197  const label maxPositions,
198  const label sampleFrequency
199 )
200 {
201  const auto availableFieldPtrs = obr.lookupClass<IOField<Type>>();
202 
203  fieldNames = availableFieldPtrs.toc();
204 
205  if (Pstream::parRun())
206  {
210  }
211 
212  const label nFields = fieldNames.size();
213 
214  if (fields.empty())
215  {
216  fields.setSize(nFields);
217  fieldNames.setSize(nFields);
218  forAll(fields, i)
219  {
220  fields[i].setSize(nTracks);
221  }
222  }
223 
224  forAll(fieldNames, fieldi)
225  {
226  const word& fieldName = fieldNames[fieldi];
227 
228  const auto* fldPtr = obr.cfindObject<IOField<Type>>(fieldName);
229 
230  createTrack
231  (
232  fldPtr ? static_cast<const Field<Type>>(*fldPtr) : Field<Type>(),
233  fields[fieldi]
234  );
235  }
236 
237  return nFields;
238 }
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 int main(int argc, char *argv[])
244 {
246  (
247  "Generate a file of particle tracks for cases that were"
248  " computed using a tracked-parcel-type cloud"
249  );
250 
252  #include "addRegionOption.H"
253 
254  #include "setRootCase.H"
255 
256  #include "createTime.H"
258  #include "createNamedMesh.H"
259  #include "createFields.H"
260 
261  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263  const fileName vtkPath(runTime.rootPath()/runTime.globalCaseName()/"VTK");
264  mkDir(vtkPath);
265 
266  Info<< "Scanning times to determine track data for cloud " << cloudName
267  << nl << endl;
268 
269  labelList maxIds(Pstream::nProcs(), -1);
270  forAll(timeDirs, timei)
271  {
272  runTime.setTime(timeDirs[timei], timei);
273  Info<< "Time = " << runTime.timeName() << endl;
274 
275  Info<< " Reading particle positions" << endl;
277 
278  Info<< " Read " << returnReduce(myCloud.size(), sumOp<label>())
279  << " particles" << endl;
280 
281  for (const passiveParticle& p : myCloud)
282  {
283  const label origId = p.origId();
284  const label origProc = p.origProc();
285 
286  // Handle case where we are processing particle data generated in
287  // parallel using more cores than when running this application.
288  if (origProc >= maxIds.size())
289  {
290  maxIds.setSize(origProc+1, -1);
291  }
292 
293  maxIds[origProc] = max(maxIds[origProc], origId);
294  }
295  }
296 
297  label maxNProcs = returnReduce(maxIds.size(), maxOp<label>());
298 
299  Info<< "Detected particles originating from " << maxNProcs
300  << " processors." << nl << endl;
301 
302  maxIds.setSize(maxNProcs, -1);
303 
306 
307  labelList numIds = maxIds + 1;
308 
309  Info<< nl << "Particle statistics:" << endl;
310  forAll(maxIds, proci)
311  {
312  Info<< " Found " << numIds[proci] << " particles originating"
313  << " from processor " << proci << endl;
314  }
315  Info<< nl << endl;
316 
317 
318  // Calculate starting ids for particles on each processor
319  labelList startIds(numIds.size(), Zero);
320  for (label i = 0; i < numIds.size()-1; ++i)
321  {
322  startIds[i+1] += startIds[i] + numIds[i];
323  }
324  label nParticle = startIds.last() + numIds[startIds.size()-1];
325 
326 
327  // Number of tracks to generate
328  const label nTracks =
329  maxTracks > 0
330  ? min(nParticle/sampleFrequency, maxTracks)
331  : nParticle/sampleFrequency;
332 
333  // Storage for all particle tracks
334  List<DynamicList<vector>> allTracks(nTracks);
335  List<DynamicList<scalar>> allTrackTimes(nTracks);
336 
337  // Lists of field, tracki, trackValues
338  //List<List<DynamicList<label>>> labelFields;
339  List<List<DynamicList<scalar>>> scalarFields;
340  List<List<DynamicList<vector>>> vectorFields;
341  List<List<DynamicList<sphericalTensor>>> sphTensorFields;
342  List<List<DynamicList<symmTensor>>> symTensorFields;
343  List<List<DynamicList<tensor>>> tensorFields;
344  //List<word> labelFieldNames;
345  List<word> scalarFieldNames;
346  List<word> vectorFieldNames;
347  List<word> sphTensorFieldNames;
348  List<word> symTensorFieldNames;
349  List<word> tensorFieldNames;
350 
351  Info<< "\nGenerating " << nTracks << " particle tracks for cloud "
352  << cloudName << nl << endl;
353 
354  forAll(timeDirs, timei)
355  {
356  runTime.setTime(timeDirs[timei], timei);
357  Info<< "Time = " << runTime.timeName() << endl;
358 
359  // Read particles. Will be size 0 if no particles.
360  Info<< " Reading particle positions" << endl;
362 
363  pointField localPositions(myCloud.size(), Zero);
364  scalarField localTimes(myCloud.size(), Zero);
365 
366  List<labelField> allOrigIds(Pstream::nProcs());
367  List<labelField> allOrigProcs(Pstream::nProcs());
368 
369  // Collect the track data on all processors that have positions
370  allOrigIds[Pstream::myProcNo()].setSize(myCloud.size(), Zero);
371  allOrigProcs[Pstream::myProcNo()].setSize(myCloud.size(), Zero);
372 
373  label i = 0;
374  for (const passiveParticle& p : myCloud)
375  {
376  allOrigIds[Pstream::myProcNo()][i] = p.origId();
377  allOrigProcs[Pstream::myProcNo()][i] = p.origProc();
378  localPositions[i] = p.position();
379  localTimes[i] = runTime.value();
380  ++i;
381  }
382 
383  // Collect the track data on the master processor
384  Pstream::gatherList(allOrigIds);
385  Pstream::gatherList(allOrigProcs);
386 
387  objectRegistry obr
388  (
389  IOobject
390  (
391  "cloudFields",
392  runTime.timeName(),
393  runTime
394  )
395  );
396 
397  myCloud.readFromFiles(obr, fieldNames);
398 
399  // Create track positions and track time fields
400  // (not registered as IOFields)
401  // Note: createTrack is a local #define to reduce arg count...
402  createTrack(localPositions, allTracks);
403  createTrack(localTimes, allTrackTimes);
404 
405  // Create the track fields
406  // Note: setFields is a local #define to reduce arg count...
407  //setFields(labelFields, labelFieldNames);
408  setFields(scalarFields, scalarFieldNames);
409  setFields(vectorFields, vectorFieldNames);
410  setFields(sphTensorFields, sphTensorFieldNames);
411  setFields(symTensorFields, symTensorFieldNames);
412  setFields(tensorFields, tensorFieldNames);
413  }
414 
415 
416  if (Pstream::master())
417  {
418  PtrList<coordSet> tracks(allTracks.size());
419  List<scalarField> times(tracks.size());
420  forAll(tracks, tracki)
421  {
422  tracks.set
423  (
424  tracki,
425  new coordSet("track" + Foam::name(tracki), "distance")
426  );
427  tracks[tracki].transfer(allTracks[tracki]);
428  times[tracki].transfer(allTrackTimes[tracki]);
429  }
430 
431  Info<< nl;
432 
433  const label Uid = vectorFieldNames.find(UName);
434  List<vectorField> dirs(nTracks);
435 
436  if (Uid != -1)
437  {
438  const auto& UTracks = vectorFields[Uid];
439  forAll(UTracks, tracki)
440  {
441  const auto& U = UTracks[tracki];
442  dirs[tracki] = U/(mag(U) + ROOTVSMALL);
443  }
444  }
445 
446  // Write track fields
447  // Note: writeFields is a local #define to reduce arg count...
448  //writeFields(allLabelFields, labelFieldNames, tracks);
449  writeFields(scalarFields, scalarFieldNames, tracks, times, dirs);
450  writeFields(vectorFields, vectorFieldNames, tracks, times, dirs);
451  writeFields(sphTensorFields, sphTensorFieldNames, tracks, times, dirs);
452  writeFields(symTensorFields, symTensorFieldNames, tracks, times, dirs);
453  writeFields(tensorFields, tensorFieldNames, tracks, times, dirs);
454  }
455 
456  Info<< nl << "End\n" << endl;
457 
458  return 0;
459 }
460 
461 
462 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:87
Foam::TimePaths::globalCaseName
const fileName & globalCaseName() const
Return global case name.
Definition: TimePathsI.H:56
Foam::maxOp
Definition: ops.H:223
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
sampleFrequency
const label sampleFrequency(propsDict.get< label >("sampleFrequency"))
p
volScalarField & p
Definition: createFieldRefs.H:8
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::OFstream::name
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< Type >
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
setFormat
const word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))
UName
const word UName(propsDict.getOrDefault< word >("U", "U"))
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::passiveParticle
Copy of base particle.
Definition: passiveParticle.H:53
Cloud.H
Foam::argList::addNote
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:412
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
maxPositions
const label maxPositions(propsDict.get< label >("maxPositions"))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (can be nullptr),.
Definition: PtrList.H:138
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
maxTracks
const label maxTracks(propsDict.getOrDefault< label >("maxTracks", -1))
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
formatOptions
const dictionary formatOptions
Definition: createFields.H:26
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::objectRegistry::lookupClass
HashTable< const Type * > lookupClass(const bool strict=false) const
Return all objects with a class satisfying isA<Type>
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
argList.H
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
addRegionOption.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:376
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::passiveParticleCloud
A Cloud of passive particles.
Definition: passiveParticleCloud.H:52
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
createNamedMesh.H
Required Variables.
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::Pstream::combineGather
static void combineGather(const List< commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:48
Foam::ListOps::uniqueEqOp
List helper to append y unique elements onto the end of x.
Definition: ListOps.H:577
Foam::maxEqOp
Definition: ops.H:80
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
IOdictionary.H
Time.H
U
U
Definition: pEqn.H:72
setRootCase.H
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::objectRegistry::cfindObject
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:390
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:102
passiveParticleCloud.H
Foam::Time::rootPath
const fileName & rootPath() const
Return root path.
Definition: TimePathsI.H:50
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::Time::setTime
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:1003
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
timeSelector.H
createTime.H
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
ListOps.H
Various functions to operate on Lists.
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
writer.H
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Definition: timeSelector.C:235
args
Foam::argList args(argc, argv)
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::Pstream::combineScatter
static void combineScatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:183
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::writeFields
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
Foam::PtrList::transfer
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition: PtrListI.H:245
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:38