MIDAPACK - MIcrowave Data Analysis PACKage  1.1b
Parallel software tools for high performance CMB DA analysis
 All Data Structures Files Functions Variables Typedefs Groups Pages
ring.c
Go to the documentation of this file.
1 
8 #ifdef W_MPI
9 
10 #include <mpi.h>
11 #include <stdlib.h>
12 #include <string.h>
13 
32 int ring_init(int *indices, int count, int **R, int *nR, int **S, int *nS, int steps, MPI_Comm comm){
33  int err, p, tag;
34  int size, rank, sp, rp;
35  int *buf, nbuf;
36  MPI_Request s_request, r_request;
37 
38  MPI_Comm_size(comm, &size);
39  MPI_Comm_rank(comm, &rank);
40  MPI_Allreduce( &count, &nbuf, 1, MPI_INT, MPI_MAX, comm); //compute the buffer size : max(count)_{comm}
41  buf = (int* ) malloc(nbuf*sizeof(int)); //allocate buffer
42  tag=0;
43  for (p=1; p < steps; p++){ //communication phase to get nb shared indices between peer of preocesses
44  sp=(rank+p)%size;
45  rp=(rank+size-p)%size;
46  MPI_Isend( &count, 1, MPI_INT, sp , 0, comm, &s_request); //send my number of indices
47  MPI_Irecv( &nbuf, 1, MPI_INT, rp, 0, comm, &r_request); //receive a number of indices
48  tag++;
49  MPI_Wait(&r_request, MPI_STATUS_IGNORE);
50  MPI_Irecv( buf, nbuf, MPI_INT, rp, tag, comm, &r_request); //receive indices tab
51  MPI_Wait(&s_request, MPI_STATUS_IGNORE);
52  MPI_Isend( indices, count, MPI_INT, sp, tag, comm, &s_request); //send indices tab
53  tag++;
54 
55  MPI_Wait(&r_request, MPI_STATUS_IGNORE);
56  nR[p] = card_and(indices, count, buf, nbuf); //compute number of shared indices
57  nS[steps-p]=nR[p];
58  R[p] = (int* ) malloc(nR[p]*sizeof(int)); //allocate receiving tab
59  S[steps-p] = (int* ) malloc(nS[steps-p]*sizeof(int)); //allocate sanding tab
60  map_and(indices, count, buf, nbuf, R[p]); //fill receiving tab
61  S[steps-p]=R[p]; //
62  }
63  free(buf);
64  nS[0]=0; //
65  nR[0]=0; //
66  return 0;
67 }
68 
82 int ring_reduce(int **R, int *nR, int nRmax, int **S, int *nS, int nSmax, double *val, double *res_val, int steps, MPI_Comm comm){
83  int tag, rank, size, p;
84  MPI_Request s_request, r_request;
85  int sp, rp;
86  double *sbuf, *rbuf;
87 
88  MPI_Comm_size(comm, &size);
89  MPI_Comm_rank(comm, &rank);
90  tag=0;
91 
92  rbuf = (double *) malloc(nRmax * sizeof(double));
93  sbuf = (double *) malloc(nSmax * sizeof(double));
94 
95  for (p=1; p < steps; p++){
96  rp=(rank+size-p)%size;
97  MPI_Irecv(rbuf, nR[p], MPI_DOUBLE, rp, tag, comm, &r_request);
98  sp=(rank+p)%size;
99  m2s(val, sbuf, S[p], nS[p]); //fill the sending buffer
100  MPI_Isend(sbuf, nS[p], MPI_DOUBLE, sp, tag, comm, &s_request);
101 
102  tag++;
103 
104  MPI_Wait(&r_request, MPI_STATUS_IGNORE);
105  s2m_sum(res_val, rbuf, R[p], nR[p]); //sum receive buffer into values
106 
107  MPI_Wait(&s_request, MPI_STATUS_IGNORE);
108  }
109  free(sbuf);
110  free(rbuf);
111  return 0;
112 }
113 
114 
126 int ring_nonblocking_reduce(int **R, int *nR, int **S, int *nS, double *val, double *res_val, int steps, MPI_Comm comm){
127  int tag, rank, size, p;
128  MPI_Request *s_request, *r_request;
129  int sp, rp;
130  double **sbuf, **rbuf;
131 
132  MPI_Comm_size(comm, &size);
133  MPI_Comm_rank(comm, &rank);
134  //printf("\n non_blocking rank %d", rank);
135 
136  s_request = (MPI_Request *) malloc((steps-1) * sizeof(MPI_Request));
137  r_request = (MPI_Request *) malloc((steps-1) * sizeof(MPI_Request));
138 
139  rbuf = (double **) malloc((steps-1) * sizeof(double *));
140  sbuf = (double **) malloc((steps-1) * sizeof(double *));
141 
142  for (p=1; p < steps; p++){
143  //printf("\n buf alloc %d", p);
144  rbuf[p-1] = (double *) malloc(nR[p] * sizeof(double));
145  sbuf[p-1] = (double *) malloc(nS[p] * sizeof(double));
146  m2s(val, sbuf[p-1], S[p], nS[p]); //fill the sending buffer
147  }
148 
149  tag=0;
150  for (p=1; p < steps; p++){
151  //printf("\n isend %d", p);
152  sp=(rank+p)%size;
153  rp=(rank+size-p)%size;
154 
155  MPI_Irecv(rbuf[p-1], nR[p], MPI_DOUBLE, rp, tag, comm, &r_request[p-1]);
156  MPI_Isend(sbuf[p-1], nS[p], MPI_DOUBLE, sp, tag, comm, &s_request[p-1]);
157  tag++;
158  }
159  MPI_Waitall(size-1, r_request, MPI_STATUSES_IGNORE);
160 
161  for (p=1; p < steps; p++){
162  s2m_sum(res_val, rbuf[p-1], R[p], nR[p]); //sum receive buffer into values
163  }
164  MPI_Waitall(size-1, s_request, MPI_STATUSES_IGNORE);
165  free(r_request);
166  free(s_request);
167  free(sbuf);
168  free(rbuf);
169  return 0;
170 }
171 
185 int ring_noempty_reduce(int **R, int *nR, int nneR, int **S, int *nS, int nneS, double *val, double *res_val, int steps, MPI_Comm comm){
186  int tag, rank, size, p;
187  MPI_Request *s_request, *r_request;
188  int sp, rp, nesi, neri;
189  double **sbuf, **rbuf;
190 
191  MPI_Comm_size(comm, &size);
192  MPI_Comm_rank(comm, &rank);
193  //printf("\n non_blocking rank %d", rank);
194 
195  s_request = (MPI_Request *) malloc(nneS * sizeof(MPI_Request));
196  r_request = (MPI_Request *) malloc(nneR * sizeof(MPI_Request));
197 
198  rbuf = (double **) malloc(nneR * sizeof(double *));
199  sbuf = (double **) malloc(nneS * sizeof(double *));
200 
201  nesi=0;
202  for (p=1; p < steps; p++){
203  if(nS[p] != 0){
204  sbuf[nesi] = (double *) malloc(nS[p] * sizeof(double));
205  m2s(val, sbuf[nesi], S[p], nS[p]); //fill the sending buffer
206  nesi++;
207  }
208  }
209 
210  tag=0;
211  nesi=0;
212  neri=0;
213  for (p=1; p < steps; p++){
214  sp=(rank+p)%size;
215  rp=(rank+size-p)%size;
216  if(nR[p] != 0){
217  rbuf[neri] = (double *) malloc(nR[p] * sizeof(double));
218  MPI_Irecv(rbuf[neri], nR[p], MPI_DOUBLE, rp, tag, comm, &r_request[neri]);
219  neri++;
220  }
221  if(nS[p] != 0){
222  MPI_Isend(sbuf[nesi], nS[p], MPI_DOUBLE, sp, tag, comm, &s_request[nesi]);
223  nesi++;
224  }
225  tag++;
226  }
227  MPI_Waitall(nneR, r_request, MPI_STATUSES_IGNORE);
228 
229  neri=0;
230  for (p=1; p < steps; p++){
231  if(nR[p] != 0){
232  s2m_sum(res_val, rbuf[neri], R[p], nR[p]); //sum receive buffer into values
233  neri++;
234  }
235  }
236  MPI_Waitall(nneS, s_request, MPI_STATUSES_IGNORE);
237  free(r_request);
238  free(s_request);
239  free(sbuf);
240  free(rbuf);
241  return 0;
242 }
243 #endif