|
ScaLAPACK
2.0.2
ScaLAPACK: Scalable Linear Algebra PACKage
|
00001 #include "redist.h" 00002 /* $Id: pdgemrdrv.c,v 1.1.1.1 2000/02/15 18:04:11 susan Exp $ 00003 * 00004 * pdgemrdrv.c : 00005 * 00006 * 00007 * PURPOSE: 00008 * 00009 * this driver is testing the PDGEMR2D routine. It calls it to obtain a new 00010 * scattered block data decomposition of a distributed DOUBLE PRECISION 00011 * (block scattered) matrix. Then it calls PDGEMR2D for the inverse 00012 * redistribution and checks the results with the initial data. 00013 * 00014 * Data are going from a Block Scattered nbrow0 x nbcol0 decomposition on the 00015 * processor grid p0 x q0, to data distributed in a BS nbrow1 x nbcol1 on the 00016 * processor grid p1 x q1, then back to the BS nbrow0 x nbcol0 decomposition 00017 * on the processor grid p0 x q0. 00018 * 00019 * See pdgemr.c file for detailed info on the PDGEMR2D function. 00020 * 00021 * 00022 * The testing parameters are read from the file GEMR2D.dat, see the file in the 00023 * distribution to have an example. 00024 * 00025 * created by Bernard Tourancheau in April 1994. 00026 * 00027 * modifications : see sccs history 00028 * 00029 * =================================== 00030 * 00031 * 00032 * NOTE : 00033 * 00034 * - the matrix elements are DOUBLE PRECISION 00035 * 00036 * - memory requirements : this procedure requires approximately 3 times the 00037 * memory space of the initial data block in grid 0 (initial block, copy for 00038 * test and second redistribution result) and 1 time the memory space of the 00039 * result data block in grid 1. with the element size = sizeof(double) 00040 * bytes, 00041 * 00042 * 00043 * - use the procedures of the files: 00044 * 00045 * pdgemr.o pdgemr2.o pdgemraux.o 00046 * 00047 * 00048 * ====================================== 00049 * 00050 * WARNING ASSUMPTIONS : 00051 * 00052 * 00053 * ======================================== 00054 * 00055 * 00056 * Planned changes: 00057 * 00058 * 00059 * 00060 * ========================================= */ 00061 #define static2 static 00062 #if defined(Add_) || defined(f77IsF2C) 00063 #define fortran_mr2d pdgemr2do_ 00064 #define fortran_mr2dnew pdgemr2d_ 00065 #elif defined(UpCase) 00066 #define fortran_mr2dnew PDGEMR2D 00067 #define fortran_mr2d PDGEMR2DO 00068 #define dcopy_ DCOPY 00069 #define dlacpy_ DLACPY 00070 #else 00071 #define fortran_mr2d pdgemr2do 00072 #define fortran_mr2dnew pdgemr2d 00073 #define dcopy_ dcopy 00074 #define dlacpy_ dlacpy 00075 #endif 00076 #define Clacpy Cdgelacpy 00077 void Clacpy(); 00078 typedef struct { 00079 int desctype; 00080 int ctxt; 00081 int m; 00082 int n; 00083 int nbrow; 00084 int nbcol; 00085 int sprow; 00086 int spcol; 00087 int lda; 00088 } MDESC; 00089 #define BLOCK_CYCLIC_2D 1 00090 typedef struct { 00091 int lstart; 00092 int len; 00093 } IDESC; 00094 #define SHIFT(row,sprow,nbrow) ((row)-(sprow)+ ((row) >= (sprow) ? 0 : (nbrow))) 00095 #define max(A,B) ((A)>(B)?(A):(B)) 00096 #define min(A,B) ((A)>(B)?(B):(A)) 00097 #define DIVUP(a,b) ( ((a)-1) /(b)+1) 00098 #define ROUNDUP(a,b) (DIVUP(a,b)*(b)) 00099 #ifdef MALLOCDEBUG 00100 #define malloc mymalloc 00101 #define free myfree 00102 #define realloc myrealloc 00103 #endif 00104 /* Cblacs */ 00105 extern void Cblacs_pcoord(); 00106 extern int Cblacs_pnum(); 00107 extern void Csetpvmtids(); 00108 extern void Cblacs_get(); 00109 extern void Cblacs_pinfo(); 00110 extern void Cblacs_gridinfo(); 00111 extern void Cblacs_gridinit(); 00112 extern void Cblacs_exit(); 00113 extern void Cblacs_gridexit(); 00114 extern void Cblacs_setup(); 00115 extern void Cigebs2d(); 00116 extern void Cigebr2d(); 00117 extern void Cigesd2d(); 00118 extern void Cigerv2d(); 00119 extern void Cigsum2d(); 00120 extern void Cigamn2d(); 00121 extern void Cigamx2d(); 00122 extern void Cdgesd2d(); 00123 extern void Cdgerv2d(); 00124 /* lapack */ 00125 void dlacpy_(); 00126 /* aux fonctions */ 00127 extern int localindice(); 00128 extern void *mr2d_malloc(); 00129 extern int ppcm(); 00130 extern int localsize(); 00131 extern int memoryblocksize(); 00132 extern int changeorigin(); 00133 extern void paramcheck(); 00134 /* tools and others function */ 00135 #define scanD0 dgescanD0 00136 #define dispmat dgedispmat 00137 #define setmemory dgesetmemory 00138 #define freememory dgefreememory 00139 #define scan_intervals dgescan_intervals 00140 extern void scanD0(); 00141 extern void dispmat(); 00142 extern void setmemory(); 00143 extern void freememory(); 00144 extern int scan_intervals(); 00145 extern void Cpdgemr2do(); 00146 extern void Cpdgemr2d(); 00147 /* some defines for Cpdgemr2do */ 00148 #define SENDBUFF 0 00149 #define RECVBUFF 1 00150 #define SIZEBUFF 2 00151 #if 0 00152 #define DEBUG 00153 #endif 00154 #ifndef DEBUG 00155 #define NDEBUG 00156 #endif 00157 #include <stdio.h> 00158 #include <stdlib.h> 00159 #include <string.h> 00160 #include <ctype.h> 00161 #include <assert.h> 00162 /* initblock: intialize the local part of a matrix with random data (well, 00163 * not very random) */ 00164 static2 void 00165 initblock(block, m, n) 00166 double *block; 00167 int m, n; 00168 { 00169 double *pdata; 00170 int i; 00171 pdata = block; 00172 for (i = 0; i < m * n; i++, pdata++) { 00173 (*pdata) = i; 00174 }; 00175 } 00176 /* getparam:read from a file a list of integer parameters, the end of the 00177 * parameters to read is given by a NULL at the end of the args list */ 00178 #ifdef __STDC__ 00179 #include <stdarg.h> 00180 static void 00181 getparam(FILE * f,...) 00182 { 00183 #else 00184 #include <varargs.h> 00185 static void 00186 getparam(va_alist) 00187 va_dcl 00188 { 00189 FILE *f; 00190 #endif 00191 va_list ap; 00192 int i; 00193 static int nbline; 00194 char *ptr, *next; 00195 int *var; 00196 static char buffer[200]; 00197 #ifdef __STDC__ 00198 va_start(ap, f); 00199 #else 00200 va_start(ap); 00201 f = va_arg(ap, FILE *); 00202 #endif 00203 do { 00204 next = fgets(buffer, 200, f); 00205 if (next == NULL) { 00206 fprintf(stderr, "bad configuration driver file:after line %d\n", nbline); 00207 exit(1); 00208 } 00209 nbline += 1; 00210 } while (buffer[0] == '#'); 00211 ptr = buffer; 00212 var = va_arg(ap, int *); 00213 while (var != NULL) { 00214 *var = strtol(ptr, &next, 10); 00215 if (ptr == next) { 00216 fprintf(stderr, "bad configuration driver file:error line %d\n", nbline); 00217 exit(1); 00218 } 00219 ptr = next; 00220 var = va_arg(ap, int *); 00221 } 00222 va_end(ap); 00223 } 00224 void 00225 initforpvm(argc, argv) 00226 int argc; 00227 char *argv[]; 00228 { 00229 int pnum, nproc; 00230 Cblacs_pinfo(&pnum, &nproc); 00231 if (nproc < 1) { /* we are with PVM */ 00232 if (pnum == 0) { 00233 if (argc < 2) { 00234 fprintf(stderr, "usage with PVM:xdgemr nbproc\n\ 00235 \t where nbproc is the number of nodes to initialize\n"); 00236 exit(1); 00237 } 00238 nproc = atoi(argv[1]); 00239 } 00240 Cblacs_setup(&pnum, &nproc); 00241 } 00242 } 00243 int 00244 main(argc, argv) 00245 int argc; 00246 char *argv[]; 00247 { 00248 /* We initialize the data-block on the current processor, then redistribute 00249 * it, and perform the inverse redistribution to compare the local memory 00250 * with the initial one. */ 00251 /* Data file */ 00252 FILE *fp; 00253 int nbre, nbremax; 00254 /* Data distribution 0 parameters */ 00255 int p0, /* # of rows in the processor grid */ 00256 q0; /* # of columns in the processor grid */ 00257 /* Data distribution 1 parameters */ 00258 int p1, q1; 00259 /* # of parameter to be read on the keyboard */ 00260 #define nbparameter 24 00261 /* General variables */ 00262 int blocksize0; 00263 int mypnum, nprocs; 00264 int parameters[nbparameter], nberrors; 00265 int i; 00266 int ia, ja, ib, jb, m, n; 00267 int gcontext, context0, context1; 00268 int myprow1, myprow0, mypcol0, mypcol1; 00269 int dummy; 00270 MDESC ma, mb; 00271 double *ptrmyblock, *ptrsavemyblock, *ptrmyblockcopy, *ptrmyblockvide; 00272 #ifdef UsingMpiBlacs 00273 MPI_Init(&argc, &argv); 00274 #endif 00275 setvbuf(stdout, NULL, _IOLBF, 0); 00276 setvbuf(stderr, NULL, _IOLBF, 0); 00277 #ifdef T3D 00278 free(malloc(14000000)); 00279 #endif 00280 initforpvm(argc, argv); 00281 /* Read physical parameters */ 00282 Cblacs_pinfo(&mypnum, &nprocs); 00283 /* initialize BLACS for the parameter communication */ 00284 Cblacs_get(0, 0, &gcontext); 00285 Cblacs_gridinit(&gcontext, "R", nprocs, 1); 00286 Cblacs_gridinfo(gcontext, &dummy, &dummy, &mypnum, &dummy); 00287 if (mypnum == 0) { 00288 if ((fp = fopen("GEMR2D.dat", "r")) == NULL) { 00289 fprintf(stderr, "Can't open GEMR2D.dat\n"); 00290 exit(1); 00291 }; 00292 printf("\n// DGEMR2D TESTER for DOUBLE PRECISION //\n"); 00293 getparam(fp, &nbre, NULL); 00294 printf("////////// %d tests \n\n", nbre); 00295 parameters[0] = nbre; 00296 Cigebs2d(gcontext, "All", "H", 1, 1, parameters, 1); 00297 } else { 00298 Cigebr2d(gcontext, "All", "H", 1, 1, parameters, 1, 0, 0); 00299 nbre = parameters[0]; 00300 }; 00301 if (mypnum == 0) { 00302 printf("\n m n m0 n0 sr0 sc0 i0 j0 p0 q0 nbr0 nbc0 \ 00303 m1 n1 sr1 sc1 i1 j1 p1 q1 nbr1 nbc1\n\n"); 00304 }; 00305 /****** TEST LOOP *****/ 00306 /* Here we are in grip 1xnprocs */ 00307 nbremax = nbre; 00308 #ifdef DEBUG 00309 fprintf(stderr, "bonjour,je suis le noeud %d\n", mypnum); 00310 #endif 00311 while (nbre-- != 0) { /* Loop on the serie of tests */ 00312 /* All the processors read the parameters so we have to be in a 1xnprocs 00313 * grid at each iteration */ 00314 /* Read processors grid and matrices parameters */ 00315 if (mypnum == 0) { 00316 int u, d; 00317 getparam(fp, 00318 &m, &n, 00319 &ma.m, &ma.n, &ma.sprow, &ma.spcol, 00320 &ia, &ja, &p0, &q0, &ma.nbrow, &ma.nbcol, 00321 &mb.m, &mb.n, &mb.sprow, &mb.spcol, 00322 &ib, &jb, &p1, &q1, &mb.nbrow, &mb.nbcol, 00323 NULL); 00324 printf("\t\t************* TEST # %d **********\n", 00325 nbremax - nbre); 00326 printf(" %3d %3d %3d %3d %3d %3d %3d %3d \ 00327 %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d", 00328 m, n, 00329 ma.m, ma.n, ma.sprow, ma.spcol, 00330 ia, ja, p0, q0, ma.nbrow, ma.nbcol, 00331 mb.m, mb.n, mb.sprow, mb.spcol, 00332 ib, jb, p1, q1, mb.nbrow, mb.nbcol); 00333 printf("\n"); 00334 if (p0 * q0 > nprocs || p1 * q1 > nprocs) { 00335 fprintf(stderr, "not enough nodes:%d processors required\n", 00336 max(p0 * q0, p1 * q1)); 00337 exit(1); 00338 } 00339 parameters[0] = p0; 00340 parameters[1] = q0; 00341 parameters[2] = ma.nbrow; 00342 parameters[3] = ma.nbcol; 00343 parameters[4] = p1; 00344 parameters[5] = q1; 00345 parameters[6] = mb.nbrow; 00346 parameters[7] = mb.nbcol; 00347 parameters[8] = ma.m; 00348 parameters[9] = ma.n; 00349 parameters[10] = ma.sprow; 00350 parameters[11] = ma.spcol; 00351 parameters[12] = mb.sprow; 00352 parameters[13] = mb.spcol; 00353 parameters[14] = ia; 00354 parameters[15] = ja; 00355 parameters[16] = ib; 00356 parameters[17] = jb; 00357 parameters[18] = m; 00358 parameters[19] = n; 00359 parameters[20] = mb.m; 00360 parameters[21] = mb.n; 00361 Cigebs2d(gcontext, "All", "H", 1, nbparameter, parameters, 1); 00362 } else { 00363 Cigebr2d(gcontext, "All", "H", 1, nbparameter, parameters, 1, 0, 0); 00364 p0 = parameters[0]; 00365 q0 = parameters[1]; 00366 ma.nbrow = parameters[2]; 00367 ma.nbcol = parameters[3]; 00368 p1 = parameters[4]; 00369 q1 = parameters[5]; 00370 mb.nbrow = parameters[6]; 00371 mb.nbcol = parameters[7]; 00372 ma.m = parameters[8]; 00373 ma.n = parameters[9]; 00374 ma.sprow = parameters[10]; 00375 ma.spcol = parameters[11]; 00376 mb.sprow = parameters[12]; 00377 mb.spcol = parameters[13]; 00378 ia = parameters[14]; 00379 ja = parameters[15]; 00380 ib = parameters[16]; 00381 jb = parameters[17]; 00382 m = parameters[18]; 00383 n = parameters[19]; 00384 mb.m = parameters[20]; 00385 mb.n = parameters[21]; 00386 ma.desctype = BLOCK_CYCLIC_2D; 00387 mb.desctype = BLOCK_CYCLIC_2D; 00388 }; 00389 Cblacs_get(0, 0, &context0); 00390 Cblacs_gridinit(&context0, "R", p0, q0); 00391 Cblacs_get(0, 0, &context1); 00392 Cblacs_gridinit(&context1, "R", p1, q1); 00393 Cblacs_gridinfo(context0, &dummy, &dummy, &myprow0, &mypcol0); 00394 if (myprow0 >= p0 || mypcol0 >= q0) 00395 myprow0 = mypcol0 = -1; 00396 Cblacs_gridinfo(context1, &dummy, &dummy, &myprow1, &mypcol1); 00397 if (myprow1 >= p1 || mypcol1 >= q1) 00398 myprow1 = mypcol1 = -1; 00399 assert((myprow0 < p0 && mypcol0 < q0) || (myprow0 == -1 && mypcol0 == -1)); 00400 assert((myprow1 < p1 && mypcol1 < q1) || (myprow1 == -1 && mypcol1 == -1)); 00401 ma.ctxt = context0; 00402 mb.ctxt = context1; 00403 /* From here, we are not assuming that only the processors working in the 00404 * redistribution are calling xxMR2D, but the ones not concerned will do 00405 * nothing. */ 00406 /* We compute the exact size of the local memory block for the memory 00407 * allocations */ 00408 if (myprow0 >= 0 && mypcol0 >= 0) { 00409 blocksize0 = memoryblocksize(&ma); 00410 ma.lda = localsize(SHIFT(myprow0, ma.sprow, p0), p0, ma.nbrow, ma.m); 00411 setmemory(&ptrmyblock, blocksize0); 00412 initblock(ptrmyblock, 1, blocksize0); 00413 setmemory(&ptrmyblockcopy, blocksize0); 00414 memcpy((char *) ptrmyblockcopy, (char *) ptrmyblock, 00415 blocksize0 * sizeof(double)); 00416 setmemory(&ptrmyblockvide, blocksize0); 00417 for (i = 0; i < blocksize0; i++) 00418 ptrmyblockvide[i] = -1; 00419 }; /* if (mypnum < p0 * q0) */ 00420 if (myprow1 >= 0 && mypcol1 >= 0) { 00421 setmemory(&ptrsavemyblock, memoryblocksize(&mb)); 00422 mb.lda = localsize(SHIFT(myprow1, mb.sprow, p1), p1, mb.nbrow, mb.m); 00423 }; /* if (mypnum < p1 * q1) */ 00424 /* Redistribute the matrix from grid 0 to grid 1 (memory location 00425 * ptrmyblock to ptrsavemyblock) */ 00426 Cpdgemr2d(m, n, 00427 ptrmyblock, ia, ja, &ma, 00428 ptrsavemyblock, ib, jb, &mb, gcontext); 00429 /* Perform the inverse redistribution of the matrix from grid 1 to grid 0 00430 * (memory location ptrsavemyblock to ptrmyblockvide) */ 00431 Cpdgemr2d(m, n, 00432 ptrsavemyblock, ib, jb, &mb, 00433 ptrmyblockvide, ia, ja, &ma, gcontext); 00434 /* Check the differences */ 00435 nberrors = 0; 00436 if (myprow0 >= 0 && mypcol0 >= 0) { 00437 /* only for the processors that do have data at the begining */ 00438 for (i = 0; i < blocksize0; i++) { 00439 int li, lj, gi, gj; 00440 int in; 00441 in = 1; 00442 li = i % ma.lda; 00443 lj = i / ma.lda; 00444 gi = (li / ma.nbrow) * p0 * ma.nbrow + 00445 SHIFT(myprow0, ma.sprow, p0) * ma.nbrow + li % ma.nbrow; 00446 gj = (lj / ma.nbcol) * q0 * ma.nbcol + 00447 SHIFT(mypcol0, ma.spcol, q0) * ma.nbcol + lj % ma.nbcol; 00448 assert(gi < ma.m && gj < ma.n); 00449 gi -= (ia - 1); 00450 gj -= (ja - 1); 00451 if (gi < 0 || gj < 0 || gi >= m || gj >= n) 00452 in = 0; 00453 if (!in) { 00454 ptrmyblockcopy[i] = -1; 00455 } 00456 if (ptrmyblockvide[i] != ptrmyblockcopy[i]) { 00457 nberrors++; 00458 printf("Proc %d : Error element number %d, value = %f , initvalue =%f \n" 00459 ,mypnum, i, 00460 ptrmyblockvide[i], ptrmyblockcopy[i]); 00461 }; 00462 }; 00463 if (nberrors > 0) { 00464 printf("Processor %d, has tested %d DOUBLE PRECISION elements,\ 00465 Number of redistribution errors = %d \n", 00466 mypnum, blocksize0, nberrors); 00467 } 00468 } 00469 /* Look at the errors on all the processors at this point. */ 00470 Cigsum2d(gcontext, "All", "H", 1, 1, &nberrors, 1, 0, 0); 00471 if (mypnum == 0) 00472 if (nberrors) 00473 printf(" => Total number of redistribution errors = %d \n", 00474 nberrors); 00475 else 00476 printf("TEST PASSED OK\n"); 00477 /* release memory for the next iteration */ 00478 if (myprow0 >= 0 && mypcol0 >= 0) { 00479 freememory((char *) ptrmyblock); 00480 freememory((char *) ptrmyblockvide); 00481 freememory((char *) ptrmyblockcopy); 00482 }; /* if (mypnum < p0 * q0) */ 00483 /* release memory for the next iteration */ 00484 if (myprow1 >= 0 && mypcol1 >= 0) { 00485 freememory((char *) ptrsavemyblock); 00486 }; 00487 if (myprow0 >= 0) 00488 Cblacs_gridexit(context0); 00489 if (myprow1 >= 0) 00490 Cblacs_gridexit(context1); 00491 }; /* while nbre != 0 */ 00492 if (mypnum == 0) { 00493 fclose(fp); 00494 }; 00495 Cblacs_exit(0); 00496 return 0; 00497 }/* main */