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-2022 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 "coordSetWriter.H"
46 #include "passiveParticleCloud.H"
47 #include "particleTracksSampler.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 using namespace Foam;
52 
53 template<class Type>
54 void writeTrackFields
55 (
57  HashTable<List<DynamicList<Type>>>& fieldTable
58 )
59 {
60  for (const word& fieldName : fieldTable.sortedToc())
61  {
62  // Steal and reshape from List<DynamicList> to List<Field>
63  auto& input = fieldTable[fieldName];
64 
65  List<Field<Type>> fields(input.size());
66  forAll(input, tracki)
67  {
68  fields[tracki].transfer(input[tracki]);
69  }
70 
71  writer.write(fieldName, fields);
72  }
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 int main(int argc, char *argv[])
79 {
81  (
82  "Generate a file of particle tracks for cases that were"
83  " computed using a tracked-parcel-type cloud"
84  );
86  #include "addRegionOption.H"
87 
88  // Less frequently used - reduce some clutter
89  argList::setAdvanced("decomposeParDict");
90 
92  (
93  "dict",
94  "file",
95  "Alternative particleTracksProperties dictionary"
96  );
98  (
99  "stride",
100  "int",
101  "Override the sample-frequency"
102  );
104  (
105  "fields",
106  "wordRes",
107  "Specify single or multiple fields to write "
108  "(default: all or 'fields' from dictionary)\n"
109  "Eg, 'T' or '( \"U.*\" )'"
110  );
112  (
113  "format",
114  "name",
115  "The writer format "
116  "(default: vtk or 'setFormat' from dictionary)"
117  );
118  argList::addVerboseOption("Additional verbosity");
119 
120  #include "setRootCase.H"
121  #include "createTime.H"
123  #include "createNamedMesh.H"
124 
125  // ------------------------------------------------------------------------
126  // Control properties
127 
128  #include "createControls.H"
129 
131  args.readIfPresent("format", setFormat);
132 
134  sampleFrequency = max(1, sampleFrequency); // sanity
135 
136  // Setup the writer
137  auto writerPtr =
139  (
140  setFormat,
142  );
143 
144  writerPtr().useTracks(true);
145 
146  if (args.verbose())
147  {
148  writerPtr().verbose(true);
149  }
150 
151  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153  particleTracksSampler trackSampler;
154 
155  Info<< "Scanning times to determine track data for cloud " << cloudName
156  << nl << endl;
157 
158  {
159  labelList maxIds(Pstream::nProcs(), -1);
160  forAll(timeDirs, timei)
161  {
162  runTime.setTime(timeDirs[timei], timei);
163  Info<< "Time = " << runTime.timeName() << endl;
164 
166 
167  Info<< " Read " << returnReduce(myCloud.size(), sumOp<label>())
168  << " particles" << endl;
169 
170  for (const passiveParticle& p : myCloud)
171  {
172  const label origId = p.origId();
173  const label origProc = p.origProc();
174 
175  // Handle case where processing particle data generated in
176  // parallel using more cores than for this application.
177  if (origProc >= maxIds.size())
178  {
179  maxIds.resize(origProc+1, -1);
180  }
181 
182  maxIds[origProc] = max(maxIds[origProc], origId);
183  }
184  }
185 
186  const label maxNProcs = returnReduce(maxIds.size(), maxOp<label>());
187  maxIds.resize(maxNProcs, -1);
188 
190 
191  // From ids to count
192  const labelList numIds = maxIds + 1;
193 
194  // Set parcel addressing
195  trackSampler.reset(numIds);
196 
197  Info<< nl
198  << "Detected particles originating from "
199  << maxNProcs << " processors." << nl
200  << "Particle statistics:" << endl;
201 
202  if (Pstream::master())
203  {
204  const globalIndex& parcelAddr = trackSampler.parcelAddr();
205 
206  for (const label proci : parcelAddr.allProcs())
207  {
208  Info<< " Found " << parcelAddr.localSize(proci)
209  << " particles originating"
210  << " from processor " << proci << nl;
211  }
212  }
213  }
214 
216 
217 
218  // Number of tracks to generate
219  const label nTracks = trackSampler.nTracks();
220 
221 
222  // Storage for all particle tracks
223  List<DynamicList<point>> allTrackPos(nTracks);
224  List<DynamicList<scalar>> allTrackTimes(nTracks);
225 
226  // Track field values by name/type
227  HashTable<List<DynamicList<label>>> labelFields; // <- mostly unused
228  HashTable<List<DynamicList<scalar>>> scalarFields;
229  HashTable<List<DynamicList<vector>>> vectorFields;
230  HashTable<List<DynamicList<sphericalTensor>>> sphericalTensorFields;
231  HashTable<List<DynamicList<symmTensor>>> symmTensorFields;
232  HashTable<List<DynamicList<tensor>>> tensorFields;
233 
234  Info<< nl
235  << "Generating " << nTracks
236  << " particle tracks for cloud " << cloudName << nl << endl;
237 
238  forAll(timeDirs, timei)
239  {
240  runTime.setTime(timeDirs[timei], timei);
241  Info<< "Time = " << runTime.timeName() << " (processing)" << endl;
242 
243  // Read particles. Will be size 0 if no particles.
245 
246  pointField localPositions(myCloud.size());
247  scalarField localTimes(myCloud.size(), runTime.value());
248 
249  // Gather track data from all processors that have positions
250  trackSampler.resetCloud(myCloud.size());
251  {
252  labelField& origIds = trackSampler.origParcelIds_;
253  labelField& origProcs = trackSampler.origProcIds_;
254 
255  label np = 0;
256  for (const passiveParticle& p : myCloud)
257  {
258  origIds[np] = p.origId();
259  origProcs[np] = p.origProc();
260  localPositions[np] = p.position();
261  ++np;
262  }
263 
264  trackSampler.gatherInplace(origIds);
265  trackSampler.gatherInplace(origProcs);
266  }
267 
268 
269  // Read cloud fields (from disk) into object registry
270  objectRegistry obr
271  (
272  IOobject
273  (
274  "cloudFields",
275  runTime.timeName(),
276  runTime
277  )
278  );
279 
280  myCloud.readFromFiles(obr, acceptFields, excludeFields);
281 
282  // Create track positions and track time fields
283  // (not registered as IOFields)
284 
285  trackSampler.createTrackField(localPositions, allTrackPos);
286  trackSampler.createTrackField(localTimes, allTrackTimes);
287 
288  // Create the track fields
290  trackSampler.setTrackFields(obr, scalarFields);
291  trackSampler.setTrackFields(obr, vectorFields);
292  trackSampler.setTrackFields(obr, sphericalTensorFields);
293  trackSampler.setTrackFields(obr, symmTensorFields);
294  trackSampler.setTrackFields(obr, tensorFields);
295  }
296 
297  const label nFields =
298  (
299  labelFields.size()
300  + scalarFields.size()
301  + vectorFields.size()
302  + sphericalTensorFields.size()
303  + symmTensorFields.size()
304  + tensorFields.size()
305  );
306 
307  Info<< nl
308  << "Extracted " << nFields << " cloud fields" << nl;
309 
310  if (nFields)
311  {
312  #undef doLocalCode
313  #define doLocalCode(FieldContent) \
314  if (!FieldContent.empty()) \
315  { \
316  Info<< " "; \
317  for (const word& fieldName : FieldContent.sortedToc()) \
318  { \
319  Info<< ' ' << fieldName; \
320  } \
321  Info<< nl; \
322  }
323 
324  doLocalCode(labelFields);
325  doLocalCode(scalarFields);
326  doLocalCode(vectorFields);
327  doLocalCode(sphericalTensorFields);
328  doLocalCode(symmTensorFields);
329  doLocalCode(tensorFields);
330  #undef doLocalCode
331  }
332 
333  Info<< nl
334  << "Writing particle tracks (" << setFormat << " format)" << nl;
335 
336  if (Pstream::master())
337  {
338  PtrList<coordSet> tracks(allTrackPos.size());
339  List<scalarField> times(allTrackPos.size());
340  forAll(tracks, tracki)
341  {
342  tracks.set
343  (
344  tracki,
345  new coordSet("track" + Foam::name(tracki), "xyz")
346  );
347  tracks[tracki].transfer(allTrackPos[tracki]);
348  times[tracki].transfer(allTrackTimes[tracki]);
349 
350  if (!tracki) tracks[0].rename("tracks");
351  }
352 
365 
366 
367  // Write tracks with fields
368  if (nFields)
369  {
370  auto& writer = *writerPtr;
371 
372  const fileName outputPath
373  (
375  / "particleTracks" / "tracks"
376  );
377 
378  Info<< " "
379  << runTime.relativePath(outputPath) << endl;
380 
381  writer.open(tracks, outputPath);
382  writer.setTrackTimes(times);
383  writer.nFields(nFields);
384 
386  writeTrackFields(writer, scalarFields);
387  writeTrackFields(writer, vectorFields);
388  writeTrackFields(writer, sphericalTensorFields);
389  writeTrackFields(writer, symmTensorFields);
390  writeTrackFields(writer, tensorFields);
391  }
392  else
393  {
394  Info<< "Warning: no fields, did not write" << endl;
395  }
396  }
397 
398  Info<< nl << "End\n" << endl;
399 
400  return 0;
401 }
402 
403 
404 // ************************************************************************* //
labelField origProcIds_
The originating processor ids.
const Type & value() const noexcept
Return const reference to value.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
A class for handling file names.
Definition: fileName.H:71
static void setAdvanced(const word &optName, bool advanced=true)
Set an existing option as being &#39;advanced&#39; or normal.
Definition: argList.C:395
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:132
labelField origParcelIds_
The originating parcel ids.
wordRes acceptFields
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Required Variables.
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: TimePathsI.H:80
void gatherInplace(List< Type > &fld) const
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:45
void reset(const labelUList &origParcelCounts)
Define the orig parcel mappings.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
A Cloud of passive particles.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
label setTrackFields(const objectRegistry &obr, HashTable< List< DynamicList< Type >>> &fieldTable) const
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:572
wordRes excludeFields
static autoPtr< coordSetWriter > New(const word &writeFormat)
Return a reference to the selected writer.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
#define doLocalCode(GeoField)
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator) is 1 for serial run.
Definition: UPstream.H:656
dynamicFvMesh & mesh
Holds list of sampling positions.
Definition: coordSet.H:49
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const word cloudName(propsDict.get< word >("cloud"))
label setSampleRate(const label sampleFreq, const label maxPositions, const label maxTracks=-1)
Set the sampling stride, upper limits.
A class for handling words, derived from Foam::string.
Definition: word.H:63
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:48
Helper class when generating particle tracks. The interface is fairly rudimentary.
Base class for writing coordSet(s) and tracks with fields.
static void addVerboseOption(const string &usage, bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:505
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:977
label maxPositions(propsDict.get< label >("maxPositions"))
String literal.
Definition: keyType.H:82
A HashTable similar to std::unordered_map.
Definition: HashTable.H:102
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
label nTracks() const noexcept
Number of tracks to generate.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
bool readListIfPresent(const word &optName, List< T > &list) const
If named option is present, get a List of values treating a single entry like a list of size 1...
Definition: argListI.H:387
virtual bool open(const fileName &file, bool parallel=UPstream::parRun())
Open file for writing (creates parent directory).
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:234
const globalIndex & parcelAddr() const noexcept
The original parcel addressing.
labelRange allProcs() const noexcept
Range of process indices for all addressed offsets (processes)
Definition: globalIndexI.H:163
label maxTracks(propsDict.getOrDefault< label >("maxTracks", -1))
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
label sampleFrequency(propsDict.get< label >("sampleFrequency"))
void createTrackField(const UList< Type > &values, List< DynamicList< Type >> &trackValues) const
static bool master(const label communicator=worldComm)
Am I the master rank.
Definition: UPstream.H:672
static word outputPrefix
Directory prefix.
void resetCloud(const label localCloudSize)
messageStream Info
Information stream (stdout output on master, null elsewhere)
word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
volScalarField & p
Registry of regIOobjects.
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
Copy of base particle.
int verbose() const noexcept
Return the verbose flag.
Definition: argListI.H:121
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.
const dictionary formatOptions(propsDict.subOrEmptyDict("formatOptions", keyType::LITERAL))
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
label localSize() const
My local size.
Definition: globalIndexI.H:225
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:93