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