source: XIOS/trunk/extern/src_ep/ep_reduce.cpp @ 1034

Last change on this file since 1034 was 1034, checked in by yushan, 7 years ago

adding src_ep into extern folder

File size: 19.2 KB
Line 
1/*!
2   \file ep_reduce.cpp
3   \since 2 may 2016
4
5   \brief Definitions of MPI collective function: MPI_Reduce, MPI_Allreduce
6 */
7
8#include "ep_lib.hpp"
9#include <mpi.h>
10#include "ep_declaration.hpp"
11
12using namespace std;
13
14
15namespace ep_lib {
16
17  template<typename T>
18  T max_op(T a, T b)
19  {
20    return max(a,b);
21  }
22
23  template<typename T>
24  T min_op(T a, T b)
25  {
26    return min(a,b);
27  }
28
29
30  int MPI_Reduce_local(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
31  {
32    if(datatype == MPI_INT)
33    {
34      Debug("datatype is INT\n");
35      return MPI_Reduce_local_int(sendbuf, recvbuf, count, op, comm);
36    }
37    else if(datatype == MPI_FLOAT)
38    {
39      Debug("datatype is FLOAT\n");
40      return MPI_Reduce_local_float(sendbuf, recvbuf, count, op, comm);
41    }
42    else if(datatype == MPI_DOUBLE)
43    {
44      Debug("datatype is DOUBLE\n");
45      return MPI_Reduce_local_double(sendbuf, recvbuf, count, op, comm);
46    }
47    else if(datatype == MPI_LONG)
48    {
49      Debug("datatype is DOUBLE\n");
50      return MPI_Reduce_local_long(sendbuf, recvbuf, count, op, comm);
51    }
52    else if(datatype == MPI_UNSIGNED_LONG)
53    {
54      Debug("datatype is DOUBLE\n");
55      return MPI_Reduce_local_ulong(sendbuf, recvbuf, count, op, comm);
56    }
57    else if(datatype == MPI_CHAR)
58    {
59      Debug("datatype is DOUBLE\n");
60      return MPI_Reduce_local_char(sendbuf, recvbuf, count, op, comm);
61    }
62    else
63    {
64      printf("MPI_Reduce Datatype not supported!\n");
65      exit(0);
66    }
67  }
68
69
70  int MPI_Reduce_local_int(const void *sendbuf, void *recvbuf, int count, MPI_Op op, MPI_Comm comm)
71  {
72    int my_rank = comm.ep_comm_ptr->size_rank_info[1].first;
73    int num_ep  = comm.ep_comm_ptr->size_rank_info[1].second;
74
75    int *buffer = comm.my_buffer->buf_int;
76    int *send_buf = static_cast<int*>(const_cast<void*>(sendbuf));
77    int *recv_buf = static_cast<int*>(const_cast<void*>(recvbuf));
78
79    for(int j=0; j<count; j+=BUFFER_SIZE)
80    {
81      if( 0 == my_rank )
82      {
83        #pragma omp critical (write_to_buffer)
84        copy(send_buf+j, send_buf+j + min(BUFFER_SIZE, count-j), buffer);
85        #pragma omp flush
86      }
87
88      MPI_Barrier_local(comm);
89
90      if(my_rank !=0 )
91      {
92        #pragma omp critical (write_to_buffer)
93        {
94          #pragma omp flush
95          if(op == MPI_SUM)
96          {
97            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, std::plus<int>());
98          }
99
100          else if (op == MPI_MAX)
101          {
102            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, max_op<int>);
103          }
104
105          else if (op == MPI_MIN)
106          {
107            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, min_op<int>);
108          }
109
110          else
111          {
112            printf("Supported operation: MPI_SUM, MPI_MAX, MPI_MIN\n");
113            exit(1);
114          }
115          #pragma omp flush
116        }
117      }
118
119      MPI_Barrier_local(comm);
120
121      if(my_rank == 0)
122      {
123        #pragma omp flush
124        copy(buffer, buffer+min(BUFFER_SIZE, count-j), recv_buf+j);
125      }
126      MPI_Barrier_local(comm);
127    }
128  }
129
130
131  int MPI_Reduce_local_float(const void *sendbuf, void *recvbuf, int count, MPI_Op op, MPI_Comm comm)
132  {
133    int my_rank = comm.ep_comm_ptr->size_rank_info[1].first;
134    int num_ep  = comm.ep_comm_ptr->size_rank_info[1].second;
135
136    float *buffer = comm.my_buffer->buf_float;
137    float *send_buf = static_cast<float*>(const_cast<void*>(sendbuf));
138    float *recv_buf = static_cast<float*>(const_cast<void*>(recvbuf));
139
140    for(int j=0; j<count; j+=BUFFER_SIZE)
141    {
142      if( 0 == my_rank )
143      {
144        #pragma omp critical (write_to_buffer)
145        copy(send_buf+j, send_buf+j + min(BUFFER_SIZE, count-j), buffer);
146        #pragma omp flush
147      }
148
149      MPI_Barrier_local(comm);
150
151      if(my_rank !=0 )
152      {
153        #pragma omp critical (write_to_buffer)
154        {
155          #pragma omp flush
156
157          if(op == MPI_SUM)
158          {
159            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, std::plus<float>());
160          }
161
162          else if (op == MPI_MAX)
163          {
164            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, max_op<float>);
165          }
166
167          else if (op == MPI_MIN)
168          {
169            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, min_op<float>);
170          }
171
172          else
173          {
174            printf("Supported operation: MPI_SUM, MPI_MAX, MPI_MIN\n");
175            exit(1);
176          }
177          #pragma omp flush
178        }
179      }
180
181      MPI_Barrier_local(comm);
182
183      if(my_rank == 0)
184      {
185        #pragma omp flush
186        copy(buffer, buffer+min(BUFFER_SIZE, count-j), recv_buf+j);
187      }
188      MPI_Barrier_local(comm);
189    }
190  }
191
192  int MPI_Reduce_local_double(const void *sendbuf, void *recvbuf, int count, MPI_Op op, MPI_Comm comm)
193  {
194    int my_rank = comm.ep_comm_ptr->size_rank_info[1].first;
195    int num_ep  = comm.ep_comm_ptr->size_rank_info[1].second;
196
197    double *buffer = comm.my_buffer->buf_double;
198    double *send_buf = static_cast<double*>(const_cast<void*>(sendbuf));
199    double *recv_buf = static_cast<double*>(const_cast<void*>(recvbuf));
200
201    for(int j=0; j<count; j+=BUFFER_SIZE)
202    {
203      if( 0 == my_rank )
204      {
205        #pragma omp critical (write_to_buffer)
206        copy(send_buf+j, send_buf+j + min(BUFFER_SIZE, count-j), buffer);
207        #pragma omp flush
208      }
209
210      MPI_Barrier_local(comm);
211
212      if(my_rank !=0 )
213      {
214        #pragma omp critical (write_to_buffer)
215        {
216          #pragma omp flush
217
218
219          if(op == MPI_SUM)
220          {
221            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, std::plus<double>());
222          }
223
224          else if (op == MPI_MAX)
225          {
226            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, max_op<double>);
227          }
228
229
230          else if (op == MPI_MIN)
231          {
232            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, min_op<double>);
233          }
234
235          else
236          {
237            printf("Supported operation: MPI_SUM, MPI_MAX, MPI_MIN\n");
238            exit(1);
239          }
240          #pragma omp flush
241        }
242      }
243
244      MPI_Barrier_local(comm);
245
246      if(my_rank == 0)
247      {
248        #pragma omp flush
249        copy(buffer, buffer+min(BUFFER_SIZE, count-j), recv_buf+j);
250      }
251      MPI_Barrier_local(comm);
252    }
253  }
254
255  int MPI_Reduce_local_long(const void *sendbuf, void *recvbuf, int count, MPI_Op op, MPI_Comm comm)
256  {
257    int my_rank = comm.ep_comm_ptr->size_rank_info[1].first;
258    int num_ep  = comm.ep_comm_ptr->size_rank_info[1].second;
259
260    long *buffer = comm.my_buffer->buf_long;
261    long *send_buf = static_cast<long*>(const_cast<void*>(sendbuf));
262    long *recv_buf = static_cast<long*>(const_cast<void*>(recvbuf));
263
264    for(int j=0; j<count; j+=BUFFER_SIZE)
265    {
266      if( 0 == my_rank )
267      {
268        #pragma omp critical (write_to_buffer)
269        copy(send_buf+j, send_buf+j + min(BUFFER_SIZE, count-j), buffer);
270        #pragma omp flush
271      }
272
273      MPI_Barrier_local(comm);
274
275      if(my_rank !=0 )
276      {
277        #pragma omp critical (write_to_buffer)
278        {
279          #pragma omp flush
280
281
282          if(op == MPI_SUM)
283          {
284            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, std::plus<long>());
285          }
286
287          else if (op == MPI_MAX)
288          {
289            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, max_op<long>);
290          }
291
292
293          else if (op == MPI_MIN)
294          {
295            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, min_op<long>);
296          }
297
298          else
299          {
300            printf("Supported operation: MPI_SUM, MPI_MAX, MPI_MIN\n");
301            exit(1);
302          }
303          #pragma omp flush
304        }
305      }
306
307      MPI_Barrier_local(comm);
308
309      if(my_rank == 0)
310      {
311        #pragma omp flush
312        copy(buffer, buffer+min(BUFFER_SIZE, count-j), recv_buf+j);
313      }
314      MPI_Barrier_local(comm);
315    }
316  }
317
318  int MPI_Reduce_local_ulong(const void *sendbuf, void *recvbuf, int count, MPI_Op op, MPI_Comm comm)
319  {
320    int my_rank = comm.ep_comm_ptr->size_rank_info[1].first;
321    int num_ep  = comm.ep_comm_ptr->size_rank_info[1].second;
322
323    unsigned long *buffer = comm.my_buffer->buf_ulong;
324    unsigned long *send_buf = static_cast<unsigned long*>(const_cast<void*>(sendbuf));
325    unsigned long *recv_buf = static_cast<unsigned long*>(const_cast<void*>(recvbuf));
326
327    for(int j=0; j<count; j+=BUFFER_SIZE)
328    {
329      if( 0 == my_rank )
330      {
331        #pragma omp critical (write_to_buffer)
332        copy(send_buf+j, send_buf+j + min(BUFFER_SIZE, count-j), buffer);
333        #pragma omp flush
334      }
335
336      MPI_Barrier_local(comm);
337
338      if(my_rank !=0 )
339      {
340        #pragma omp critical (write_to_buffer)
341        {
342          #pragma omp flush
343
344
345          if(op == MPI_SUM)
346          {
347            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, std::plus<unsigned long>());
348          }
349
350          else if (op == MPI_MAX)
351          {
352            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, max_op<unsigned long>);
353          }
354
355
356          else if (op == MPI_MIN)
357          {
358            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, min_op<unsigned long>);
359          }
360
361          else
362          {
363            printf("Supported operation: MPI_SUM, MPI_MAX, MPI_MIN\n");
364            exit(1);
365          }
366          #pragma omp flush
367        }
368      }
369
370      MPI_Barrier_local(comm);
371
372      if(my_rank == 0)
373      {
374        #pragma omp flush
375        copy(buffer, buffer+min(BUFFER_SIZE, count-j), recv_buf+j);
376      }
377      MPI_Barrier_local(comm);
378    }
379  }
380
381  int MPI_Reduce_local_char(const void *sendbuf, void *recvbuf, int count, MPI_Op op, MPI_Comm comm)
382  {
383    int my_rank = comm.ep_comm_ptr->size_rank_info[1].first;
384    int num_ep  = comm.ep_comm_ptr->size_rank_info[1].second;
385
386    char *buffer = comm.my_buffer->buf_char;
387    char *send_buf = static_cast<char*>(const_cast<void*>(sendbuf));
388    char *recv_buf = static_cast<char*>(const_cast<void*>(recvbuf));
389
390    for(int j=0; j<count; j+=BUFFER_SIZE)
391    {
392      if( 0 == my_rank )
393      {
394        #pragma omp critical (write_to_buffer)
395        copy(send_buf+j, send_buf+j + min(BUFFER_SIZE, count-j), buffer);
396        #pragma omp flush
397      }
398
399      MPI_Barrier_local(comm);
400
401      if(my_rank !=0 )
402      {
403        #pragma omp critical (write_to_buffer)
404        {
405          #pragma omp flush
406
407
408          if(op == MPI_SUM)
409          {
410            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, std::plus<char>());
411          }
412
413          else if (op == MPI_MAX)
414          {
415            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, max_op<char>);
416          }
417
418
419          else if (op == MPI_MIN)
420          {
421            transform(buffer, buffer+min(BUFFER_SIZE, count-j), send_buf+j, buffer, min_op<char>);
422          }
423
424          else
425          {
426            printf("Supported operation: MPI_SUM, MPI_MAX, MPI_MIN\n");
427            exit(1);
428          }
429          #pragma omp flush
430        }
431      }
432
433      MPI_Barrier_local(comm);
434
435      if(my_rank == 0)
436      {
437        #pragma omp flush
438        copy(buffer, buffer+min(BUFFER_SIZE, count-j), recv_buf+j);
439      }
440      MPI_Barrier_local(comm);
441    }
442  }
443
444
445  int MPI_Reduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
446  {
447    if(!comm.is_ep && comm.mpi_comm)
448    {
449      #ifdef _serialized
450      #pragma omp critical (_mpi_call)
451      #endif // _serialized
452      ::MPI_Reduce(sendbuf, recvbuf, count, static_cast< ::MPI_Datatype>(datatype), static_cast< ::MPI_Op>(op), root,
453                   static_cast< ::MPI_Comm>(comm.mpi_comm));
454      return 0;
455    }
456
457
458    if(!comm.mpi_comm) return 0;
459
460    int root_mpi_rank = comm.rank_map->at(root).second;
461    int root_ep_loc = comm.rank_map->at(root).first;
462
463    int ep_rank, ep_rank_loc, mpi_rank;
464    int ep_size, num_ep, mpi_size;
465
466    ep_rank = comm.ep_comm_ptr->size_rank_info[0].first;
467    ep_rank_loc = comm.ep_comm_ptr->size_rank_info[1].first;
468    mpi_rank = comm.ep_comm_ptr->size_rank_info[2].first;
469    ep_size = comm.ep_comm_ptr->size_rank_info[0].second;
470    num_ep = comm.ep_comm_ptr->size_rank_info[1].second;
471    mpi_size = comm.ep_comm_ptr->size_rank_info[2].second;
472
473
474    ::MPI_Aint recvsize, lb;
475    #ifdef _serialized
476    #pragma omp critical (_mpi_call)
477    #endif // _serialized
478    ::MPI_Type_get_extent(static_cast< ::MPI_Datatype>(datatype), &lb, &recvsize);
479
480    void *local_recvbuf;
481    if(ep_rank_loc==0)
482    {
483      local_recvbuf = new void*[recvsize*count];
484    }
485
486    MPI_Reduce_local(sendbuf, local_recvbuf, count, datatype, op, comm);
487
488
489    if(ep_rank_loc==0)
490    {
491      #ifdef _serialized
492      #pragma omp critical (_mpi_call)
493      #endif // _serialized
494      ::MPI_Reduce(local_recvbuf, recvbuf, count, static_cast< ::MPI_Datatype>(datatype), static_cast< ::MPI_Op>(op), root_mpi_rank, static_cast< ::MPI_Comm>(comm.mpi_comm));
495    }
496
497    if(root_ep_loc != 0 && mpi_rank == root_mpi_rank) // root is not master, master send to root and root receive from master
498    {
499      innode_memcpy(0, recvbuf, root_ep_loc, recvbuf, count, datatype, comm);
500    }
501
502    if(ep_rank_loc==0)
503    {
504      if(datatype == MPI_INT) delete[] static_cast<int*>(local_recvbuf);
505      else if(datatype == MPI_FLOAT) delete[] static_cast<float*>(local_recvbuf);
506      else if(datatype == MPI_DOUBLE) delete[] static_cast<double*>(local_recvbuf);
507      else if(datatype == MPI_LONG) delete[] static_cast<long*>(local_recvbuf);
508      else if(datatype == MPI_UNSIGNED_LONG) delete[] static_cast<unsigned long*>(local_recvbuf);
509      else delete[] static_cast<char*>(local_recvbuf);
510    }
511
512    Message_Check(comm);
513
514    return 0;
515  }
516
517
518
519
520  int MPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
521  {
522    if(!comm.is_ep && comm.mpi_comm)
523    {
524      #ifdef _serialized
525      #pragma omp critical (_mpi_call)
526      #endif // _serialized
527      ::MPI_Allreduce(sendbuf, recvbuf, count, static_cast< ::MPI_Datatype>(datatype), static_cast< ::MPI_Op>(op),
528                      static_cast< ::MPI_Comm>(comm.mpi_comm));
529      return 0;
530    }
531
532    if(!comm.mpi_comm) return 0;
533
534
535    int ep_rank, ep_rank_loc, mpi_rank;
536    int ep_size, num_ep, mpi_size;
537
538    ep_rank = comm.ep_comm_ptr->size_rank_info[0].first;
539    ep_rank_loc = comm.ep_comm_ptr->size_rank_info[1].first;
540    mpi_rank = comm.ep_comm_ptr->size_rank_info[2].first;
541    ep_size = comm.ep_comm_ptr->size_rank_info[0].second;
542    num_ep = comm.ep_comm_ptr->size_rank_info[1].second;
543    mpi_size = comm.ep_comm_ptr->size_rank_info[2].second;
544
545
546    ::MPI_Aint recvsize, lb;
547    #ifdef _serialized
548    #pragma omp critical (_mpi_call)
549    #endif // _serialized
550    ::MPI_Type_get_extent(static_cast< ::MPI_Datatype>(datatype), &lb, &recvsize);
551
552    void *local_recvbuf;
553    if(ep_rank_loc==0)
554    {
555      local_recvbuf = new void*[recvsize*count];
556    }
557
558    MPI_Reduce_local(sendbuf, local_recvbuf, count, datatype, op, comm);
559
560
561    if(ep_rank_loc==0)
562    {
563      #ifdef _serialized
564      #pragma omp critical (_mpi_call)
565      #endif // _serialized
566      ::MPI_Allreduce(local_recvbuf, recvbuf, count, static_cast< ::MPI_Datatype>(datatype), static_cast< ::MPI_Op>(op), static_cast< ::MPI_Comm>(comm.mpi_comm));
567    }
568
569    MPI_Bcast_local(recvbuf, count, datatype, comm);
570
571    if(ep_rank_loc==0)
572    {
573      if(datatype == MPI_INT) delete[] static_cast<int*>(local_recvbuf);
574      else if(datatype == MPI_FLOAT) delete[] static_cast<float*>(local_recvbuf);
575      else if(datatype == MPI_DOUBLE) delete[] static_cast<double*>(local_recvbuf);
576      else if(datatype == MPI_LONG) delete[] static_cast<long*>(local_recvbuf);
577      else if(datatype == MPI_UNSIGNED_LONG) delete[] static_cast<unsigned long*>(local_recvbuf);
578      else delete[] static_cast<char*>(local_recvbuf);
579    }
580
581    Message_Check(comm);
582
583    return 0;
584  }
585
586
587  int MPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int recvcounts[], MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
588  {
589
590    if(!comm.is_ep && comm.mpi_comm)
591    {
592      #ifdef _serialized
593      #pragma omp critical (_mpi_call)
594      #endif // _serialized
595      ::MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, static_cast< ::MPI_Datatype>(datatype), static_cast< ::MPI_Op>(op),
596                           static_cast< ::MPI_Comm>(comm.mpi_comm));
597      return 0;
598    }
599
600    if(!comm.mpi_comm) return 0;
601
602    int ep_rank, ep_rank_loc, mpi_rank;
603    int ep_size, num_ep, mpi_size;
604
605    ep_rank = comm.ep_comm_ptr->size_rank_info[0].first;
606    ep_rank_loc = comm.ep_comm_ptr->size_rank_info[1].first;
607    mpi_rank = comm.ep_comm_ptr->size_rank_info[2].first;
608    ep_size = comm.ep_comm_ptr->size_rank_info[0].second;
609    num_ep = comm.ep_comm_ptr->size_rank_info[1].second;
610    mpi_size = comm.ep_comm_ptr->size_rank_info[2].second;
611
612    void* local_buf;
613    void* local_buf2;
614    int local_buf_size = accumulate(recvcounts, recvcounts+ep_size, 0);
615    int local_buf2_size = accumulate(recvcounts+ep_rank-ep_rank_loc, recvcounts+ep_rank-ep_rank_loc+num_ep, 0);
616
617    ::MPI_Aint datasize, lb;
618    #ifdef _serialized
619    #pragma omp critical (_mpi_call)
620    #endif // _serialized
621    ::MPI_Type_get_extent(static_cast< ::MPI_Datatype>(datatype), &lb, &datasize);
622
623    if(ep_rank_loc == 0)
624    {
625      local_buf = new void*[local_buf_size*datasize];
626      local_buf2 = new void*[local_buf2_size*datasize];
627    }
628    MPI_Reduce_local(sendbuf, local_buf, local_buf_size, MPI_INT, op, comm);
629
630
631    if(ep_rank_loc == 0)
632    {
633      int local_recvcnt[mpi_size];
634      for(int i=0; i<mpi_size; i++)
635      {
636        local_recvcnt[i] = accumulate(recvcounts+ep_rank, recvcounts+ep_rank+num_ep, 0);
637      }
638      #ifdef _serialized
639      #pragma omp critical (_mpi_call)
640      #endif // _serialized
641      ::MPI_Reduce_scatter(local_buf, local_buf2, local_recvcnt, static_cast< ::MPI_Datatype>(datatype),
642                         static_cast< ::MPI_Op>(op), static_cast< ::MPI_Comm>(comm.mpi_comm));
643    }
644
645
646    int displs[num_ep];
647    displs[0] = 0;
648    for(int i=1; i<num_ep; i++)
649    {
650      displs[i] = displs[i-1] + recvcounts[ep_rank-ep_rank_loc+i-1];
651    }
652
653    MPI_Scatterv_local(local_buf2, recvcounts+ep_rank-ep_rank_loc, displs, datatype, recvbuf, comm);
654
655    if(ep_rank_loc == 0)
656    {
657      if(datatype == MPI_INT)
658      {
659        delete[] static_cast<int*>(local_buf);
660        delete[] static_cast<int*>(local_buf2);
661      }
662      else if(datatype == MPI_FLOAT)
663      {
664        delete[] static_cast<float*>(local_buf);
665        delete[] static_cast<float*>(local_buf2);
666      }
667      else if(datatype == MPI_DOUBLE)
668      {
669        delete[] static_cast<double*>(local_buf);
670        delete[] static_cast<double*>(local_buf2);
671      }
672      else if(datatype == MPI_LONG)
673      {
674        delete[] static_cast<long*>(local_buf);
675        delete[] static_cast<long*>(local_buf2);
676      }
677      else if(datatype == MPI_UNSIGNED_LONG)
678      {
679        delete[] static_cast<unsigned long*>(local_buf);
680        delete[] static_cast<unsigned long*>(local_buf2);
681      }
682      else // if(datatype == MPI_DOUBLE)
683      {
684        delete[] static_cast<char*>(local_buf);
685        delete[] static_cast<char*>(local_buf2);
686      }
687    }
688
689    Message_Check(comm);
690    return 0;
691  }
692}
693
Note: See TracBrowser for help on using the repository browser.