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