Actual source code: mpiutils.h
1: #pragma once
3: #include <petscsys.h>
5: PETSC_EXTERN PetscErrorCode PetscGatherNumberOfMessages_Private(MPI_Comm, const PetscMPIInt[], const PetscInt[], PetscMPIInt *);
6: PETSC_EXTERN PetscErrorCode PetscGatherMessageLengths_Private(MPI_Comm, PetscMPIInt, PetscMPIInt, const PetscInt[], PetscMPIInt **, PetscInt **);
8: #if !defined(PETSC_HAVE_MPI_LARGE_COUNT) /* No matter PetscInt is 32-bit or 64-bit, without MPI large count we always do casting before MPI calls */
9: /* Cast PetscInt <a> to PetscMPIInt <b>, where <a> is likely used for the 'count' argument in MPI routines.
10: It is similar to PetscMPIIntCast() except that here it returns an MPI error code.
11: */
12: static inline PetscMPIInt PetscMPIIntCast_Internal(PetscInt a, PetscMPIInt *b)
13: {
14: *b = (PetscMPIInt)(a);
15: if (PetscDefined(USE_64BIT_INDICIES) && PetscUnlikely(a > PETSC_MPI_INT_MAX)) return MPI_ERR_COUNT;
16: return MPI_SUCCESS;
17: }
19: static inline PetscMPIInt MPIU_Send(const void *buf, PetscInt count, MPI_Datatype datatype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm)
20: {
21: PetscMPIInt count2;
23: PetscFunctionBegin;
24: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
25: PetscCallMPI(MPI_Send(buf, count2, datatype, dest, tag, comm));
26: PetscFunctionReturn(MPI_SUCCESS);
27: }
29: static inline PetscMPIInt MPIU_Send_init(const void *buf, PetscInt count, MPI_Datatype datatype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
30: {
31: PetscMPIInt count2;
33: PetscFunctionBegin;
34: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
35: PetscCallMPI(MPI_Send_init(buf, count2, datatype, dest, tag, comm, request));
36: PetscFunctionReturn(MPI_SUCCESS);
37: }
39: static inline PetscMPIInt MPIU_Isend(const void *buf, PetscInt count, MPI_Datatype datatype, PetscMPIInt dest, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
40: {
41: PetscMPIInt count2;
43: PetscFunctionBegin;
44: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
45: PetscCallMPI(MPI_Isend(buf, count2, datatype, dest, tag, comm, request));
46: PetscFunctionReturn(MPI_SUCCESS);
47: }
49: static inline PetscMPIInt MPIU_Recv(void *buf, PetscInt count, MPI_Datatype datatype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Status *status)
50: {
51: PetscMPIInt count2;
53: PetscFunctionBegin;
54: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
55: PetscCallMPI(MPI_Recv(buf, count2, datatype, source, tag, comm, status));
56: PetscFunctionReturn(MPI_SUCCESS);
57: }
59: static inline PetscMPIInt MPIU_Recv_init(void *buf, PetscInt count, MPI_Datatype datatype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
60: {
61: PetscMPIInt count2;
63: PetscFunctionBegin;
64: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
65: PetscCallMPI(MPI_Recv_init(buf, count2, datatype, source, tag, comm, request));
66: PetscFunctionReturn(MPI_SUCCESS);
67: }
69: static inline PetscMPIInt MPIU_Irecv(void *buf, PetscInt count, MPI_Datatype datatype, PetscMPIInt source, PetscMPIInt tag, MPI_Comm comm, MPI_Request *request)
70: {
71: PetscMPIInt count2;
73: PetscFunctionBegin;
74: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
75: PetscCallMPI(MPI_Irecv(buf, count2, datatype, source, tag, comm, request));
76: PetscFunctionReturn(MPI_SUCCESS);
77: }
78: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
79: static inline PetscMPIInt MPIU_Reduce_local(const void *inbuf, void *inoutbuf, PetscInt count, MPI_Datatype datatype, MPI_Op op)
80: {
81: PetscMPIInt count2;
83: PetscFunctionBegin;
84: PetscCallMPI(PetscMPIIntCast_Internal(count, &count2));
85: PetscCallMPI(MPI_Reduce_local(inbuf, inoutbuf, count, datatype, op));
86: PetscFunctionReturn(MPI_SUCCESS);
87: }
88: #endif
90: #elif defined(PETSC_USE_64BIT_INDICES)
91: #define MPIU_Send(buf, count, datatype, dest, tag, comm) MPI_Send_c(buf, count, datatype, dest, tag, comm)
92: #define MPIU_Send_init(buf, count, datatype, dest, tag, comm, request) MPI_Send_init_c(buf, count, datatype, dest, tag, comm, request)
93: #define MPIU_Isend(buf, count, datatype, dest, tag, comm, request) MPI_Isend_c(buf, count, datatype, dest, tag, comm, request)
94: #define MPIU_Recv(buf, count, datatype, source, tag, comm, status) MPI_Recv_c(buf, count, datatype, source, tag, comm, status)
95: #define MPIU_Recv_init(buf, count, datatype, source, tag, comm, request) MPI_Recv_init_c(buf, count, datatype, source, tag, comm, request)
96: #define MPIU_Irecv(buf, count, datatype, source, tag, comm, request) MPI_Irecv_c(buf, count, datatype, source, tag, comm, request)
97: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
98: #define MPIU_Reduce_local(inbuf, inoutbuf, count, datatype, op) MPI_Reduce_local_c(inbuf, inoutbuf, count, datatype, op)
99: #endif
100: #else
101: #define MPIU_Send(buf, count, datatype, dest, tag, comm) MPI_Send(buf, count, datatype, dest, tag, comm)
102: #define MPIU_Send_init(buf, count, datatype, dest, tag, comm, request) MPI_Send_init(buf, count, datatype, dest, tag, comm, request)
103: #define MPIU_Isend(buf, count, datatype, dest, tag, comm, request) MPI_Isend(buf, count, datatype, dest, tag, comm, request)
104: #define MPIU_Recv(buf, count, datatype, source, tag, comm, status) MPI_Recv(buf, count, datatype, source, tag, comm, status)
105: #define MPIU_Recv_init(buf, count, datatype, source, tag, comm, request) MPI_Recv_init(buf, count, datatype, source, tag, comm, request)
106: #define MPIU_Irecv(buf, count, datatype, source, tag, comm, request) MPI_Irecv(buf, count, datatype, source, tag, comm, request)
107: #if defined(PETSC_HAVE_MPI_REDUCE_LOCAL)
108: #define MPIU_Reduce_local(inbuf, inoutbuf, count, datatype, op) MPI_Reduce_local(inbuf, inoutbuf, count, datatype, op)
109: #endif
110: #endif
112: /* These APIs use arrays of MPI_Count/MPI_Aint */
113: #if defined(PETSC_HAVE_MPI_LARGE_COUNT) && defined(PETSC_USE_64BIT_INDICES)
114: #define MPIU_Neighbor_alltoallv(a, b, c, d, e, f, g, h, i) MPI_Neighbor_alltoallv_c(a, b, c, d, e, f, g, h, i)
115: #define MPIU_Neighbor_alltoallv_init(a, b, c, d, e, f, g, h, i, j, k) MPI_Neighbor_alltoallv_init_c(a, b, c, d, e, f, g, h, i, j, k)
116: #define MPIU_Ineighbor_alltoallv(a, b, c, d, e, f, g, h, i, j) MPI_Ineighbor_alltoallv_c(a, b, c, d, e, f, g, h, i, j)
117: #else
118: #define MPIU_Neighbor_alltoallv(a, b, c, d, e, f, g, h, i) MPI_Neighbor_alltoallv(a, b, c, d, e, f, g, h, i)
119: #define MPIU_Neighbor_alltoallv_init(a, b, c, d, e, f, g, h, i, j, k) MPI_Neighbor_alltoallv_init(a, b, c, d, e, f, g, h, i, j, k)
120: #define MPIU_Ineighbor_alltoallv(a, b, c, d, e, f, g, h, i, j) MPI_Ineighbor_alltoallv(a, b, c, d, e, f, g, h, i, j)
121: #endif