COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
BitMapCarousel.h
Go to the documentation of this file.
1/****************************************************************/
2/* Parallel Combinatorial BLAS Library (for Graph Computations) */
3/* version 1.4 -------------------------------------------------*/
4/* date: 1/17/2014 ---------------------------------------------*/
5/* authors: Aydin Buluc (abuluc@lbl.gov), Adam Lugowski --------*/
6/* this file contributed by Scott Beamer of UC Berkeley --------*/
7/****************************************************************/
8/*
9 Copyright (c) 2010-2014, The Regents of the University of California
10
11 Permission is hereby granted, free of charge, to any person obtaining a copy
12 of this software and associated documentation files (the "Software"), to deal
13 in the Software without restriction, including without limitation the rights
14 to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
15 copies of the Software, and to permit persons to whom the Software is
16 furnished to do so, subject to the following conditions:
17
18 The above copyright notice and this permission notice shall be included in
19 all copies or substantial portions of the Software.
20
21 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
22 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
23 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
24 AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
25 LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
26 OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
27 THE SOFTWARE.
28 */
29
30
31#ifndef BITMAPCAROUSEL_H
32#define BITMAPCAROUSEL_H
33
34// #include <algorithm>
35#include "BitMap.h"
36#include "BitMapFringe.h"
37#include "CommGrid.h"
38
39namespace combblas {
40
41template <class IT, class NT>
42class BitMapCarousel {
43 public:
44 BitMapCarousel(std::shared_ptr<CommGrid> grid, IT glen, int local_subword_disp) {
45 commGrid.reset(new CommGrid(*grid));
46 rotation_index = 0;
47 global_size = glen;
48 my_procrow = commGrid->GetRankInProcCol();
49 my_proccol = commGrid->GetRankInProcRow();
50 procrows = commGrid->GetGridRows();
51 proccols = commGrid->GetGridCols();
52 local_size = SizeOfChunk();
53 rotlenuntil = RotLengthUntil();
54 IT biggest_size = global_size - RotLengthUntil(procrows-1, proccols-1) + 63;
55 bm = new BitMap(biggest_size);
56 recv_buff = new BitMap(biggest_size);
57 old_bm = new BitMap(biggest_size);
58 sub_disps = new int[proccols];
59 sub_disps[my_proccol] = local_subword_disp;
60 MPI_Allgather(MPI_IN_PLACE, 1, MPI_INT, sub_disps, 1, MPI_INT, commGrid->GetRowWorld());
61 curr_subword_disp = local_subword_disp;
62 }
63
65 delete bm;
66 delete recv_buff;
67 delete old_bm;
68 delete[] sub_disps;
69 }
70
71 // Return the global index of the start of the local chunk (inclusive)
73 return rotlenuntil;
74 }
75
77 return RotLengthUntil(my_procrow,col);
78 }
79
80 // Return the global index of the end of the local chunk (exclusive)
82 return rotlenuntil + local_size;
83 }
84
86 return GetGlobalStartOfLocal(col) + local_size;
87 }
88
89 bool GetBit(IT index) const {
90 return bm->get_bit(index - rotlenuntil + curr_subword_disp);
91 }
92
93 void SetBit(IT index) {
94 bm->set_bit(index - rotlenuntil + curr_subword_disp);
95 }
96
97 template <class NT1>
99 bm->reset();
100 local_size = x.LocArrSize();
101 for (DenseVectorLocalIterator<IT,NT> it(x); it.HasNext(); it.Next())
102 if (it.GetValue() != -1)
103 bm->set_bit(it.GetLocIndex() + curr_subword_disp);
104 }
105
107 int curr_col = (my_proccol + rotation_index) % proccols;
108 return RotLengthUntil(my_procrow, curr_col);
109 }
110
111 IT RotLengthUntil(int my_procrow, int my_proccol) const {
112 IT n_perprocrow = global_size / procrows;
114 if(my_procrow == procrows-1)
115 n_thisrow = global_size - (n_perprocrow*(procrows-1));
116 else
118 IT n_perproc = n_thisrow / proccols;
119 return ((n_perprocrow * my_procrow)+(n_perproc*my_proccol));
120 }
121
122 IT SizeOfChunk() const {
123 int curr_col = (my_proccol + rotation_index) % proccols;
124 return SizeOfChunk(my_procrow, curr_col);
125 }
126
127 IT SizeOfChunk(int my_procrow, int my_proccol) const {
129 if (my_proccol == (proccols-1)) {
130 if (my_procrow == (procrows-1)) {
131 my_upper_limit = global_size;
132 }else {
133 my_upper_limit = RotLengthUntil(my_procrow+1, 0);
134 }
135 } else {
136 my_upper_limit = RotLengthUntil(my_procrow, (my_proccol+1) % proccols);
137 }
138 return my_upper_limit - RotLengthUntil(my_procrow,my_proccol);
139 }
140
142 int source = (my_proccol+1) % proccols;
143 int dest = (my_proccol+(proccols-1)) % proccols;
144 rotation_index = (rotation_index + 1) % proccols;
145 long send_words = (local_size + 63 + curr_subword_disp)>>6;
146 long recv_words = (SizeOfChunk() + 63 + sub_disps[(my_proccol+rotation_index)%proccols])>>6;
147
148 MPI_Comm RowWorld = commGrid->GetRowWorld();
152
153 local_size = SizeOfChunk();
154 rotlenuntil = RotLengthUntil();
155 std::swap(bm, recv_buff);
156 curr_subword_disp = sub_disps[(my_proccol+rotation_index)%proccols];
157 }
158
160 uint64_t* dest = bm_fringe.AccessBM()->data();
161 uint64_t* curr = bm->data();
162 uint64_t* old = old_bm->data();
163 IT num_words = (local_size + 63 + curr_subword_disp) / 64;
164 for (IT i=0; i<num_words; i++) {
165 dest[i] = curr[i] ^ old[i];
166 }
167 }
168
169 void SaveOld() {
170 old_bm->copy_from(bm);
171 }
172
173private:
174 std::shared_ptr<CommGrid> commGrid;
175 int rotation_index;
176 int my_procrow;
177 int my_proccol;
178 int procrows;
179 int proccols;
180 int curr_subword_disp;
181 IT rotlenuntil;
182 IT global_size;
183 IT local_size;
184 BitMap* bm;
185 BitMap* recv_buff;
186 BitMap* old_bm;
187 int* sub_disps;
188};
189
190}
191
192#endif // BITMAPCAROUSEL_H
int64_t IT
IT SizeOfChunk(int my_procrow, int my_proccol) const
IT GetGlobalEndOfLocal(int col) const
IT GetGlobalStartOfLocal(int col) const
BitMapCarousel(std::shared_ptr< CommGrid > grid, IT glen, int local_subword_disp)
bool GetBit(IT index) const
IT RotLengthUntil(int my_procrow, int my_proccol) const
void UpdateFringe(BitMapFringe< IT, NT > &bm_fringe)
void LoadVec(FullyDistVec< IT, NT1 > &x)
void reset()
Definition BitMap.h:79
bool get_bit(uint64_t pos)
Definition BitMap.h:106
void set_bit(uint64_t pos)
Definition BitMap.h:85
uint64_t * data()
Definition BitMap.h:145
void copy_from(const BitMap *other)
Definition BitMap.h:149
#define ROTATE
Definition SpDefs.h:109