root / src / pgesv / HPL_spreadN.c @ 9
Historique | Voir | Annoter | Télécharger (12,54 ko)
1 | 1 | equemene | /*
|
---|---|---|---|
2 | 1 | equemene | * -- High Performance Computing Linpack Benchmark (HPL)
|
3 | 1 | equemene | * HPL - 2.0 - September 10, 2008
|
4 | 1 | equemene | * Antoine P. Petitet
|
5 | 1 | equemene | * University of Tennessee, Knoxville
|
6 | 1 | equemene | * Innovative Computing Laboratory
|
7 | 1 | equemene | * (C) Copyright 2000-2008 All Rights Reserved
|
8 | 1 | equemene | *
|
9 | 1 | equemene | * -- Copyright notice and Licensing terms:
|
10 | 1 | equemene | *
|
11 | 1 | equemene | * Redistribution and use in source and binary forms, with or without
|
12 | 1 | equemene | * modification, are permitted provided that the following conditions
|
13 | 1 | equemene | * are met:
|
14 | 1 | equemene | *
|
15 | 1 | equemene | * 1. Redistributions of source code must retain the above copyright
|
16 | 1 | equemene | * notice, this list of conditions and the following disclaimer.
|
17 | 1 | equemene | *
|
18 | 1 | equemene | * 2. Redistributions in binary form must reproduce the above copyright
|
19 | 1 | equemene | * notice, this list of conditions, and the following disclaimer in the
|
20 | 1 | equemene | * documentation and/or other materials provided with the distribution.
|
21 | 1 | equemene | *
|
22 | 1 | equemene | * 3. All advertising materials mentioning features or use of this
|
23 | 1 | equemene | * software must display the following acknowledgement:
|
24 | 1 | equemene | * This product includes software developed at the University of
|
25 | 1 | equemene | * Tennessee, Knoxville, Innovative Computing Laboratory.
|
26 | 1 | equemene | *
|
27 | 1 | equemene | * 4. The name of the University, the name of the Laboratory, or the
|
28 | 1 | equemene | * names of its contributors may not be used to endorse or promote
|
29 | 1 | equemene | * products derived from this software without specific written
|
30 | 1 | equemene | * permission.
|
31 | 1 | equemene | *
|
32 | 1 | equemene | * -- Disclaimer:
|
33 | 1 | equemene | *
|
34 | 1 | equemene | * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
|
35 | 1 | equemene | * ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
|
36 | 1 | equemene | * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
|
37 | 1 | equemene | * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE UNIVERSITY
|
38 | 1 | equemene | * OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
|
39 | 1 | equemene | * SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
|
40 | 1 | equemene | * LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
|
41 | 1 | equemene | * DATA OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
|
42 | 1 | equemene | * THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
|
43 | 1 | equemene | * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
|
44 | 1 | equemene | * OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
|
45 | 1 | equemene | * ---------------------------------------------------------------------
|
46 | 1 | equemene | */
|
47 | 1 | equemene | /*
|
48 | 1 | equemene | * Include files
|
49 | 1 | equemene | */
|
50 | 1 | equemene | #include "hpl.h" |
51 | 1 | equemene | |
52 | 1 | equemene | #ifdef STDC_HEADERS
|
53 | 1 | equemene | void HPL_spreadN
|
54 | 1 | equemene | ( |
55 | 1 | equemene | HPL_T_panel * PBCST, |
56 | 1 | equemene | int * IFLAG,
|
57 | 1 | equemene | HPL_T_panel * PANEL, |
58 | 1 | equemene | const enum HPL_SIDE SIDE, |
59 | 1 | equemene | const int N, |
60 | 1 | equemene | double * U,
|
61 | 1 | equemene | const int LDU, |
62 | 1 | equemene | const int SRCDIST, |
63 | 1 | equemene | const int * IPLEN, |
64 | 1 | equemene | const int * IPMAP, |
65 | 1 | equemene | const int * IPMAPM1 |
66 | 1 | equemene | ) |
67 | 1 | equemene | #else
|
68 | 1 | equemene | void HPL_spreadN
|
69 | 1 | equemene | ( PBCST, IFLAG, PANEL, SIDE, N, U, LDU, SRCDIST, IPLEN, IPMAP, IPMAPM1 ) |
70 | 1 | equemene | HPL_T_panel * PBCST; |
71 | 1 | equemene | int * IFLAG;
|
72 | 1 | equemene | HPL_T_panel * PANEL; |
73 | 1 | equemene | const enum HPL_SIDE SIDE; |
74 | 1 | equemene | const int N; |
75 | 1 | equemene | double * U;
|
76 | 1 | equemene | const int LDU; |
77 | 1 | equemene | const int SRCDIST; |
78 | 1 | equemene | const int * IPLEN; |
79 | 1 | equemene | const int * IPMAP; |
80 | 1 | equemene | const int * IPMAPM1; |
81 | 1 | equemene | #endif
|
82 | 1 | equemene | { |
83 | 1 | equemene | /*
|
84 | 1 | equemene | * Purpose
|
85 | 1 | equemene | * =======
|
86 | 1 | equemene | *
|
87 | 1 | equemene | * HPL_spreadN spreads the local array containing local pieces of U, so
|
88 | 1 | equemene | * that on exit to this function, a piece of U is contained in every
|
89 | 1 | equemene | * process row. The array IPLEN contains the number of rows of U, that
|
90 | 1 | equemene | * should be spread on any given process row. This function also probes
|
91 | 1 | equemene | * for the presence of the column panel PBCST. In case of success, this
|
92 | 1 | equemene | * panel will be forwarded. If PBCST is NULL on input, this probing
|
93 | 1 | equemene | * mechanism will be disabled.
|
94 | 1 | equemene | *
|
95 | 1 | equemene | * Arguments
|
96 | 1 | equemene | * =========
|
97 | 1 | equemene | *
|
98 | 1 | equemene | * PBCST (local input/output) HPL_T_panel *
|
99 | 1 | equemene | * On entry, PBCST points to the data structure containing the
|
100 | 1 | equemene | * panel (to be broadcast) information.
|
101 | 1 | equemene | *
|
102 | 1 | equemene | * IFLAG (local input/output) int *
|
103 | 1 | equemene | * On entry, IFLAG indicates whether or not the broadcast has
|
104 | 1 | equemene | * already been completed. If not, probing will occur, and the
|
105 | 1 | equemene | * outcome will be contained in IFLAG on exit.
|
106 | 1 | equemene | *
|
107 | 1 | equemene | * PANEL (local input/output) HPL_T_panel *
|
108 | 1 | equemene | * On entry, PANEL points to the data structure containing the
|
109 | 1 | equemene | * panel (to be spread) information.
|
110 | 1 | equemene | *
|
111 | 1 | equemene | * SIDE (global input) const enum HPL_SIDE
|
112 | 1 | equemene | * On entry, SIDE specifies whether the local piece of U located
|
113 | 1 | equemene | * in process IPMAP[SRCDIST] should be spread to the right or to
|
114 | 1 | equemene | * the left. This feature is used by the equilibration process.
|
115 | 1 | equemene | *
|
116 | 1 | equemene | * N (global input) const int
|
117 | 1 | equemene | * On entry, N specifies the local number of columns of U. N
|
118 | 1 | equemene | * must be at least zero.
|
119 | 1 | equemene | *
|
120 | 1 | equemene | * U (local input/output) double *
|
121 | 1 | equemene | * On entry, U is an array of dimension (LDU,*) containing the
|
122 | 1 | equemene | * local pieces of U.
|
123 | 1 | equemene | *
|
124 | 1 | equemene | * LDU (local input) const int
|
125 | 1 | equemene | * On entry, LDU specifies the local leading dimension of U. LDU
|
126 | 1 | equemene | * should be at least MAX(1,IPLEN[nprow]).
|
127 | 1 | equemene | *
|
128 | 1 | equemene | * SRCDIST (local input) const int
|
129 | 1 | equemene | * On entry, SRCDIST specifies the source process that spreads
|
130 | 1 | equemene | * its piece of U.
|
131 | 1 | equemene | *
|
132 | 1 | equemene | * IPLEN (global input) const int *
|
133 | 1 | equemene | * On entry, IPLEN is an array of dimension NPROW+1. This array
|
134 | 1 | equemene | * is such that IPLEN[i+1] - IPLEN[i] is the number of rows of U
|
135 | 1 | equemene | * in each process before process IPMAP[i], with the convention
|
136 | 1 | equemene | * that IPLEN[nprow] is the total number of rows. In other words
|
137 | 1 | equemene | * IPLEN[i+1] - IPLEN[i] is the local number of rows of U that
|
138 | 1 | equemene | * should be moved to process IPMAP[i].
|
139 | 1 | equemene | *
|
140 | 1 | equemene | * IPMAP (global input) const int *
|
141 | 1 | equemene | * On entry, IPMAP is an array of dimension NPROW. This array
|
142 | 1 | equemene | * contains the logarithmic mapping of the processes. In other
|
143 | 1 | equemene | * words, IPMAP[myrow] is the absolute coordinate of the sorted
|
144 | 1 | equemene | * process.
|
145 | 1 | equemene | *
|
146 | 1 | equemene | * IPMAPM1 (global input) const int *
|
147 | 1 | equemene | * On entry, IPMAPM1 is an array of dimension NPROW. This array
|
148 | 1 | equemene | * contains the inverse of the logarithmic mapping contained in
|
149 | 1 | equemene | * IPMAP: For i in [0.. NPROW) IPMAPM1[IPMAP[i]] = i.
|
150 | 1 | equemene | *
|
151 | 1 | equemene | * ---------------------------------------------------------------------
|
152 | 1 | equemene | */
|
153 | 1 | equemene | /*
|
154 | 1 | equemene | * .. Local Variables ..
|
155 | 1 | equemene | */
|
156 | 1 | equemene | MPI_Datatype type; |
157 | 1 | equemene | MPI_Status status; |
158 | 1 | equemene | MPI_Comm comm; |
159 | 1 | equemene | unsigned int ip2=1, mask=1, mydist, mydist2; |
160 | 1 | equemene | int Cmsgid=MSGID_BEGIN_PFACT, ibuf,
|
161 | 1 | equemene | ierr=MPI_SUCCESS, il, k, lbuf, lgth, myrow, |
162 | 1 | equemene | npm1, nprow, partner; |
163 | 1 | equemene | /* ..
|
164 | 1 | equemene | * .. Executable Statements ..
|
165 | 1 | equemene | */
|
166 | 1 | equemene | myrow = PANEL->grid->myrow; nprow = PANEL->grid->nprow; |
167 | 1 | equemene | comm = PANEL->grid->col_comm; |
168 | 1 | equemene | /*
|
169 | 1 | equemene | * Spread U to the left
|
170 | 1 | equemene | */
|
171 | 1 | equemene | if( SIDE == HplLeft )
|
172 | 1 | equemene | { |
173 | 1 | equemene | nprow = ( npm1 = SRCDIST ) + 1;
|
174 | 1 | equemene | if( ( ( mydist = (unsigned int)(IPMAPM1[myrow]) ) > |
175 | 1 | equemene | (unsigned int)(SRCDIST) ) || ( npm1 == 0 ) ) return; |
176 | 1 | equemene | |
177 | 1 | equemene | k = npm1; while( k > 1 ) { k >>= 1; ip2 <<= 1; mask <<= 1; mask++; } |
178 | 1 | equemene | mydist2 = ( mydist = npm1 - mydist ); il = npm1 - ip2; |
179 | 1 | equemene | lgth = IPLEN[nprow]; |
180 | 1 | equemene | |
181 | 1 | equemene | do
|
182 | 1 | equemene | { |
183 | 1 | equemene | mask ^= ip2; |
184 | 1 | equemene | |
185 | 1 | equemene | if( ( mydist & mask ) == 0 ) |
186 | 1 | equemene | { |
187 | 1 | equemene | lbuf = IPLEN[il+1] - ( ibuf = IPLEN[il-Mmin(il, (int)(ip2))] ); |
188 | 1 | equemene | |
189 | 1 | equemene | if( lbuf > 0 ) |
190 | 1 | equemene | { |
191 | 1 | equemene | partner = mydist ^ ip2; |
192 | 1 | equemene | |
193 | 1 | equemene | if( mydist & ip2 )
|
194 | 1 | equemene | { |
195 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
196 | 1 | equemene | ierr = MPI_Type_vector( N, lbuf, LDU, MPI_DOUBLE, |
197 | 1 | equemene | &type ); |
198 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
199 | 1 | equemene | ierr = MPI_Type_commit( &type ); |
200 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
201 | 1 | equemene | ierr = MPI_Recv( Mptr( U, ibuf, 0, LDU ), 1, type, |
202 | 1 | equemene | IPMAP[npm1-partner], Cmsgid, comm, |
203 | 1 | equemene | &status ); |
204 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
205 | 1 | equemene | ierr = MPI_Type_free( &type ); |
206 | 1 | equemene | } |
207 | 1 | equemene | else if( partner < nprow ) |
208 | 1 | equemene | { |
209 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
210 | 1 | equemene | ierr = MPI_Type_vector( N, lbuf, LDU, MPI_DOUBLE, |
211 | 1 | equemene | &type ); |
212 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
213 | 1 | equemene | ierr = MPI_Type_commit( &type ); |
214 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
215 | 1 | equemene | ierr = MPI_Send( Mptr( U, ibuf, 0, LDU ), 1, type, |
216 | 1 | equemene | IPMAP[npm1-partner], Cmsgid, comm ); |
217 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
218 | 1 | equemene | ierr = MPI_Type_free( &type ); |
219 | 1 | equemene | } |
220 | 1 | equemene | } |
221 | 1 | equemene | } |
222 | 1 | equemene | |
223 | 1 | equemene | if( mydist2 < ip2 ) { ip2 >>= 1; il += ip2; } |
224 | 1 | equemene | else { mydist2 -= ip2; ip2 >>= 1; il -= ip2; } |
225 | 1 | equemene | /*
|
226 | 1 | equemene | * Probe for column panel - forward it when available
|
227 | 1 | equemene | */
|
228 | 1 | equemene | if( *IFLAG == HPL_KEEP_TESTING ) (void) HPL_bcast( PBCST, IFLAG ); |
229 | 1 | equemene | |
230 | 1 | equemene | } while( ip2 > 0 ); |
231 | 1 | equemene | } |
232 | 1 | equemene | else
|
233 | 1 | equemene | { |
234 | 1 | equemene | npm1 = ( nprow -= SRCDIST ) - 1;
|
235 | 1 | equemene | if( ( ( mydist = (unsigned int)(IPMAPM1[myrow]) ) < |
236 | 1 | equemene | (unsigned int)(SRCDIST) ) || ( npm1 == 0 ) ) return; |
237 | 1 | equemene | |
238 | 1 | equemene | k = npm1; while( k > 1 ) { k >>= 1; ip2 <<= 1; mask <<= 1; mask++; } |
239 | 1 | equemene | mydist2 = ( mydist -= SRCDIST ); il = ip2; |
240 | 1 | equemene | lgth = IPLEN[SRCDIST+nprow]; |
241 | 1 | equemene | /*
|
242 | 1 | equemene | * Spread U to the right - offset the IPLEN, and IPMAP arrays
|
243 | 1 | equemene | */
|
244 | 1 | equemene | do
|
245 | 1 | equemene | { |
246 | 1 | equemene | mask ^= ip2; |
247 | 1 | equemene | |
248 | 1 | equemene | if( ( mydist & mask ) == 0 ) |
249 | 1 | equemene | { |
250 | 1 | equemene | k = il + ip2; ibuf = IPLEN[SRCDIST+il]; |
251 | 1 | equemene | lbuf = ( k >= nprow ? lgth : IPLEN[SRCDIST+k] ) - ibuf; |
252 | 1 | equemene | |
253 | 1 | equemene | if( lbuf > 0 ) |
254 | 1 | equemene | { |
255 | 1 | equemene | partner = mydist ^ ip2; |
256 | 1 | equemene | |
257 | 1 | equemene | if( mydist & ip2 )
|
258 | 1 | equemene | { |
259 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
260 | 1 | equemene | ierr = MPI_Type_vector( N, lbuf, LDU, MPI_DOUBLE, |
261 | 1 | equemene | &type ); |
262 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
263 | 1 | equemene | ierr = MPI_Type_commit( &type ); |
264 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
265 | 1 | equemene | ierr = MPI_Recv( Mptr( U, ibuf, 0, LDU ), 1, type, |
266 | 1 | equemene | IPMAP[SRCDIST+partner], Cmsgid, |
267 | 1 | equemene | comm, &status ); |
268 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
269 | 1 | equemene | ierr = MPI_Type_free( &type ); |
270 | 1 | equemene | } |
271 | 1 | equemene | else if( partner < nprow ) |
272 | 1 | equemene | { |
273 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
274 | 1 | equemene | ierr = MPI_Type_vector( N, lbuf, LDU, MPI_DOUBLE, |
275 | 1 | equemene | &type ); |
276 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
277 | 1 | equemene | ierr = MPI_Type_commit( &type ); |
278 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
279 | 1 | equemene | ierr = MPI_Send( Mptr( U, ibuf, 0, LDU ), 1, type, |
280 | 1 | equemene | IPMAP[SRCDIST+partner], Cmsgid, |
281 | 1 | equemene | comm ); |
282 | 1 | equemene | if( ierr == MPI_SUCCESS )
|
283 | 1 | equemene | ierr = MPI_Type_free( &type ); |
284 | 1 | equemene | } |
285 | 1 | equemene | } |
286 | 1 | equemene | } |
287 | 1 | equemene | |
288 | 1 | equemene | if( mydist2 < ip2 ) { ip2 >>= 1; il -= ip2; } |
289 | 1 | equemene | else { mydist2 -= ip2; ip2 >>= 1; il += ip2; } |
290 | 1 | equemene | /*
|
291 | 1 | equemene | * Probe for column panel - forward it when available
|
292 | 1 | equemene | */
|
293 | 1 | equemene | if( *IFLAG == HPL_KEEP_TESTING ) (void) HPL_bcast( PBCST, IFLAG ); |
294 | 1 | equemene | |
295 | 1 | equemene | } while( ip2 > 0 ); |
296 | 1 | equemene | } |
297 | 1 | equemene | |
298 | 1 | equemene | if( ierr != MPI_SUCCESS )
|
299 | 1 | equemene | { HPL_pabort( __LINE__, "HPL_spreadN", "MPI call failed" ); } |
300 | 1 | equemene | /*
|
301 | 1 | equemene | * End of HPL_spreadN
|
302 | 1 | equemene | */
|
303 | 1 | equemene | } |