00001 #ifndef BR_TREE_ALGORITHMS_HPP_
00002 #define BR_TREE_ALGORITHMS_HPP_
00003
00004 #include "libmpib_coll.h"
00005 #include "mpib_coll.h"
00006 #include "comm_tree.hpp"
00007 #include <boost/graph/graphviz.hpp>
00008
00013 template <typename Builder>
00014 int MPIB_Bcast_tree_algorithm(Builder builder, MPIB_child_traverse_order order,
00015 void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
00016 {
00017 using namespace MPIB::comm_tree;
00018
00019 int size;
00020 MPI_Comm_size(comm, &size);
00021 int rank;
00022 MPI_Comm_rank(comm, &rank);
00023 MPI_Aint extent;
00024 MPI_Type_extent(datatype, &extent);
00025 Graph graph;
00026 Vertex r, u, v;
00027 builder.build(size, root, rank, count * extent, graph, r, u, v);
00028 if ((rank == root) && mpib_coll_verbose)
00029 write_graphviz(cout, graph);
00030
00031 if (rank != root)
00032 MPI_Recv(buffer, count, datatype, get(vertex_index, graph, u), 0, comm, MPI_STATUS_IGNORE);
00033
00034 MPI_Request* reqs = (MPI_Request*)malloc(sizeof(MPI_Request) * size);
00035 int child_counter;
00036 Adjacency_iterator ai, ai_end;
00037 tie(ai, ai_end) = adjacent_vertices(v, graph);
00038 for(child_counter = 0; ai != ai_end; child_counter++)
00039 MPI_Isend(buffer, count, datatype, get(vertex_index, graph, (order == R2L) ? *(--ai_end) : *(ai++)), 0, comm, &reqs[child_counter]);
00040 MPI_Waitall(child_counter, reqs, MPI_STATUSES_IGNORE);
00041 free(reqs);
00042 return MPI_SUCCESS;
00043 }
00044
00065 template <typename Builder>
00066 int MPIB_Reduce_tree_algorithm(Builder builder, MPIB_child_traverse_order order,
00067 void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype,
00068 MPI_Op op, int root, MPI_Comm comm)
00069 {
00070 using namespace MPIB::comm_tree;
00071
00072 int size;
00073 MPI_Comm_size(comm, &size);
00074 int rank;
00075 MPI_Comm_rank(comm, &rank);
00076 MPI_Aint extent;
00077 MPI_Type_extent(datatype, &extent);
00078 Graph graph;
00079 Vertex r, u, v;
00080 builder.build(size, root, rank, count * extent, graph, r, u, v);
00081 if ((rank == root) && mpib_coll_verbose)
00082 write_graphviz(cout, graph);
00083
00084 void* buffer = rank == root ? recvbuf : malloc(count * extent);
00085
00086 MPI_Request* reqs = (MPI_Request*)malloc(sizeof(MPI_Request) * size);
00087 int child_counter;
00088 Adjacency_iterator ai, ai_end;
00089 tie(ai, ai_end) = adjacent_vertices(v, graph);
00090 for(child_counter = 0; ai != ai_end; child_counter++) {
00091 MPI_Irecv(buffer, count, datatype, get(vertex_index, graph, (order == R2L) ? *(--ai_end) : *(ai++)), 0, comm, &reqs[child_counter]);
00092
00093 }
00094 MPI_Waitall(child_counter, reqs, MPI_STATUSES_IGNORE);
00095 free(reqs);
00096
00097 if (rank != root) {
00098 MPI_Send(buffer, count, datatype, get(vertex_index, graph, u), 0, comm);
00099 free(buffer);
00100 }
00101 return MPI_SUCCESS;
00102 }
00103
00104 #endif