]> git.donarmstrong.com Git - libparallel-mpi-simple-perl.git/blob - Simple.xs
require libopenpmi-dev
[libparallel-mpi-simple-perl.git] / Simple.xs
1 #include "EXTERN.h"
2 #include "perl.h"
3 #include "XSUB.h"
4
5 /* The Inline::C source for this is available in the __DATA__ section
6    of Simple.pm */
7
8 #include <mpi.h> 
9 #define GATHER_TAG 11001 /* used to be unlikely to upset other sends */
10
11 /*
12   root process first broadcasts length of stored data then broadcasts
13   the data.  Non-root processes receive length (via bcast), allocate
14   space to take incomming data from root
15
16   Both root and non-root processes then create and return a new scalar
17   with contents identical to those root started with.
18 */
19
20 SV* _Bcast (SV* data, int root, SV* comm) {
21   int buf_len[1];
22   int rank;
23   SV* rval;
24   MPI_Comm_rank((MPI_Comm)SvIVX(comm), &rank);
25   if (rank == root) {
26     buf_len[0] = sv_len(data);
27     MPI_Bcast(buf_len, 1, MPI_INT, root, (MPI_Comm)SvIVX(comm));
28     MPI_Bcast(SvPVX(data), buf_len[0], MPI_CHAR, root, (MPI_Comm)SvIVX(comm));
29     rval = newSVpvn(SvPVX(data), buf_len[0]);
30   }
31   else {
32     char *recv_buf;
33     MPI_Bcast(buf_len, 1, MPI_INT, root, (MPI_Comm)SvIVX(comm));
34     recv_buf = (char*)malloc((buf_len[0]+1)*sizeof(char));
35     if (recv_buf == NULL) croak("Allocation error in _Bcast");
36     MPI_Bcast(recv_buf, buf_len[0], MPI_CHAR, root, (MPI_Comm)SvIVX(comm));
37     rval = newSVpvn(recv_buf, buf_len[0]);
38     free(recv_buf);
39   }
40   return rval;
41 }
42
43 /*
44   Finds length of data in stor_ref, sends this to receiver, then
45   sends actual data, uses same tag for each message.
46 */
47
48 int _Send(SV* stor_ref, int dest, int tag, SV*comm) {
49   int str_len[1];
50   str_len[0] = sv_len(stor_ref);
51   MPI_Send(str_len, 1, MPI_INT, dest, tag, (MPI_Comm)SvIVX(comm));
52   MPI_Send(SvPVX(stor_ref), sv_len(stor_ref),MPI_CHAR,
53            dest, tag, (MPI_Comm)SvIVX(comm));
54   return 0;
55 }
56
57 /*
58   Receives int for length of data it should then expect, allocates space
59   then receives data into that space.  Creates a new SV and returns it.
60 */
61
62 SV* _Recv (int source, int tag, SV*comm, SV*status) {
63   MPI_Status tstatus;
64   SV* rval;
65   int len_buf[1];
66   char *recv_buf;
67
68   MPI_Recv(len_buf, 1, MPI_INT, source, tag, (MPI_Comm)SvIVX(comm), &tstatus);
69   recv_buf = (char*)malloc((len_buf[0]+1)*sizeof(char));
70   if (recv_buf == NULL) croak("Allocation error in _Recv");
71   MPI_Recv(recv_buf, len_buf[0], MPI_CHAR, source, tag,
72             (MPI_Comm)SvIVX(comm), &tstatus);
73   rval = newSVpvn(recv_buf, len_buf[0]);
74   sv_setiv(status, tstatus.MPI_SOURCE);
75   free(recv_buf);
76   return rval;
77 }
78
79 /* Calls MPI_Init with dummy arguments, a bit dodgy but sort of ok */
80 int Init () {
81   MPI_Init(&PL_origargc, &PL_origargv);
82 }
83
84 /* Returns rank of process within comm */
85 int _Comm_rank (SV* comm) {
86   int trank;
87   MPI_Comm_rank((MPI_Comm)SvIVX(comm),&trank);
88   return trank;
89 }
90
91 /* returns total number of processes within comm */
92 int _Comm_size (SV* comm) {
93   int tsize;
94   MPI_Comm_size((MPI_Comm)SvIVX(comm), &tsize);
95   return tsize;
96 }
97
98 /* returns SV whose IV slot is a cast pointer to the MPI_COMM_WORLD object */
99 SV* COMM_WORLD () {
100   return newSViv((IV)MPI_COMM_WORLD);
101 }
102
103 /* returns SV whose IV slot is a cast pointer to the MPI_ANY_SOURCE value */
104 SV* ANY_SOURCE () {
105   return newSViv((IV)MPI_ANY_SOURCE);
106 }
107
108 /* calls MPI_Barrier for comm */
109 int Barrier (SV*comm) {
110   MPI_Barrier((MPI_Comm)SvIVX(comm));
111 }
112
113 /* ends MPI participation */
114 int Finalize () {
115   MPI_Finalize();
116 }
117
118 /*
119   If non-root:  participates in Gather so that root finds length of data
120                 to expect from this process.  Then send (using MPI_Send)
121                 data to root.
122
123   If root: receives array of ints detailing length of scalars held by
124    other processes, then receives from each in turn (using MPI_Recv)
125    returns an array ref to root process only.
126   
127  */
128 SV* _Gather (SV* data, int root, SV* comm) {
129   int rank, size, *buf_lens, i, max;
130   char* recv_buf;
131   int my_buf[1];
132   AV* ret_arr;
133   MPI_Status tstatus;
134
135   /* find out how long data is */
136   ret_arr = av_make(0,(SV**)NULL);
137   my_buf[0] = sv_len(data);
138   if (_Comm_rank(comm) == root) {
139     MPI_Comm_size((MPI_Comm)SvIVX(comm), &size);
140     buf_lens = malloc(size*sizeof(int));
141     if (buf_lens == NULL) croak("Allocation error (lens) in _Gather");
142     /* gather all scalar length data */
143     MPI_Gather(my_buf, 1, MPI_INT, buf_lens, 1,
144                MPI_INT, root, (MPI_Comm)SvIVX(comm));
145     max = 0; // keep buffer allocation calls to minimum
146     for (i=0;i<size;i++) {
147       max = max < buf_lens[i] ? buf_lens[i] : max;
148     }
149     recv_buf = malloc(max * sizeof(char));
150     if (recv_buf == NULL) croak("Allocation error (recv) in _Gather");
151     for (i=0;i<size;i++) {
152       if (i == root) {
153         av_push(ret_arr, data);
154         continue; /* me, no point sending */
155       }
156       MPI_Recv(recv_buf, buf_lens[i], MPI_CHAR, i, GATHER_TAG,
157                (MPI_Comm)SvIVX(comm), &tstatus );
158       av_push(ret_arr, sv_2mortal( newSVpvn(recv_buf, buf_lens[i]) ) );
159     }
160     free(recv_buf);
161     free(buf_lens);
162   }
163   else {
164     /* send out how long my scalar data is */ 
165       MPI_Gather(my_buf, 1, MPI_INT, buf_lens, 1,
166                MPI_INT, root, (MPI_Comm)SvIVX(comm) );
167     /* send out my scalar data as normal send with tag of ???? */
168       MPI_Send(SvPVX(data), my_buf[0], MPI_CHAR,
169                root, GATHER_TAG,(MPI_Comm)SvIVX(comm));
170   }
171
172   return newRV_inc((SV*)ret_arr);
173 }
174
175 /* compares two communicators, translates MPI constants into something I
176    can use as constants in the module interface */
177 int _Comm_compare(SV* comm1, SV* comm2) {
178     int result = 0;
179     MPI_Comm_compare((MPI_Comm)SvIVX(comm1), (MPI_Comm)SvIVX(comm2), &result);
180     switch (result) {
181         case MPI_IDENT:
182                        return(1);
183         case MPI_CONGRUENT:
184                        return(2);
185         case MPI_SIMILAR:
186                        return(3);
187         case MPI_UNEQUAL:
188                        return(0);
189         default:
190                        return(0);
191     }
192 }
193
194 /* frees a communicator, once all pending communication has taken place */
195 void _Comm_free (SV* comm) {
196     MPI_Comm_free((MPI_Comm*)&SvIVX(comm));
197     if ((MPI_Comm)SvIVX(comm) != MPI_COMM_NULL)
198         croak("Communicator not freed properly\n");
199 }
200
201 SV* _Comm_dup (SV*comm) {
202     MPI_Comm newcomm;
203     MPI_Comm_dup((MPI_Comm)SvIVX(comm), &newcomm);
204     return newSViv((IV)newcomm);
205 }
206
207 SV* _Comm_split (SV* comm, int colour, int key) {
208     MPI_Comm newcomm;
209     int realcolour;
210     MPI_Comm_split((MPI_Comm)SvIVX(comm),
211                     (colour < 0 ? MPI_UNDEFINED : colour),
212                     key, &newcomm);
213     return newSViv((IV)newcomm);
214 }
215
216
217
218
219 MODULE = Parallel::MPI::Simple  PACKAGE = Parallel::MPI::Simple 
220
221 PROTOTYPES: DISABLE
222
223 SV *
224 _Bcast (data, root, comm)
225         SV *    data
226         int     root
227         SV *    comm
228
229 int
230 _Send (stor_ref, dest, tag, comm)
231         SV *    stor_ref
232         int     dest
233         int     tag
234         SV *    comm
235
236 SV *
237 _Recv (source, tag, comm, status)
238         int     source
239         int     tag
240         SV *    comm
241         SV *    status
242
243 int
244 Init ()
245
246 int
247 _Comm_rank (comm)
248         SV *    comm
249
250 int
251 _Comm_size (comm)
252         SV *    comm
253
254 SV *
255 COMM_WORLD ()
256
257 SV *
258 ANY_SOURCE ()
259
260 int
261 Barrier (comm)
262         SV *    comm
263
264 int
265 Finalize ()
266
267 SV *
268 _Gather (data, root, comm)
269         SV *    data
270         int     root
271         SV *    comm
272
273 int
274 _Comm_compare (comm1, comm2)
275         SV *    comm1
276         SV *    comm2
277
278 void
279 _Comm_free (comm)
280         SV *    comm
281         PREINIT:
282         I32* temp;
283         PPCODE:
284         temp = PL_markstack_ptr++;
285         _Comm_free(comm);
286         if (PL_markstack_ptr != temp) {
287           /* truly void, because dXSARGS not invoked */
288           PL_markstack_ptr = temp;
289           XSRETURN_EMPTY; /* return empty stack */
290         }
291         /* must have used dXSARGS; list context implied */
292         return; /* assume stack size is correct */
293
294 SV *
295 _Comm_dup (comm)
296         SV *    comm
297
298 SV *
299 _Comm_split (comm, colour, key)
300         SV *    comm
301         int     colour
302         int     key
303