1 #define PETSCKSP_DLL 2 3 /*************************************xyt.c************************************ 4 Module Name: xyt 5 Module Info: 6 7 author: Henry M. Tufo III 8 e-mail: hmt@asci.uchicago.edu 9 contact: 10 +--------------------------------+--------------------------------+ 11 |MCS Division - Building 221 |Department of Computer Science | 12 |Argonne National Laboratory |Ryerson 152 | 13 |9700 S. Cass Avenue |The University of Chicago | 14 |Argonne, IL 60439 |Chicago, IL 60637 | 15 |(630) 252-5354/5986 ph/fx |(773) 702-6019/8487 ph/fx | 16 +--------------------------------+--------------------------------+ 17 18 Last Modification: 3.20.01 19 **************************************xyt.c***********************************/ 20 #include "src/ksp/pc/impls/tfs/tfs.h" 21 22 #define LEFT -1 23 #define RIGHT 1 24 #define BOTH 0 25 26 typedef struct xyt_solver_info { 27 PetscInt n, m, n_global, m_global; 28 PetscInt nnz, max_nnz, msg_buf_sz; 29 PetscInt *nsep, *lnsep, *fo, nfo, *stages; 30 PetscInt *xcol_sz, *xcol_indices; 31 PetscScalar **xcol_vals, *x, *solve_uu, *solve_w; 32 PetscInt *ycol_sz, *ycol_indices; 33 PetscScalar **ycol_vals, *y; 34 PetscInt nsolves; 35 PetscScalar tot_solve_time; 36 } xyt_info; 37 38 39 typedef struct matvec_info { 40 PetscInt n, m, n_global, m_global; 41 PetscInt *local2global; 42 gs_ADT gs_handle; 43 PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*); 44 void *grid_data; 45 } mv_info; 46 47 struct xyt_CDT{ 48 PetscInt id; 49 PetscInt ns; 50 PetscInt level; 51 xyt_info *info; 52 mv_info *mvi; 53 }; 54 55 static PetscInt n_xyt=0; 56 static PetscInt n_xyt_handles=0; 57 58 /* prototypes */ 59 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs); 60 static PetscErrorCode check_handle(xyt_ADT xyt_handle); 61 static PetscErrorCode det_separators(xyt_ADT xyt_handle); 62 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 63 static PetscInt xyt_generate(xyt_ADT xyt_handle); 64 static PetscInt do_xyt_factor(xyt_ADT xyt_handle); 65 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data); 66 67 /**************************************xyt.c***********************************/ 68 xyt_ADT XYT_new(void) 69 { 70 xyt_ADT xyt_handle; 71 72 /* rolling count on n_xyt ... pot. problem here */ 73 n_xyt_handles++; 74 xyt_handle = (xyt_ADT)malloc(sizeof(struct xyt_CDT)); 75 xyt_handle->id = ++n_xyt; 76 xyt_handle->info = NULL; 77 xyt_handle->mvi = NULL; 78 79 return(xyt_handle); 80 } 81 82 /**************************************xyt.c***********************************/ 83 PetscInt XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt handle */ 84 PetscInt *local2global, /* global column mapping */ 85 PetscInt n, /* local num rows */ 86 PetscInt m, /* local num cols */ 87 void *matvec, /* b_loc=A_local.x_loc */ 88 void *grid_data /* grid data for matvec */ 89 ) 90 { 91 92 comm_init(); 93 check_handle(xyt_handle); 94 95 /* only 2^k for now and all nodes participating */ 96 if ((1<<(xyt_handle->level=i_log2_num_nodes))!=num_nodes) 97 {SETERRQ2(PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %d != %d\n",1<<i_log2_num_nodes,num_nodes);} 98 99 /* space for X info */ 100 xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info)); 101 102 /* set up matvec handles */ 103 xyt_handle->mvi = set_mvi(local2global, n, m, matvec, grid_data); 104 105 /* matrix is assumed to be of full rank */ 106 /* LATER we can reset to indicate rank def. */ 107 xyt_handle->ns=0; 108 109 /* determine separators and generate firing order - NB xyt info set here */ 110 det_separators(xyt_handle); 111 112 return(do_xyt_factor(xyt_handle)); 113 } 114 115 /**************************************xyt.c***********************************/ 116 PetscInt XYT_solve(xyt_ADT xyt_handle, double *x, double *b) 117 { 118 comm_init(); 119 check_handle(xyt_handle); 120 121 /* need to copy b into x? */ 122 if (b) 123 {rvec_copy(x,b,xyt_handle->mvi->n);} 124 do_xyt_solve(xyt_handle,x); 125 126 return(0); 127 } 128 129 /**************************************xyt.c***********************************/ 130 PetscInt XYT_free(xyt_ADT xyt_handle) 131 { 132 comm_init(); 133 check_handle(xyt_handle); 134 n_xyt_handles--; 135 136 free(xyt_handle->info->nsep); 137 free(xyt_handle->info->lnsep); 138 free(xyt_handle->info->fo); 139 free(xyt_handle->info->stages); 140 free(xyt_handle->info->solve_uu); 141 free(xyt_handle->info->solve_w); 142 free(xyt_handle->info->x); 143 free(xyt_handle->info->xcol_vals); 144 free(xyt_handle->info->xcol_sz); 145 free(xyt_handle->info->xcol_indices); 146 free(xyt_handle->info->y); 147 free(xyt_handle->info->ycol_vals); 148 free(xyt_handle->info->ycol_sz); 149 free(xyt_handle->info->ycol_indices); 150 free(xyt_handle->info); 151 free(xyt_handle->mvi->local2global); 152 gs_free(xyt_handle->mvi->gs_handle); 153 free(xyt_handle->mvi); 154 free(xyt_handle); 155 156 157 /* if the check fails we nuke */ 158 /* if NULL pointer passed to free we nuke */ 159 /* if the calls to free fail that's not my problem */ 160 return(0); 161 } 162 163 /**************************************xyt.c***********************************/ 164 PetscInt XYT_stats(xyt_ADT xyt_handle) 165 { 166 PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 167 PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 168 PetscInt vals[9], work[9]; 169 PetscScalar fvals[3], fwork[3]; 170 171 172 comm_init(); 173 check_handle(xyt_handle); 174 175 /* if factorization not done there are no stats */ 176 if (!xyt_handle->info||!xyt_handle->mvi) 177 { 178 if (!my_id) 179 {printf("XYT_stats() :: no stats available!\n");} 180 return 1; 181 } 182 183 vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz; 184 vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n; 185 vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz; 186 giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 187 188 fvals[0]=fvals[1]=fvals[2] 189 =xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++; 190 grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop); 191 192 if (!my_id) 193 { 194 printf("%d :: min xyt_nnz=%d\n",my_id,vals[0]); 195 printf("%d :: max xyt_nnz=%d\n",my_id,vals[1]); 196 printf("%d :: avg xyt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes); 197 printf("%d :: tot xyt_nnz=%d\n",my_id,vals[2]); 198 printf("%d :: xyt C(2d) =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5))); 199 printf("%d :: xyt C(3d) =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667))); 200 printf("%d :: min xyt_n =%d\n",my_id,vals[3]); 201 printf("%d :: max xyt_n =%d\n",my_id,vals[4]); 202 printf("%d :: avg xyt_n =%g\n",my_id,1.0*vals[5]/num_nodes); 203 printf("%d :: tot xyt_n =%d\n",my_id,vals[5]); 204 printf("%d :: min xyt_buf=%d\n",my_id,vals[6]); 205 printf("%d :: max xyt_buf=%d\n",my_id,vals[7]); 206 printf("%d :: avg xyt_buf=%g\n",my_id,1.0*vals[8]/num_nodes); 207 printf("%d :: min xyt_slv=%g\n",my_id,fvals[0]); 208 printf("%d :: max xyt_slv=%g\n",my_id,fvals[1]); 209 printf("%d :: avg xyt_slv=%g\n",my_id,fvals[2]/num_nodes); 210 } 211 212 return(0); 213 } 214 215 216 /*************************************xyt.c************************************ 217 218 Description: get A_local, local portion of global coarse matrix which 219 is a row dist. nxm matrix w/ n<m. 220 o my_ml holds address of ML struct associated w/A_local and coarse grid 221 o local2global holds global number of column i (i=0,...,m-1) 222 o local2global holds global number of row i (i=0,...,n-1) 223 o mylocmatvec performs A_local . vec_local (note that gs is performed using 224 gs_init/gop). 225 226 mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 227 mylocmatvec (void :: void *data, double *in, double *out) 228 **************************************xyt.c***********************************/ 229 static PetscInt do_xyt_factor(xyt_ADT xyt_handle) 230 { 231 return xyt_generate(xyt_handle); 232 } 233 234 /**************************************xyt.c***********************************/ 235 static PetscInt xyt_generate(xyt_ADT xyt_handle) 236 { 237 PetscInt i,j,k,idx; 238 PetscInt dim, col; 239 PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 240 PetscInt *segs; 241 PetscInt op[] = {GL_ADD,0}; 242 PetscInt off, len; 243 PetscScalar *x_ptr, *y_ptr; 244 PetscInt *iptr, flag; 245 PetscInt start=0, end, work; 246 PetscInt op2[] = {GL_MIN,0}; 247 gs_ADT gs_handle; 248 PetscInt *nsep, *lnsep, *fo; 249 PetscInt a_n=xyt_handle->mvi->n; 250 PetscInt a_m=xyt_handle->mvi->m; 251 PetscInt *a_local2global=xyt_handle->mvi->local2global; 252 PetscInt level; 253 PetscInt n, m; 254 PetscInt *xcol_sz, *xcol_indices, *stages; 255 PetscScalar **xcol_vals, *x; 256 PetscInt *ycol_sz, *ycol_indices; 257 PetscScalar **ycol_vals, *y; 258 PetscInt n_global; 259 PetscInt xt_nnz=0, xt_max_nnz=0; 260 PetscInt yt_nnz=0, yt_max_nnz=0; 261 PetscInt xt_zero_nnz =0; 262 PetscInt xt_zero_nnz_0=0; 263 PetscInt yt_zero_nnz =0; 264 PetscInt yt_zero_nnz_0=0; 265 PetscBLASInt i1 = 1; 266 PetscScalar dm1 = -1.0; 267 PetscErrorCode ierr; 268 269 n=xyt_handle->mvi->n; 270 nsep=xyt_handle->info->nsep; 271 lnsep=xyt_handle->info->lnsep; 272 fo=xyt_handle->info->fo; 273 end=lnsep[0]; 274 level=xyt_handle->level; 275 gs_handle=xyt_handle->mvi->gs_handle; 276 277 /* is there a null space? */ 278 /* LATER add in ability to detect null space by checking alpha */ 279 for (i=0, j=0; i<=level; i++) 280 {j+=nsep[i];} 281 282 m = j-xyt_handle->ns; 283 if (m!=j) 284 {printf("xyt_generate() :: null space exists %d %d %d\n",m,j,xyt_handle->ns);} 285 286 ierr = PetscInfo2(0,"xyt_generate() :: X(%d,%d)\n",n,m);CHKERRQ(ierr); 287 288 /* get and initialize storage for x local */ 289 /* note that x local is nxm and stored by columns */ 290 xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 291 xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 292 xcol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *)); 293 for (i=j=0; i<m; i++, j+=2) 294 { 295 xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1; 296 xcol_vals[i] = NULL; 297 } 298 xcol_indices[j]=-1; 299 300 /* get and initialize storage for y local */ 301 /* note that y local is nxm and stored by columns */ 302 ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 303 ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 304 ycol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *)); 305 for (i=j=0; i<m; i++, j+=2) 306 { 307 ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1; 308 ycol_vals[i] = NULL; 309 } 310 ycol_indices[j]=-1; 311 312 /* size of separators for each sub-hc working from bottom of tree to top */ 313 /* this looks like nsep[]=segments */ 314 stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 315 segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 316 ivec_zero(stages,level+1); 317 ivec_copy(segs,nsep,level+1); 318 for (i=0; i<level; i++) 319 {segs[i+1] += segs[i];} 320 stages[0] = segs[0]; 321 322 /* temporary vectors */ 323 u = (PetscScalar *) malloc(n*sizeof(PetscScalar)); 324 z = (PetscScalar *) malloc(n*sizeof(PetscScalar)); 325 v = (PetscScalar *) malloc(a_m*sizeof(PetscScalar)); 326 uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 327 w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 328 329 /* extra nnz due to replication of vertices across separators */ 330 for (i=1, j=0; i<=level; i++) 331 {j+=nsep[i];} 332 333 /* storage for sparse x values */ 334 n_global = xyt_handle->info->n_global; 335 xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes; 336 x = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar)); 337 y = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar)); 338 339 /* LATER - can embed next sep to fire in gs */ 340 /* time to make the donuts - generate X factor */ 341 for (dim=i=j=0;i<m;i++) 342 { 343 /* time to move to the next level? */ 344 while (i==segs[dim]) 345 { 346 if (dim==level) 347 {SETERRQ(PETSC_ERR_PLIB,"dim about to exceed level\n"); break;} 348 349 stages[dim++]=i; 350 end+=lnsep[dim]; 351 } 352 stages[dim]=i; 353 354 /* which column are we firing? */ 355 /* i.e. set v_l */ 356 /* use new seps and do global min across hc to determine which one to fire */ 357 (start<end) ? (col=fo[start]) : (col=INT_MAX); 358 giop_hc(&col,&work,1,op2,dim); 359 360 /* shouldn't need this */ 361 if (col==INT_MAX) 362 { 363 ierr = PetscInfo(0,"hey ... col==INT_MAX??\n"); 364 continue; 365 } 366 367 /* do I own it? I should */ 368 rvec_zero(v ,a_m); 369 if (col==fo[start]) 370 { 371 start++; 372 idx=ivec_linear_search(col, a_local2global, a_n); 373 if (idx!=-1) 374 {v[idx] = 1.0; j++;} 375 else 376 {SETERRQ(PETSC_ERR_PLIB,"NOT FOUND!\n");} 377 } 378 else 379 { 380 idx=ivec_linear_search(col, a_local2global, a_m); 381 if (idx!=-1) 382 {v[idx] = 1.0;} 383 } 384 385 /* perform u = A.v_l */ 386 rvec_zero(u,n); 387 do_matvec(xyt_handle->mvi,v,u); 388 389 /* uu = X^T.u_l (local portion) */ 390 /* technically only need to zero out first i entries */ 391 /* later turn this into an XYT_solve call ? */ 392 rvec_zero(uu,m); 393 y_ptr=y; 394 iptr = ycol_indices; 395 for (k=0; k<i; k++) 396 { 397 off = *iptr++; 398 len = *iptr++; 399 400 uu[k] = BLASdot_(&len,u+off,&i1,y_ptr,&i1); 401 y_ptr+=len; 402 } 403 404 /* uu = X^T.u_l (comm portion) */ 405 ssgl_radd (uu, w, dim, stages); 406 407 /* z = X.uu */ 408 rvec_zero(z,n); 409 x_ptr=x; 410 iptr = xcol_indices; 411 for (k=0; k<i; k++) 412 { 413 off = *iptr++; 414 len = *iptr++; 415 416 BLASaxpy_(&len,&uu[k],x_ptr,&i1,z+off,&i1); 417 x_ptr+=len; 418 } 419 420 /* compute v_l = v_l - z */ 421 rvec_zero(v+a_n,a_m-a_n); 422 BLASaxpy_(&n,&dm1,z,&i1,v,&i1); 423 424 /* compute u_l = A.v_l */ 425 if (a_n!=a_m) 426 {gs_gop_hc(gs_handle,v,"+\0",dim);} 427 rvec_zero(u,n); 428 do_matvec(xyt_handle->mvi,v,u); 429 430 /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 431 alpha = BLASdot_(&n,u,&i1,u,&i1); 432 /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 433 grop_hc(&alpha, &alpha_w, 1, op, dim); 434 435 alpha = (PetscScalar) sqrt((double)alpha); 436 437 /* check for small alpha */ 438 /* LATER use this to detect and determine null space */ 439 if (fabs(alpha)<1.0e-14) 440 {SETERRQ1(PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);} 441 442 /* compute v_l = v_l/sqrt(alpha) */ 443 rvec_scale(v,1.0/alpha,n); 444 rvec_scale(u,1.0/alpha,n); 445 446 /* add newly generated column, v_l, to X */ 447 flag = 1; 448 off=len=0; 449 for (k=0; k<n; k++) 450 { 451 if (v[k]!=0.0) 452 { 453 len=k; 454 if (flag) 455 {off=k; flag=0;} 456 } 457 } 458 459 len -= (off-1); 460 461 if (len>0) 462 { 463 if ((xt_nnz+len)>xt_max_nnz) 464 { 465 ierr = PetscInfo(0,"increasing space for X by 2x!\n"); 466 xt_max_nnz *= 2; 467 x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar)); 468 rvec_copy(x_ptr,x,xt_nnz); 469 free(x); 470 x = x_ptr; 471 x_ptr+=xt_nnz; 472 } 473 xt_nnz += len; 474 rvec_copy(x_ptr,v+off,len); 475 476 /* keep track of number of zeros */ 477 if (dim) 478 { 479 for (k=0; k<len; k++) 480 { 481 if (x_ptr[k]==0.0) 482 {xt_zero_nnz++;} 483 } 484 } 485 else 486 { 487 for (k=0; k<len; k++) 488 { 489 if (x_ptr[k]==0.0) 490 {xt_zero_nnz_0++;} 491 } 492 } 493 xcol_indices[2*i] = off; 494 xcol_sz[i] = xcol_indices[2*i+1] = len; 495 xcol_vals[i] = x_ptr; 496 } 497 else 498 { 499 xcol_indices[2*i] = 0; 500 xcol_sz[i] = xcol_indices[2*i+1] = 0; 501 xcol_vals[i] = x_ptr; 502 } 503 504 505 /* add newly generated column, u_l, to Y */ 506 flag = 1; 507 off=len=0; 508 for (k=0; k<n; k++) 509 { 510 if (u[k]!=0.0) 511 { 512 len=k; 513 if (flag) 514 {off=k; flag=0;} 515 } 516 } 517 518 len -= (off-1); 519 520 if (len>0) 521 { 522 if ((yt_nnz+len)>yt_max_nnz) 523 { 524 ierr = PetscInfo(0,"increasing space for Y by 2x!\n"); 525 yt_max_nnz *= 2; 526 y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar)); 527 rvec_copy(y_ptr,y,yt_nnz); 528 free(y); 529 y = y_ptr; 530 y_ptr+=yt_nnz; 531 } 532 yt_nnz += len; 533 rvec_copy(y_ptr,u+off,len); 534 535 /* keep track of number of zeros */ 536 if (dim) 537 { 538 for (k=0; k<len; k++) 539 { 540 if (y_ptr[k]==0.0) 541 {yt_zero_nnz++;} 542 } 543 } 544 else 545 { 546 for (k=0; k<len; k++) 547 { 548 if (y_ptr[k]==0.0) 549 {yt_zero_nnz_0++;} 550 } 551 } 552 ycol_indices[2*i] = off; 553 ycol_sz[i] = ycol_indices[2*i+1] = len; 554 ycol_vals[i] = y_ptr; 555 } 556 else 557 { 558 ycol_indices[2*i] = 0; 559 ycol_sz[i] = ycol_indices[2*i+1] = 0; 560 ycol_vals[i] = y_ptr; 561 } 562 } 563 564 /* close off stages for execution phase */ 565 while (dim!=level) 566 { 567 stages[dim++]=i; 568 ierr = PetscInfo2(0,"disconnected!!! dim(%d)!=level(%d)\n",dim,level);CHKERRQ(ierr); 569 } 570 stages[dim]=i; 571 572 xyt_handle->info->n=xyt_handle->mvi->n; 573 xyt_handle->info->m=m; 574 xyt_handle->info->nnz=xt_nnz + yt_nnz; 575 xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz; 576 xyt_handle->info->msg_buf_sz=stages[level]-stages[0]; 577 xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 578 xyt_handle->info->solve_w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 579 xyt_handle->info->x=x; 580 xyt_handle->info->xcol_vals=xcol_vals; 581 xyt_handle->info->xcol_sz=xcol_sz; 582 xyt_handle->info->xcol_indices=xcol_indices; 583 xyt_handle->info->stages=stages; 584 xyt_handle->info->y=y; 585 xyt_handle->info->ycol_vals=ycol_vals; 586 xyt_handle->info->ycol_sz=ycol_sz; 587 xyt_handle->info->ycol_indices=ycol_indices; 588 589 free(segs); 590 free(u); 591 free(v); 592 free(uu); 593 free(z); 594 free(w); 595 596 return(0); 597 } 598 599 /**************************************xyt.c***********************************/ 600 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 601 { 602 PetscInt off, len, *iptr; 603 PetscInt level =xyt_handle->level; 604 PetscInt n =xyt_handle->info->n; 605 PetscInt m =xyt_handle->info->m; 606 PetscInt *stages =xyt_handle->info->stages; 607 PetscInt *xcol_indices=xyt_handle->info->xcol_indices; 608 PetscInt *ycol_indices=xyt_handle->info->ycol_indices; 609 PetscScalar *x_ptr, *y_ptr, *uu_ptr; 610 PetscScalar *solve_uu=xyt_handle->info->solve_uu; 611 PetscScalar *solve_w =xyt_handle->info->solve_w; 612 PetscScalar *x =xyt_handle->info->x; 613 PetscScalar *y =xyt_handle->info->y; 614 PetscBLASInt i1 = 1; 615 616 PetscFunctionBegin; 617 uu_ptr=solve_uu; 618 rvec_zero(uu_ptr,m); 619 620 /* x = X.Y^T.b */ 621 /* uu = Y^T.b */ 622 for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) 623 { 624 off=*iptr++; len=*iptr++; 625 *uu_ptr++ = BLASdot_(&len,uc+off,&i1,y_ptr,&i1); 626 } 627 628 /* comunication of beta */ 629 uu_ptr=solve_uu; 630 if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);} 631 632 rvec_zero(uc,n); 633 634 /* x = X.uu */ 635 for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) 636 { 637 off=*iptr++; len=*iptr++; 638 BLASaxpy_(&len,uu_ptr++,x_ptr,&i1,uc+off,&i1); 639 } 640 PetscFunctionReturn(0); 641 } 642 643 /**************************************xyt.c***********************************/ 644 static PetscErrorCode check_handle(xyt_ADT xyt_handle) 645 { 646 PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 647 648 PetscFunctionBegin; 649 if (xyt_handle==NULL) 650 {SETERRQ1(PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %d\n",xyt_handle);} 651 652 vals[0]=vals[1]=xyt_handle->id; 653 giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 654 if ((vals[0]!=vals[1])||(xyt_handle->id<=0)) 655 {SETERRQ3(PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %d/%d %d\n", vals[0],vals[1], xyt_handle->id);} 656 PetscFunctionReturn(0); 657 } 658 659 /**************************************xyt.c***********************************/ 660 static PetscErrorCode det_separators(xyt_ADT xyt_handle) 661 { 662 PetscInt i, ct, id; 663 PetscInt mask, edge, *iptr; 664 PetscInt *dir, *used; 665 PetscInt sum[4], w[4]; 666 PetscScalar rsum[4], rw[4]; 667 PetscInt op[] = {GL_ADD,0}; 668 PetscScalar *lhs, *rhs; 669 PetscInt *nsep, *lnsep, *fo, nfo=0; 670 gs_ADT gs_handle=xyt_handle->mvi->gs_handle; 671 PetscInt *local2global=xyt_handle->mvi->local2global; 672 PetscInt n=xyt_handle->mvi->n; 673 PetscInt m=xyt_handle->mvi->m; 674 PetscInt level=xyt_handle->level; 675 PetscInt shared=FALSE; 676 PetscErrorCode ierr; 677 678 PetscFunctionBegin; 679 dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 680 nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 681 lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 682 fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 683 used = (PetscInt*)malloc(sizeof(PetscInt)*n); 684 685 ivec_zero(dir ,level+1); 686 ivec_zero(nsep ,level+1); 687 ivec_zero(lnsep,level+1); 688 ivec_set (fo ,-1,n+1); 689 ivec_zero(used,n); 690 691 lhs = (double*)malloc(sizeof(PetscScalar)*m); 692 rhs = (double*)malloc(sizeof(PetscScalar)*m); 693 694 /* determine the # of unique dof */ 695 rvec_zero(lhs,m); 696 rvec_set(lhs,1.0,n); 697 gs_gop_hc(gs_handle,lhs,"+\0",level); 698 ierr = PetscInfo(0,"done first gs_gop_hc\n");CHKERRQ(ierr); 699 rvec_zero(rsum,2); 700 for (ct=i=0;i<n;i++) 701 { 702 if (lhs[i]!=0.0) 703 {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];} 704 705 if (lhs[i]!=1.0) 706 { 707 shared=TRUE; 708 } 709 } 710 711 grop_hc(rsum,rw,2,op,level); 712 rsum[0]+=0.1; 713 rsum[1]+=0.1; 714 715 /* 716 if (!my_id) 717 { 718 printf("xyt n unique = %d (%g)\n",(int) rsum[0], rsum[0]); 719 printf("xyt n shared = %d (%g)\n",(int) rsum[1], rsum[1]); 720 } 721 */ 722 723 xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0]; 724 xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0]; 725 726 /* determine separator sets top down */ 727 if (shared) 728 { 729 /* solution is to do as in the symmetric shared case but then */ 730 /* pick the sub-hc with the most free dofs and do a mat-vec */ 731 /* and pick up the responses on the other sub-hc from the */ 732 /* initial separator set obtained from the symm. shared case */ 733 SETERRQ(PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n"); 734 for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) 735 { 736 /* set rsh of hc, fire, and collect lhs responses */ 737 (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m); 738 gs_gop_hc(gs_handle,lhs,"+\0",edge); 739 740 /* set lsh of hc, fire, and collect rhs responses */ 741 (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m); 742 gs_gop_hc(gs_handle,rhs,"+\0",edge); 743 744 for (i=0;i<n;i++) 745 { 746 if (id< mask) 747 { 748 if (lhs[i]!=0.0) 749 {lhs[i]=1.0;} 750 } 751 if (id>=mask) 752 { 753 if (rhs[i]!=0.0) 754 {rhs[i]=1.0;} 755 } 756 } 757 758 if (id< mask) 759 {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);} 760 else 761 {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);} 762 763 /* count number of dofs I own that have signal and not in sep set */ 764 rvec_zero(rsum,4); 765 for (ivec_zero(sum,4),ct=i=0;i<n;i++) 766 { 767 if (!used[i]) 768 { 769 /* number of unmarked dofs on node */ 770 ct++; 771 /* number of dofs to be marked on lhs hc */ 772 if (id< mask) 773 { 774 if (lhs[i]!=0.0) 775 {sum[0]++; rsum[0]+=1.0/lhs[i];} 776 } 777 /* number of dofs to be marked on rhs hc */ 778 if (id>=mask) 779 { 780 if (rhs[i]!=0.0) 781 {sum[1]++; rsum[1]+=1.0/rhs[i];} 782 } 783 } 784 } 785 786 /* go for load balance - choose half with most unmarked dofs, bias LHS */ 787 (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 788 (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct); 789 giop_hc(sum,w,4,op,edge); 790 grop_hc(rsum,rw,4,op,edge); 791 rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1; 792 793 if (id<mask) 794 { 795 /* mark dofs I own that have signal and not in sep set */ 796 for (ct=i=0;i<n;i++) 797 { 798 if ((!used[i])&&(lhs[i]!=0.0)) 799 { 800 ct++; nfo++; 801 802 if (nfo>n) 803 {SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");} 804 805 *--iptr = local2global[i]; 806 used[i]=edge; 807 } 808 } 809 if (ct>1) {ivec_sort(iptr,ct);} 810 811 lnsep[edge]=ct; 812 nsep[edge]=(PetscInt) rsum[0]; 813 dir [edge]=LEFT; 814 } 815 816 if (id>=mask) 817 { 818 /* mark dofs I own that have signal and not in sep set */ 819 for (ct=i=0;i<n;i++) 820 { 821 if ((!used[i])&&(rhs[i]!=0.0)) 822 { 823 ct++; nfo++; 824 825 if (nfo>n) 826 {SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");} 827 828 *--iptr = local2global[i]; 829 used[i]=edge; 830 } 831 } 832 if (ct>1) {ivec_sort(iptr,ct);} 833 834 lnsep[edge]=ct; 835 nsep[edge]= (PetscInt) rsum[1]; 836 dir [edge]=RIGHT; 837 } 838 839 /* LATER or we can recur on these to order seps at this level */ 840 /* do we need full set of separators for this? */ 841 842 /* fold rhs hc into lower */ 843 if (id>=mask) 844 {id-=mask;} 845 } 846 } 847 else 848 { 849 for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) 850 { 851 /* set rsh of hc, fire, and collect lhs responses */ 852 (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m); 853 gs_gop_hc(gs_handle,lhs,"+\0",edge); 854 855 /* set lsh of hc, fire, and collect rhs responses */ 856 (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m); 857 gs_gop_hc(gs_handle,rhs,"+\0",edge); 858 859 /* count number of dofs I own that have signal and not in sep set */ 860 for (ivec_zero(sum,4),ct=i=0;i<n;i++) 861 { 862 if (!used[i]) 863 { 864 /* number of unmarked dofs on node */ 865 ct++; 866 /* number of dofs to be marked on lhs hc */ 867 if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;} 868 /* number of dofs to be marked on rhs hc */ 869 if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;} 870 } 871 } 872 873 /* for the non-symmetric case we need separators of width 2 */ 874 /* so take both sides */ 875 (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 876 giop_hc(sum,w,4,op,edge); 877 878 ct=0; 879 if (id<mask) 880 { 881 /* mark dofs I own that have signal and not in sep set */ 882 for (i=0;i<n;i++) 883 { 884 if ((!used[i])&&(lhs[i]!=0.0)) 885 { 886 ct++; nfo++; 887 *--iptr = local2global[i]; 888 used[i]=edge; 889 } 890 } 891 /* LSH hc summation of ct should be sum[0] */ 892 } 893 else 894 { 895 /* mark dofs I own that have signal and not in sep set */ 896 for (i=0;i<n;i++) 897 { 898 if ((!used[i])&&(rhs[i]!=0.0)) 899 { 900 ct++; nfo++; 901 *--iptr = local2global[i]; 902 used[i]=edge; 903 } 904 } 905 /* RSH hc summation of ct should be sum[1] */ 906 } 907 908 if (ct>1) {ivec_sort(iptr,ct);} 909 lnsep[edge]=ct; 910 nsep[edge]=sum[0]+sum[1]; 911 dir [edge]=BOTH; 912 913 /* LATER or we can recur on these to order seps at this level */ 914 /* do we need full set of separators for this? */ 915 916 /* fold rhs hc into lower */ 917 if (id>=mask) 918 {id-=mask;} 919 } 920 } 921 922 /* level 0 is on processor case - so mark the remainder */ 923 for (ct=i=0;i<n;i++) 924 { 925 if (!used[i]) 926 { 927 ct++; nfo++; 928 *--iptr = local2global[i]; 929 used[i]=edge; 930 } 931 } 932 if (ct>1) {ivec_sort(iptr,ct);} 933 lnsep[edge]=ct; 934 nsep [edge]=ct; 935 dir [edge]=BOTH; 936 937 xyt_handle->info->nsep=nsep; 938 xyt_handle->info->lnsep=lnsep; 939 xyt_handle->info->fo=fo; 940 xyt_handle->info->nfo=nfo; 941 942 free(dir); 943 free(lhs); 944 free(rhs); 945 free(used); 946 PetscFunctionReturn(0); 947 } 948 949 /**************************************xyt.c***********************************/ 950 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data) 951 { 952 mv_info *mvi; 953 954 955 mvi = (mv_info*)malloc(sizeof(mv_info)); 956 mvi->n=n; 957 mvi->m=m; 958 mvi->n_global=-1; 959 mvi->m_global=-1; 960 mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt)); 961 ivec_copy(mvi->local2global,local2global,m); 962 mvi->local2global[m] = INT_MAX; 963 mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec; 964 mvi->grid_data=grid_data; 965 966 /* set xyt communication handle to perform restricted matvec */ 967 mvi->gs_handle = gs_init(local2global, m, num_nodes); 968 969 return(mvi); 970 } 971 972 /**************************************xyt.c***********************************/ 973 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 974 { 975 PetscFunctionBegin; 976 A->matvec((mv_info*)A->grid_data,v,u); 977 PetscFunctionReturn(0); 978 } 979 980 981 982