Changeset 17627 for issm/trunk-jpl/src/c/toolkits/mpi/issmmpi.h
- Timestamp:
- 04/02/14 15:36:30 (11 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/toolkits/mpi/issmmpi.h
r16144 r17627 12 12 #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!" 13 13 #endif 14 15 14 #include <cstddef> 15 #include <cassert> 16 #include "../../shared/Numerics/types.h" 16 17 17 18 #ifdef _HAVE_MPI_ … … 83 84 /*}}}*/ 84 85 # endif 86 87 template <class T> ISSM_MPI_Datatype TypeToMPIType(){assert(false);}; 88 template <> inline ISSM_MPI_Datatype TypeToMPIType<IssmDouble>(){return ISSM_MPI_DOUBLE;}; 89 #if _HAVE_ADOLC_ 90 template <> inline ISSM_MPI_Datatype TypeToMPIType<IssmPDouble>(){return ISSM_MPI_PDOUBLE;}; 91 #endif 92 template <> inline ISSM_MPI_Datatype TypeToMPIType<int>(){return ISSM_MPI_INT;}; 93 template <> inline ISSM_MPI_Datatype TypeToMPIType<char>(){return ISSM_MPI_CHAR;}; 94 template <class T> ISSM_MPI_Datatype ISSM_MPI_Reduce(T *sendbuf,T *recvbuf,int count,ISSM_MPI_Op op, int root, ISSM_MPI_Comm comm){ /*{{{*/ 95 96 /*Get MPI type*/ 97 ISSM_MPI_Datatype datatype = TypeToMPIType<T>(); 98 99 int rc=0; 100 #ifdef _HAVE_MPI_ 101 # ifdef _HAVE_AMPI_ 102 rc=AMPI_Reduce(sendbuf, 103 recvbuf, 104 count, 105 datatype, 106 op, 107 root, 108 comm); 109 # else 110 rc=MPI_Reduce(sendbuf, 111 recvbuf, 112 count, 113 datatype, 114 op, 115 root, 116 comm); 117 # endif 118 #else 119 # ifdef _HAVE_ADOLC_ 120 if (datatype==ISSM_MPI_DOUBLE) { 121 IssmDouble* activeSendBuf=(IssmDouble*)sendbuf; 122 IssmDouble* activeRecvBuf=(IssmDouble*)recvbuf; 123 for(int i=0;i<count;++i) activeRecvBuf[i]=activeSendBuf[i]; 124 } 125 else 126 # endif 127 memcpy(recvbuf,sendbuf,sizeHelper(datatype)*count); 128 #endif 129 return rc; 130 }/*}}}*/ 131 template <class T> int ISSM_MPI_Bcast(T *buffer, int count,int root, ISSM_MPI_Comm comm){ /*{{{*/ 132 133 int rc=0; 134 135 /*Get MPI type*/ 136 ISSM_MPI_Datatype datatype = TypeToMPIType<T>(); 137 138 #ifdef _HAVE_MPI_ 139 # ifdef _HAVE_AMPI_ 140 rc=AMPI_Bcast(buffer, 141 count, 142 datatype, 143 root, 144 comm); 145 # else 146 rc=MPI_Bcast(buffer, 147 count, 148 datatype, 149 root, 150 comm); 151 # endif 152 #else 153 // nothing to be done here 154 #endif 155 return rc; 156 }/*}}}*/ 157 85 158 /* interfaces {{{*/ 86 159 int ISSM_MPI_Allgather(void *sendbuf, int sendcount, ISSM_MPI_Datatype sendtype, void *recvbuf, int recvcount, ISSM_MPI_Datatype recvtype, ISSM_MPI_Comm comm);
Note:
See TracChangeset
for help on using the changeset viewer.