00001 #ifndef SG_TREE_ALGORITHMS_HPP_
00002 #define SG_TREE_ALGORITHMS_HPP_
00003
00004 #include "libmpib_coll.h"
00005 #include "comm_tree.hpp"
00006 #include <boost/graph/graphviz.hpp>
00007
00008 namespace MPIB {
00010 namespace SG {
00011 using namespace comm_tree;
00012
00017 class Assembler {
00018 private:
00020 const map<int, int>& rank2index;
00022 const int rscount;
00024 int& count;
00026 int* lengths;
00028 int* indices;
00029 public:
00030 Assembler(const map<int, int>& _rank2index, const int _rscount,
00031 int& _count, int* _lengths, int* _indices):
00032 rank2index(_rank2index), rscount(_rscount),
00033 count(_count), lengths(_lengths), indices(_indices) {
00034 count = 0;
00035 }
00036 void preorder(Vertex vertex, Tree& tree) {
00037 lengths[count] = rscount;
00038 indices[count] = rscount * rank2index.find(get(vertex_index, tree, vertex))->second;
00039 count++;
00040 }
00041 void inorder(Vertex vertex, Tree& tree) {}
00042 void postorder(Vertex vertex, Tree& tree) {}
00043 };
00044
00048 class Indexer {
00049 private:
00051 int& index;
00053 map<int, int>& rank2index;
00054 public:
00055 Indexer(int& _index, map<int, int>& _rank2index):
00056 index(_index), rank2index(_rank2index) {
00057 index = 0;
00058 }
00059 void preorder(Vertex vertex, Tree& tree) {
00060 rank2index[get(vertex_index, tree, vertex)] = index++;
00061 }
00062 void inorder(Vertex vertex, Tree& tree) {}
00063 void postorder(Vertex vertex, Tree& tree) {}
00064 };
00065
00067 class Edge_writer {
00068 private:
00069 Graph& graph;
00070 const int rscount;
00071
00073 class Indexer {
00074 private:
00075 const int rscount;
00076 int& count;
00077 public:
00078 Indexer(const int rscount, int& _count): rscount(rscount), count(_count) {
00079 count = 0;
00080 }
00081 void preorder(Vertex vertex, Tree& tree) {
00082 count += rscount;
00083 }
00084 void inorder(Vertex vertex, Tree& tree) {}
00085 void postorder(Vertex vertex, Tree& tree) {}
00086 };
00087
00088 public:
00089 Edge_writer(Graph& _graph, const int _rscount):
00090 graph(_graph), rscount(_rscount) {}
00091 void operator()(std::ostream& out, const Edge& e) const {
00092 Vertex v = target(e, graph);
00093 Tree tree = Tree(graph, v);
00094 int count;
00095 traverse_tree(v, tree, Indexer(rscount, count));
00096 out << "[label = \"" << count << "\"]";
00097 }
00098 };
00099 }
00100 }
00101
00106 template <typename Builder>
00107 int MPIB_Scatter_tree_algorithm(Builder builder, MPIB_child_traverse_order order,
00108 void* sendbuf, int sendcount, MPI_Datatype sendtype,
00109 void* recvbuf, int recvcount, MPI_Datatype recvtype,
00110 int root, MPI_Comm comm)
00111 {
00112 using namespace MPIB::SG;
00113
00114 int size;
00115 MPI_Comm_size(comm, &size);
00116 int rank;
00117 MPI_Comm_rank(comm, &rank);
00118 MPI_Aint recvext;
00119 MPI_Type_extent(recvtype, &recvext);
00120 Graph graph;
00121 Vertex r, u, v;
00122 builder.build(size, root, rank, recvcount * recvext, graph, r, u, v);
00123 if ((rank == root) && mpib_coll_verbose)
00124 write_graphviz(cout, graph, default_writer(), Edge_writer(graph, recvcount * recvext));
00125
00126 Tree tree = Tree(graph, r);
00127
00128 map<int, int> rank2index;
00129 char* buffer = NULL;
00130 if (rank == root) {
00131
00132 for (int i = 0; i < size; i++)
00133 rank2index[i] = i;
00134 buffer = (char*)sendbuf;
00135 } else {
00136
00137 int index;
00138 traverse_tree(v, tree, Indexer(index, rank2index));
00139 int count = rank2index.size() * recvcount;
00140 buffer = (char*)malloc(sizeof(char) * count * recvext);
00141
00142 MPI_Recv(buffer, count, recvtype, get(vertex_index, graph, u), 0, comm, MPI_STATUS_IGNORE);
00143 }
00144
00145 memcpy(recvbuf, buffer + rank2index[rank] * recvcount * recvext, recvcount * recvext);
00146
00147 int subtree_size = rank2index.size() - 1;
00148 int count;
00149 int* lengths = (int*)malloc(sizeof(int) * subtree_size);
00150 int* indices = (int*)malloc(sizeof(int) * subtree_size);
00151 MPI_Datatype* dts = (MPI_Datatype*)malloc(sizeof(MPI_Datatype) * subtree_size);
00152 MPI_Request* reqs = (MPI_Request*)malloc(sizeof(MPI_Request) * subtree_size);
00153 int child_counter;
00154 Adjacency_iterator ai, ai_end;
00155 tie(ai, ai_end) = adjacent_vertices(v, graph);
00156 for(child_counter = 0; ai != ai_end; child_counter++) {
00157 Vertex t = (order == R2L) ? *(--ai_end) : *(ai++);
00158
00159 traverse_tree(t, tree,
00160 Assembler(rank2index, recvcount, count, lengths, indices));
00161
00162 MPI_Type_indexed(count, lengths, indices, recvtype, &dts[child_counter]);
00163 MPI_Type_commit(&dts[child_counter]);
00164 MPI_Isend(buffer, 1, dts[child_counter], get(vertex_index, graph, t), 0, comm, &reqs[child_counter]);
00165 }
00166 free(lengths);
00167 free(indices);
00168 MPI_Waitall(child_counter, reqs, MPI_STATUSES_IGNORE);
00169 free(reqs);
00170 for (int i = 0; i < child_counter; i++)
00171 MPI_Type_free(&dts[i]);
00172 free(dts);
00173 if (rank != root)
00174 free(buffer);
00175 return MPI_SUCCESS;
00176 }
00177
00182 template <typename Builder>
00183 int MPIB_Gather_tree_algorithm(Builder builder, MPIB_child_traverse_order order,
00184 void* sendbuf, int sendcount, MPI_Datatype sendtype,
00185 void* recvbuf, int recvcount, MPI_Datatype recvtype,
00186 int root, MPI_Comm comm)
00187 {
00188 using namespace MPIB::SG;
00189
00190 int size;
00191 MPI_Comm_size(comm, &size);
00192 int rank;
00193 MPI_Comm_rank(comm, &rank);
00194 MPI_Aint sendext;
00195 MPI_Type_extent(sendtype, &sendext);
00196 Graph graph;
00197 Vertex r, u, v;
00198 builder.build(size, root, rank, sendcount * sendext, graph, r, u, v);
00199 if ((rank == root) && mpib_coll_verbose)
00200 write_graphviz(cout, graph, default_writer(), Edge_writer(graph, sendcount * sendext));
00201
00202 Tree tree = Tree(graph, r);
00203
00204 map<int, int> rank2index;
00205 char* buffer = NULL;
00206 if (rank == root) {
00207
00208 for (int i = 0; i < size; i++)
00209 rank2index[i] = i;
00210 buffer = (char*)recvbuf;
00211 } else {
00212
00213 int index;
00214 traverse_tree(v, tree, Indexer(index, rank2index));
00215 buffer = (char*)malloc(sizeof(char) * rank2index.size() * sendcount * sendext);
00216 }
00217
00218 memcpy(buffer + rank2index[rank] * sendcount * sendext, sendbuf, sendcount * sendext);
00219
00220 int subtree_size = rank2index.size() - 1;
00221 int count;
00222 int* lengths = (int*)malloc(sizeof(int) * subtree_size);
00223 int* indices = (int*)malloc(sizeof(int) * subtree_size);
00224 MPI_Datatype* dts = (MPI_Datatype*)malloc(sizeof(MPI_Datatype) * subtree_size);
00225 MPI_Request* reqs = (MPI_Request*)malloc(sizeof(MPI_Request) * subtree_size);
00226 int child_counter;
00227 Adjacency_iterator ai, ai_end;
00228 tie(ai, ai_end) = adjacent_vertices(v, graph);
00229 for(child_counter = 0; ai != ai_end; child_counter++) {
00230 Vertex t = (order == R2L) ? *(--ai_end) : *(ai++);
00231
00232 traverse_tree(t, tree,
00233 Assembler(rank2index, sendcount, count, lengths, indices));
00234
00235 MPI_Type_indexed(count, lengths, indices, sendtype, &dts[child_counter]);
00236 MPI_Type_commit(&dts[child_counter]);
00237 MPI_Irecv(buffer, 1, dts[child_counter], get(vertex_index, graph, t), 0, comm, &reqs[child_counter]);
00238 }
00239 free(lengths);
00240 free(indices);
00241 MPI_Waitall(child_counter, reqs, MPI_STATUSES_IGNORE);
00242 free(reqs);
00243 for (int i = 0; i < child_counter; i++)
00244 MPI_Type_free(&dts[i]);
00245 free(dts);
00246
00247 if (rank != root)
00248 MPI_Send(buffer, rank2index.size() * sendcount, sendtype, get(vertex_index, graph, u), 0, comm);
00249 if (rank != root)
00250 free(buffer);
00251 return MPI_SUCCESS;
00252 }
00253
00254 #endif