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