1 /* 2 Copyright (c) 2009, David Cheng, Viral B. Shah. 3 4 Permission is hereby granted, free of charge, to any person obtaining a copy 5 of this software and associated documentation files (the "Software"), to deal 6 in the Software without restriction, including without limitation the rights 7 to use, copy, modify, merge, publish, distribute, sublicense, and/or sell 8 copies of the Software, and to permit persons to whom the Software is 9 furnished to do so, subject to the following conditions: 10 11 The above copyright notice and this permission notice shall be included in 12 all copies or substantial portions of the Software. 13 14 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR 15 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, 16 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE 17 AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER 18 LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, 19 OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN 20 THE SOFTWARE. 21 */ 22 23 #ifndef PSORT_H 24 #define PSORT_H 25 26 #include "psort_util.h" 27 28 namespace vpsort { 29 using namespace std; 30 31 /* 32 SeqSort can be STLSort, STLStableSort 33 Split can be MedianSplit, SampleSplit 34 Merge can be FlatMerge, TreeMerge, OOPTreeMerge, FunnelMerge2, FunnelMerge4 35 */ 36 template<typename _RandomAccessIter, typename _Compare, 37 typename _SeqSortType, typename _SplitType, typename _MergeType> parallel_sort(_RandomAccessIter first,_RandomAccessIter last,_Compare comp,long * dist_in,SeqSort<_SeqSortType> & mysort,Split<_SplitType> & mysplit,Merge<_MergeType> & mymerge,MPI_Comm comm)38 void parallel_sort (_RandomAccessIter first, _RandomAccessIter last, 39 _Compare comp, long *dist_in, 40 SeqSort<_SeqSortType> &mysort, 41 Split<_SplitType> &mysplit, 42 Merge<_MergeType> &mymerge, 43 MPI_Comm comm) { 44 45 typedef typename iterator_traits<_RandomAccessIter>::value_type _ValueType; 46 typedef typename iterator_traits<_RandomAccessIter>::difference_type _Distance; 47 48 int nproc, rank; 49 MPI_Comm_size (comm, &nproc); 50 MPI_Comm_rank (comm, &rank); 51 52 MPI_Datatype MPI_valueType, MPI_distanceType; 53 MPI_Type_contiguous (sizeof(_ValueType), MPI_CHAR, &MPI_valueType); 54 MPI_Type_commit (&MPI_valueType); // ABAB: Any type committed needs to be freed to claim storage 55 MPI_Type_contiguous (sizeof(_Distance), MPI_CHAR, &MPI_distanceType); 56 MPI_Type_commit (&MPI_distanceType); 57 58 _Distance *dist = new _Distance[nproc]; 59 for (int i=0; i<nproc; ++i) dist[i] = (_Distance) dist_in[i]; 60 61 // Sort the data locally 62 // ABAB: Progress calls use MPI::COMM_WORLD, instead of the passed communicator 63 // Since they are just debug prints, avoid them 64 progress (rank, 0, mysort.description(), comm); 65 mysort.seqsort (first, last, comp); 66 67 if (nproc == 1) 68 { 69 MPI_Type_free(&MPI_valueType); 70 MPI_Type_free(&MPI_distanceType); 71 return; 72 } 73 74 // Find splitters 75 progress (rank, 1, mysplit.description(), comm); 76 // explicit vector ( size_type n, const T& value= T(), const Allocator& = Allocator() ); 77 // repetitive sequence constructor: Initializes the vector with its content set to a repetition, n times, of copies of value. 78 vector< vector<_Distance> > right_ends(nproc + 1, vector<_Distance>(nproc, 0)); 79 // ABAB: The following hangs if some dist entries are zeros, i.e. first = last for some processors 80 mysplit.split (first, last, dist, comp, right_ends, 81 MPI_valueType, MPI_distanceType, comm); 82 83 // Communicate to destination 84 progress (rank, 2, const_cast<char *>(string("alltoall").c_str()), comm); 85 _Distance n_loc = last - first; 86 _ValueType *trans_data = new _ValueType[n_loc]; 87 _Distance *boundaries = new _Distance[nproc+1]; 88 alltoall (right_ends, first, last, 89 trans_data, boundaries, 90 MPI_valueType, MPI_distanceType, comm); 91 92 // Merge streams from all processors 93 // progress (rank, 3, mymerge.description()); 94 mymerge.merge (trans_data, first, boundaries, nproc, comp); 95 96 delete [] boundaries; 97 delete [] dist; 98 delete [] trans_data; 99 MPI_Type_free (&MPI_valueType); 100 MPI_Type_free (&MPI_distanceType); 101 102 // Finish 103 progress (rank, 4, const_cast<char *>(string("finish").c_str()), comm); 104 105 return; 106 } 107 108 template<typename _RandomAccessIter, typename _Compare> parallel_sort(_RandomAccessIter first,_RandomAccessIter last,_Compare comp,long * dist,MPI_Comm comm)109 void parallel_sort (_RandomAccessIter first, 110 _RandomAccessIter last, 111 _Compare comp, long *dist, MPI_Comm comm) { 112 113 STLSort mysort; 114 MedianSplit mysplit; 115 OOPTreeMerge mymerge; 116 117 parallel_sort (first, last, comp, dist, mysort, mysplit, mymerge, comm); 118 } 119 120 template<typename _RandomAccessIter> parallel_sort(_RandomAccessIter first,_RandomAccessIter last,long * dist,MPI_Comm comm)121 void parallel_sort (_RandomAccessIter first, 122 _RandomAccessIter last, 123 long *dist, MPI_Comm comm) { 124 125 typedef typename iterator_traits<_RandomAccessIter>::value_type _ValueType; 126 127 STLSort mysort; 128 MedianSplit mysplit; 129 OOPTreeMerge mymerge; 130 131 parallel_sort (first, last, less<_ValueType>(), 132 dist, mysort, mysplit, mymerge, comm); 133 } 134 135 } 136 137 138 #endif /* PSORT_H */ 139