Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright notices
[simgrid.git] / src / smpi / colls / allgather-loosely-lr.c
1 /* Copyright (c) 2013-2014. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "colls_private.h"
8
9 int smpi_coll_tuned_allgather_loosely_lr(void *sbuf, int scount,
10                                          MPI_Datatype stype, void *rbuf,
11                                          int rcount, MPI_Datatype rtype,
12                                          MPI_Comm comm)
13 {
14   int comm_size, rank;
15   int tag = COLL_TAG_ALLGATHER;
16   int i, j, send_offset, recv_offset;
17   int intra_rank, inter_rank, inter_comm_size, intra_comm_size;
18   int inter_dst, inter_src;
19
20   comm_size = smpi_comm_size(comm);
21
22 if(smpi_comm_get_leaders_comm(comm)==MPI_COMM_NULL){
23     smpi_comm_init_smp(comm);
24   }
25   int num_core=1;
26   if (smpi_comm_is_uniform(comm)){
27     num_core = smpi_comm_size(smpi_comm_get_intra_comm(comm));
28   }
29
30   if(comm_size%num_core)
31     THROWF(arg_error,0, "allgather loosely lr algorithm can't be used with non multiple of NUM_CORE=%d number of processes ! ",num_core);
32
33   rank = smpi_comm_rank(comm);
34   MPI_Aint rextent, sextent;
35   rextent = smpi_datatype_get_extent(rtype);
36   sextent = smpi_datatype_get_extent(stype);
37   MPI_Request inter_rrequest;
38   MPI_Request rrequest_array[128];
39   MPI_Request srequest_array[128];
40   MPI_Request inter_srequest_array[128];
41
42
43   int rrequest_count = 0;
44   int srequest_count = 0;
45   int inter_srequest_count = 0;
46
47   MPI_Status status;
48
49   intra_rank = rank % num_core;
50   inter_rank = rank / num_core;
51   inter_comm_size = (comm_size + num_core - 1) / num_core;
52   intra_comm_size = num_core;
53
54   int src_seg, dst_seg;
55
56   //copy corresponding message from sbuf to rbuf
57   recv_offset = rank * rextent * rcount;
58   smpi_mpi_sendrecv(sbuf, scount, stype, rank, tag,
59                (char *)rbuf + recv_offset, rcount, rtype, rank, tag, comm, &status);
60
61   int dst, src;
62   int inter_send_offset, inter_recv_offset;
63
64   rrequest_count = 0;
65   srequest_count = 0;
66   inter_srequest_count = 0;
67
68   for (i = 0; i < inter_comm_size; i++) {
69
70     // inter_communication
71
72     inter_dst = (rank + intra_comm_size) % comm_size;
73     inter_src = (rank - intra_comm_size + comm_size) % comm_size;
74
75     src_seg =
76         ((inter_rank - 1 - i +
77           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
78     dst_seg =
79         ((inter_rank - i +
80           inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
81
82     inter_send_offset = dst_seg * sextent * scount;
83     inter_recv_offset = src_seg * rextent * rcount;
84
85     for (j = 0; j < intra_comm_size; j++) {
86
87       // inter communication
88       if (intra_rank == j) {
89         if (i != inter_comm_size - 1) {
90
91           inter_rrequest = smpi_mpi_irecv((char *)rbuf + inter_recv_offset, rcount, rtype,
92                                           inter_src, tag, comm);
93           inter_srequest_array[inter_srequest_count++] = smpi_mpi_isend((char *)rbuf + inter_send_offset, scount, stype,
94                                                                         inter_dst, tag, comm);
95         }
96       }
97       //intra_communication
98       src = inter_rank * intra_comm_size + j;
99       dst = inter_rank * intra_comm_size + j;
100
101       src_seg =
102           ((inter_rank - i +
103             inter_comm_size) % inter_comm_size) * intra_comm_size + j;
104       dst_seg =
105           ((inter_rank - i +
106             inter_comm_size) % inter_comm_size) * intra_comm_size + intra_rank;
107
108       send_offset = dst_seg * sextent * scount;
109       recv_offset = src_seg * rextent * rcount;
110
111
112       if (j != intra_rank) {
113
114         rrequest_array[rrequest_count++] = smpi_mpi_irecv((char *)rbuf + recv_offset, rcount, rtype, src, tag, comm);
115         srequest_array[srequest_count++] = smpi_mpi_isend((char *)rbuf + send_offset, scount, stype, dst, tag, comm);
116
117       }
118     }                           // intra loop
119
120
121     // wait for inter communication to finish for these rounds (# of round equals num_core)
122     if (i != inter_comm_size - 1) {
123       smpi_mpi_wait(&inter_rrequest, &status);
124     }
125
126   }                             //inter loop
127
128   smpi_mpi_waitall(rrequest_count, rrequest_array, MPI_STATUSES_IGNORE);
129   smpi_mpi_waitall(srequest_count, srequest_array, MPI_STATUSES_IGNORE);
130   smpi_mpi_waitall(inter_srequest_count, inter_srequest_array, MPI_STATUSES_IGNORE);
131
132   return MPI_SUCCESS;
133 }