40 #define MPI_SCALAR MPI_FLOAT 42 #define MPI_SCALAR MPI_DOUBLE 53 validParOptions.insert(
"np",
"");
54 validParOptions.insert(
"p4pg",
"PI file");
55 validParOptions.insert(
"p4wd",
"directory");
56 validParOptions.insert(
"p4amslave",
"");
57 validParOptions.insert(
"p4yourname",
"hostname");
58 validParOptions.insert(
"machinefile",
"machine file");
65 int provided_thread_support;
71 &provided_thread_support
75 MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
77 MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
81 Pout<<
"UPstream::init : initialised with numProcs:" << numprocs
82 <<
" myRank:" << myRank <<
endl;
88 <<
"bool IPstream::init(int& argc, char**& argv) : " 89 "attempt to run parallel on 1 processor" 97 if (
Pstream::master() && provided_thread_support != MPI_THREAD_MULTIPLE)
100 <<
"mpi does not seem to have thread support." 101 <<
" There might be issues with e.g. threaded IO" 107 string bufferSizeName =
getEnv(
"MPI_BUFFER_SIZE");
109 if (bufferSizeName.size())
111 int bufferSize = atoi(bufferSizeName.c_str());
115 MPI_Buffer_attach(
new char[bufferSize], bufferSize);
121 <<
"UPstream::init(int& argc, char**& argv) : " 122 <<
"environment variable MPI_BUFFER_SIZE not defined" 148 MPI_Buffer_detach(&buff, &size);
158 <<
"There are still " << n <<
" outstanding MPI_Requests." <<
endl 159 <<
"This means that your code exited before doing a" 160 <<
" UPstream::waitRequests()." <<
endl 161 <<
"This should not happen for a normal code exit." 166 forAll(myProcNo_, communicator)
168 if (myProcNo_[communicator] != -1)
170 freePstreamCommunicator(communicator);
181 MPI_Abort(MPI_COMM_WORLD, errnum);
188 MPI_Abort(MPI_COMM_WORLD, 1);
195 const sumOp<scalar>& bop,
197 const label communicator
202 Pout<<
"** reducing:" << Value <<
" with comm:" << communicator
207 allReduce(Value, 1, MPI_SCALAR, MPI_SUM, bop, tag, communicator);
214 const minOp<scalar>& bop,
216 const label communicator
221 Pout<<
"** reducing:" << Value <<
" with comm:" << communicator
226 allReduce(Value, 1, MPI_SCALAR, MPI_MIN, bop, tag, communicator);
233 const sumOp<vector2D>& bop,
235 const label communicator
240 Pout<<
"** reducing:" << Value <<
" with comm:" << communicator
245 allReduce(Value, 2, MPI_SCALAR, MPI_SUM, bop, tag, communicator);
254 const label communicator
259 Pout<<
"** reducing:" << Value <<
" with comm:" << communicator
264 vector2D twoScalars(Value, scalar(Count));
265 reduce(twoScalars, sumOp<vector2D>(), tag, communicator);
267 Value = twoScalars.x();
268 Count = twoScalars.y();
275 const sumOp<scalar>& bop,
277 const label communicator,
281 #ifdef MPIX_COMM_TYPE_SHARED 303 Pout<<
"UPstream::allocateRequest for non-blocking reduce" 304 <<
" : request:" << requestID
309 reduce(Value, bop, tag, communicator);
319 const label communicator
324 if (sendData.size() != np || recvData.size() != np)
327 <<
"Size of sendData " << sendData.size()
328 <<
" or size of recvData " << recvData.size()
329 <<
" is not equal to the number of processors in the domain " 336 recvData.deepCopy(sendData);
346 const_cast<label*>(sendData.begin()),
357 <<
"MPI_Alltoall failed for " << sendData
358 <<
" on communicator " << communicator
365 void Foam::UPstream::allocatePstreamCommunicator
367 const label parentIndex,
374 MPI_Group newGroup = MPI_GROUP_NULL;
376 MPI_Comm newComm = MPI_COMM_NULL;
382 <<
"PstreamGlobals out of sync with UPstream data. Problem." 387 if (parentIndex == -1)
394 <<
"world communicator should always be index " 411 procIDs_[index].setSize(numProcs);
412 forAll(procIDs_[index], i)
414 procIDs_[index][i] = i;
423 procIDs_[index].size(),
424 procIDs_[index].begin(),
438 myProcNo_[index] = -1;
453 <<
" when allocating communicator at " << index
454 <<
" from ranks " << procIDs_[index]
455 <<
" of parent " << parentIndex
456 <<
" cannot find my own rank" 464 void Foam::UPstream::freePstreamCommunicator(
const label communicator)
501 Pout<<
"UPstream::waitRequests : starting wait for " 503 <<
" outstanding requests starting at " << start <<
endl;
526 <<
"MPI_Waitall returned with error" <<
Foam::endl;
534 Pout<<
"UPstream::waitRequests : finished wait." <<
endl;
543 Pout<<
"UPstream::waitRequest : starting wait for request:" << i
551 <<
" outstanding send requests and you are asking for i=" << i
553 <<
"Maybe you are mixing blocking/non-blocking comms?" 567 <<
"MPI_Wait returned with error" <<
Foam::endl;
572 Pout<<
"UPstream::waitRequest : finished wait for request:" << i
582 Pout<<
"UPstream::finishedRequest : checking request:" << i
590 <<
" outstanding send requests and you are asking for i=" << i
592 <<
"Maybe you are mixing blocking/non-blocking comms?" 606 Pout<<
"UPstream::finishedRequest : finished request:" << i
634 Pout<<
"UPstream::allocateTag " << s
663 Pout<<
"UPstream::allocateTag " << s
682 Pout<<
"UPstream::freeTag " << s <<
" tag:" << tag <<
endl;
698 Pout<<
"UPstream::freeTag " << s <<
" tag:" << tag <<
endl;
DynamicList< int > freedTags_
string getEnv(const word &)
Return environment variable of given name.
static void printStack(Ostream &)
Helper function to print a stack.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Inter-processor communication reduction functions.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Various functions to wrap MPI_Allreduce.
Vector2D< scalar > vector2D
vector2D obtained from generic Vector2D
Ostream & endl(Ostream &os)
Add newline and flush stream.
static bool master(const label communicator=0)
Am I the master process.
static bool finishedRequest(const label i)
Non-blocking comms: has request i finished?
DynamicList< MPI_Request > outstandingRequests_
static label nRequests()
Get number of outstanding requests.
static label worldComm
Default communicator (all processors)
DynamicList< MPI_Group > MPIGroups_
UList< label > labelUList
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static int allocateTag(const char *)
static void freeTag(const char *, const int tag)
A class for handling words, derived from string.
void sumReduce(T &Value, label &Count, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
static void resetRequests(const label sz)
Truncate number of outstanding requests.
static bool init(int &argc, char **&argv)
Initialisation function called from main.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
static label warnComm
Debugging: warn for use of any communicator differing from warnComm.
errorManip< error > abort(error &err)
static void exit(int errnum=1)
Exit program.
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T remove()
Remove and return the top element.
static void abort()
Abort program.
static bool & parRun()
Is this a parallel run?
static label nProcs(const label communicator=0)
Number of processes in parallel run.
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
DynamicList< MPI_Comm > MPICommunicators_
static void allToAll(const labelUList &sendData, labelUList &recvData, const label communicator=0)
Exchange label with all processors (in the communicator).
static void waitRequest(const label i)
Wait until request i has finished.
static void addValidParOptions(HashTable< string > &validParOptions)
Add the valid option this type of communications library.
void allReduce(Type &Value, int count, MPI_Datatype MPIType, MPI_Op op, const BinaryOp &bop, const int tag, const label communicator)