Go to the documentation of this file.
47 int main(
int argc,
char *argv[])
57 #include "createFields.H"
61 fileName vtkPath(runTime.path()/
"VTK");
64 Info<<
"Scanning times to determine track data for cloud " <<
cloudName
70 runTime.setTime(
timeDirs[timeI], timeI);
71 Info<<
"Time = " << runTime.timeName() <<
endl;
73 Info<<
" Reading particle positions" <<
endl;
76 <<
" particles" <<
endl;
80 label origId = iter().origId();
81 label origProc = iter().origProc();
83 if (origProc >= maxIds.size())
86 maxIds.setSize(origProc+1, -1);
89 maxIds[origProc] =
max(maxIds[origProc], origId);
95 Info<<
"Detected particles originating from " << maxNProcs
96 <<
" processors." <<
nl <<
endl;
98 maxIds.setSize(maxNProcs, -1);
108 Info<<
" Found " << numIds[procI] <<
" particles originating"
109 <<
" from processor " << procI <<
endl;
116 for (
label i = 0; i < numIds.
size()-1; i++)
118 startIds[i+1] += startIds[i] + numIds[i];
120 label nParticle = startIds.last() + numIds[startIds.size()-1];
130 Info<<
"\nGenerating " << nTracks <<
" particle tracks for cloud "
135 runTime.setTime(
timeDirs[timeI], timeI);
136 Info<<
"Time = " << runTime.timeName() <<
endl;
143 Info<<
" Reading particle positions" <<
endl;
172 forAll(allPositions, procI)
174 forAll(allPositions[procI], i)
177 startIds[allOrigProcs[procI][i]]
178 + allOrigIds[procI][i];
185 allTracks[trackId].append
187 allPositions[procI][i]
211 tracks[trackI].transfer(allTracks[trackI]);
222 scalarFormatterPtr().getFileName
232 /
"particleTracks." + vtkFile.ext()
236 <<
" format to " << vtkTracks.name()
240 scalarFormatterPtr().write
word setFormat(propsDict.lookupOrDefault< word >("setFormat", "vtk"))
label sampleFrequency(readLabel(propsDict.lookup("sampleFrequency")))
A class for handling file names.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
#define forAll(list, i)
Loop across all elements in list.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Ostream & endl(Ostream &os)
Add newline and flush stream.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
List< word > wordList
A List of words.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
label maxPositions(readLabel(propsDict.lookup("maxPositions")))
int main(int argc, char *argv[])
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A Cloud of passive particles.
Holds list of sampling positions.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
static bool master(const label communicator=0)
Am I the master process.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
static instantList timeDirs
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const word cloudName(propsDict.lookup("cloudName"))
void size(const label)
Override size to be inconsistent with allocated storage.
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Foam::argList args(argc, argv)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
word name(const complex &)
Return a string representation of a complex.
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.