particleTracks.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 "setWriter.H"
42 #include "systemDict.H"
43 
44 using namespace Foam;
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
51  #include "addRegionOption.H"
52 
53  #include "setRootCase.H"
54 
55  #include "createTime.H"
57  #include "createNamedMesh.H"
58  #include "createFields.H"
59 
60  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62  const fileName vtkPath(runTime.rootPath()/runTime.globalCaseName()/"VTK");
63  mkDir(vtkPath);
64 
65  Info<< "Scanning times to determine track data for cloud " << cloudName
66  << nl << endl;
67 
68  labelList maxIds(Pstream::nProcs(), -1);
69  forAll(timeDirs, timeI)
70  {
71  runTime.setTime(timeDirs[timeI], timeI);
72  Info<< "Time = " << runTime.userTimeName() << endl;
73 
74  mesh.readUpdate();
75 
76  Info<< " Reading particle positions" << endl;
77  passiveParticleCloud myCloud(mesh, cloudName);
78 
79  Info<< " Read " << returnReduce(myCloud.size(), sumOp<label>())
80  << " particles" << endl;
81 
83  {
84  label origId = iter().origId();
85  label origProc = iter().origProc();
86 
87  if (origProc >= maxIds.size())
88  {
89  maxIds.setSize(origProc+1, -1);
90  }
91 
92  maxIds[origProc] = max(maxIds[origProc], origId);
93  }
94  }
95 
96  label maxNProcs = returnReduce(maxIds.size(), maxOp<label>());
97 
98  Info<< "Detected particles originating from " << maxNProcs
99  << " processors." << nl << endl;
100 
101  maxIds.setSize(maxNProcs, -1);
102 
105 
106  labelList numIds = maxIds + 1;
107 
108  Info<< nl << "Particle statistics:" << endl;
109  forAll(maxIds, proci)
110  {
111  Info<< " Found " << numIds[proci] << " particles originating"
112  << " from processor " << proci << endl;
113  }
114  Info<< nl << endl;
115 
116 
117  // Calculate starting ids for particles on each processor
118  List<label> startIds(numIds.size(), 0);
119  for (label i = 0; i < numIds.size()-1; i++)
120  {
121  startIds[i+1] += startIds[i] + numIds[i];
122  }
123  label nParticle = startIds.last() + numIds[startIds.size()-1];
124 
125 
126  // Number of tracks to generate
127  label nTracks = nParticle/sampleFrequency;
128 
129  // Storage for all particle tracks
130  List<DynamicList<vector>> allTracks(nTracks);
131 
132  Info<< "\nGenerating " << nTracks << " particle tracks for cloud "
133  << cloudName << nl << endl;
134 
135  forAll(timeDirs, timeI)
136  {
137  runTime.setTime(timeDirs[timeI], timeI);
138  Info<< "Time = " << runTime.userTimeName() << endl;
139 
140  mesh.readUpdate();
141 
142  List<pointField> allPositions(Pstream::nProcs());
143  List<labelField> allOrigIds(Pstream::nProcs());
144  List<labelField> allOrigProcs(Pstream::nProcs());
145 
146  // Read particles. Will be size 0 if no particles.
147  Info<< " Reading particle positions" << endl;
148  passiveParticleCloud myCloud(mesh, cloudName);
149 
150  // Collect the track data on all processors that have positions
151  allPositions[Pstream::myProcNo()].setSize
152  (
153  myCloud.size(),
155  );
156  allOrigIds[Pstream::myProcNo()].setSize(myCloud.size(), 0);
157  allOrigProcs[Pstream::myProcNo()].setSize(myCloud.size(), 0);
158 
159  label i = 0;
160  forAllConstIter(passiveParticleCloud, myCloud, iter)
161  {
162  allPositions[Pstream::myProcNo()][i] = iter().position(mesh);
163  allOrigIds[Pstream::myProcNo()][i] = iter().origId();
164  allOrigProcs[Pstream::myProcNo()][i] = iter().origProc();
165  i++;
166  }
167 
168  // Collect the track data on the master processor
169  Pstream::gatherList(allPositions);
170  Pstream::gatherList(allOrigIds);
171  Pstream::gatherList(allOrigProcs);
172 
173  Info<< " Constructing tracks" << nl << endl;
174  if (Pstream::master())
175  {
176  forAll(allPositions, proci)
177  {
178  forAll(allPositions[proci], i)
179  {
180  label globalId =
181  startIds[allOrigProcs[proci][i]]
182  + allOrigIds[proci][i];
183 
184  if (globalId % sampleFrequency == 0)
185  {
186  label trackId = globalId/sampleFrequency;
187  if (allTracks[trackId].size() < maxPositions)
188  {
189  allTracks[trackId].append
190  (
191  allPositions[proci][i]
192  );
193  }
194  }
195  }
196  }
197  }
198  }
199 
200  if (allTracks.size() && Pstream::master())
201  {
202  DynamicList<point> allTrack;
203  DynamicList<label> allTrackIDs;
204  forAll(allTracks, trackI)
205  {
206  allTrack.append(allTracks[trackI]);
207  allTrackIDs.append(labelList(allTracks[trackI].size(), trackI));
208  }
209 
210  Info<< "\nWriting particle tracks in " << setFormat
211  << " format to " << vtkPath << nl << endl;
212 
214  (
215  vtkPath,
216  "tracks",
217  coordSet(allTrackIDs, word::null, pointField(allTrack))
218  );
219  }
220 
221  Info<< "End\n" << endl;
222 
223  return 0;
224 }
225 
226 
227 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
T & last()
Return the last element of the list.
Definition: UListI.H:128
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
static const Form zero
Definition: VectorSpace.H:113
Holds list of sampling positions.
Definition: coordSet.H:51
A class for handling file names.
Definition: fileName.H:82
A Cloud of passive particles.
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
Definition: setWriter.C:286
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
static const word null
An empty word.
Definition: word.H:77
int main(int argc, char *argv[])
Definition: financialFoam.C:44
static instantList timeDirs
Definition: globalFoam.H:44
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static const char nl
Definition: Ostream.H:260
Foam::argList args(argc, argv)
IOdictionary propsDict(systemDict("particleTracksDict", args, runTime))
label sampleFrequency(propsDict.lookup< label >("sampleFrequency"))
word setFormat(propsDict.lookup< word >("setFormat"))
label maxPositions(propsDict.lookup< label >("maxPositions"))
const word cloudName(propsDict.lookup("cloudName"))