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_COMM_SELF,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, PetscScalar *x, PetscScalar *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 PetscErrorCode ierr; 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 {ierr = PetscPrintf(PETSC_COMM_WORLD,"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 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_nnz=%D\n",my_id,vals[0]); 195 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_nnz=%D\n",my_id,vals[1]); 196 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes); 197 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_nnz=%D\n",my_id,vals[2]); 198 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(2d) =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5))); 199 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(3d) =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667))); 200 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_n =%D\n",my_id,vals[3]); 201 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_n =%D\n",my_id,vals[4]); 202 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_n =%g\n",my_id,1.0*vals[5]/num_nodes); 203 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_n =%D\n",my_id,vals[5]); 204 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_buf=%D\n",my_id,vals[6]); 205 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_buf=%D\n",my_id,vals[7]); 206 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_buf=%g\n",my_id,1.0*vals[8]/num_nodes); 207 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_slv=%g\n",my_id,fvals[0]); 208 ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_slv=%g\n",my_id,fvals[1]); 209 ierr = PetscPrintf(PETSC_COMM_WORLD,"%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,dlen; 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 {ierr = PetscPrintf(PETSC_COMM_WORLD,"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 if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n"); 346 stages[dim++]=i; 347 end+=lnsep[dim]; 348 } 349 stages[dim]=i; 350 351 /* which column are we firing? */ 352 /* i.e. set v_l */ 353 /* use new seps and do global min across hc to determine which one to fire */ 354 (start<end) ? (col=fo[start]) : (col=INT_MAX); 355 giop_hc(&col,&work,1,op2,dim); 356 357 /* shouldn't need this */ 358 if (col==INT_MAX) 359 { 360 ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr); 361 continue; 362 } 363 364 /* do I own it? I should */ 365 rvec_zero(v ,a_m); 366 if (col==fo[start]) 367 { 368 start++; 369 idx=ivec_linear_search(col, a_local2global, a_n); 370 if (idx!=-1) 371 {v[idx] = 1.0; j++;} 372 else 373 {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");} 374 } 375 else 376 { 377 idx=ivec_linear_search(col, a_local2global, a_m); 378 if (idx!=-1) 379 {v[idx] = 1.0;} 380 } 381 382 /* perform u = A.v_l */ 383 rvec_zero(u,n); 384 do_matvec(xyt_handle->mvi,v,u); 385 386 /* uu = X^T.u_l (local portion) */ 387 /* technically only need to zero out first i entries */ 388 /* later turn this into an XYT_solve call ? */ 389 rvec_zero(uu,m); 390 y_ptr=y; 391 iptr = ycol_indices; 392 for (k=0; k<i; k++) 393 { 394 off = *iptr++; 395 len = *iptr++; 396 dlen = PetscBLASIntCast(len); 397 uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1); 398 y_ptr+=len; 399 } 400 401 /* uu = X^T.u_l (comm portion) */ 402 ssgl_radd (uu, w, dim, stages); 403 404 /* z = X.uu */ 405 rvec_zero(z,n); 406 x_ptr=x; 407 iptr = xcol_indices; 408 for (k=0; k<i; k++) 409 { 410 off = *iptr++; 411 len = *iptr++; 412 dlen = PetscBLASIntCast(len); 413 BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1); 414 x_ptr+=len; 415 } 416 417 /* compute v_l = v_l - z */ 418 rvec_zero(v+a_n,a_m-a_n); 419 dlen = PetscBLASIntCast(n); 420 BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1); 421 422 /* compute u_l = A.v_l */ 423 if (a_n!=a_m) 424 {gs_gop_hc(gs_handle,v,"+\0",dim);} 425 rvec_zero(u,n); 426 do_matvec(xyt_handle->mvi,v,u); 427 428 /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 429 dlen = PetscBLASIntCast(n); 430 alpha = BLASdot_(&dlen,u,&i1,u,&i1); 431 /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 432 grop_hc(&alpha, &alpha_w, 1, op, dim); 433 434 alpha = (PetscScalar) sqrt((double)alpha); 435 436 /* check for small alpha */ 437 /* LATER use this to detect and determine null space */ 438 if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha); 439 440 /* compute v_l = v_l/sqrt(alpha) */ 441 rvec_scale(v,1.0/alpha,n); 442 rvec_scale(u,1.0/alpha,n); 443 444 /* add newly generated column, v_l, to X */ 445 flag = 1; 446 off=len=0; 447 for (k=0; k<n; k++) 448 { 449 if (v[k]!=0.0) 450 { 451 len=k; 452 if (flag) 453 {off=k; flag=0;} 454 } 455 } 456 457 len -= (off-1); 458 459 if (len>0) 460 { 461 if ((xt_nnz+len)>xt_max_nnz) 462 { 463 ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr); 464 xt_max_nnz *= 2; 465 x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar)); 466 rvec_copy(x_ptr,x,xt_nnz); 467 free(x); 468 x = x_ptr; 469 x_ptr+=xt_nnz; 470 } 471 xt_nnz += len; 472 rvec_copy(x_ptr,v+off,len); 473 474 /* keep track of number of zeros */ 475 if (dim) 476 { 477 for (k=0; k<len; k++) 478 { 479 if (x_ptr[k]==0.0) 480 {xt_zero_nnz++;} 481 } 482 } 483 else 484 { 485 for (k=0; k<len; k++) 486 { 487 if (x_ptr[k]==0.0) 488 {xt_zero_nnz_0++;} 489 } 490 } 491 xcol_indices[2*i] = off; 492 xcol_sz[i] = xcol_indices[2*i+1] = len; 493 xcol_vals[i] = x_ptr; 494 } 495 else 496 { 497 xcol_indices[2*i] = 0; 498 xcol_sz[i] = xcol_indices[2*i+1] = 0; 499 xcol_vals[i] = x_ptr; 500 } 501 502 503 /* add newly generated column, u_l, to Y */ 504 flag = 1; 505 off=len=0; 506 for (k=0; k<n; k++) 507 { 508 if (u[k]!=0.0) 509 { 510 len=k; 511 if (flag) 512 {off=k; flag=0;} 513 } 514 } 515 516 len -= (off-1); 517 518 if (len>0) 519 { 520 if ((yt_nnz+len)>yt_max_nnz) 521 { 522 ierr = PetscInfo(0,"increasing space for Y by 2x!\n");CHKERRQ(ierr); 523 yt_max_nnz *= 2; 524 y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar)); 525 rvec_copy(y_ptr,y,yt_nnz); 526 free(y); 527 y = y_ptr; 528 y_ptr+=yt_nnz; 529 } 530 yt_nnz += len; 531 rvec_copy(y_ptr,u+off,len); 532 533 /* keep track of number of zeros */ 534 if (dim) 535 { 536 for (k=0; k<len; k++) 537 { 538 if (y_ptr[k]==0.0) 539 {yt_zero_nnz++;} 540 } 541 } 542 else 543 { 544 for (k=0; k<len; k++) 545 { 546 if (y_ptr[k]==0.0) 547 {yt_zero_nnz_0++;} 548 } 549 } 550 ycol_indices[2*i] = off; 551 ycol_sz[i] = ycol_indices[2*i+1] = len; 552 ycol_vals[i] = y_ptr; 553 } 554 else 555 { 556 ycol_indices[2*i] = 0; 557 ycol_sz[i] = ycol_indices[2*i+1] = 0; 558 ycol_vals[i] = y_ptr; 559 } 560 } 561 562 /* close off stages for execution phase */ 563 while (dim!=level) 564 { 565 stages[dim++]=i; 566 ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr); 567 } 568 stages[dim]=i; 569 570 xyt_handle->info->n=xyt_handle->mvi->n; 571 xyt_handle->info->m=m; 572 xyt_handle->info->nnz=xt_nnz + yt_nnz; 573 xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz; 574 xyt_handle->info->msg_buf_sz=stages[level]-stages[0]; 575 xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 576 xyt_handle->info->solve_w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 577 xyt_handle->info->x=x; 578 xyt_handle->info->xcol_vals=xcol_vals; 579 xyt_handle->info->xcol_sz=xcol_sz; 580 xyt_handle->info->xcol_indices=xcol_indices; 581 xyt_handle->info->stages=stages; 582 xyt_handle->info->y=y; 583 xyt_handle->info->ycol_vals=ycol_vals; 584 xyt_handle->info->ycol_sz=ycol_sz; 585 xyt_handle->info->ycol_indices=ycol_indices; 586 587 free(segs); 588 free(u); 589 free(v); 590 free(uu); 591 free(z); 592 free(w); 593 594 return(0); 595 } 596 597 /**************************************xyt.c***********************************/ 598 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 599 { 600 PetscInt off, len, *iptr; 601 PetscInt level =xyt_handle->level; 602 PetscInt n =xyt_handle->info->n; 603 PetscInt m =xyt_handle->info->m; 604 PetscInt *stages =xyt_handle->info->stages; 605 PetscInt *xcol_indices=xyt_handle->info->xcol_indices; 606 PetscInt *ycol_indices=xyt_handle->info->ycol_indices; 607 PetscScalar *x_ptr, *y_ptr, *uu_ptr; 608 PetscScalar *solve_uu=xyt_handle->info->solve_uu; 609 PetscScalar *solve_w =xyt_handle->info->solve_w; 610 PetscScalar *x =xyt_handle->info->x; 611 PetscScalar *y =xyt_handle->info->y; 612 PetscBLASInt i1 = 1,dlen; 613 614 PetscFunctionBegin; 615 uu_ptr=solve_uu; 616 rvec_zero(uu_ptr,m); 617 618 /* x = X.Y^T.b */ 619 /* uu = Y^T.b */ 620 for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) 621 { 622 off=*iptr++; len=*iptr++; 623 dlen = PetscBLASIntCast(len); 624 *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1); 625 } 626 627 /* comunication of beta */ 628 uu_ptr=solve_uu; 629 if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);} 630 631 rvec_zero(uc,n); 632 633 /* x = X.uu */ 634 for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) 635 { 636 off=*iptr++; len=*iptr++; 637 dlen = PetscBLASIntCast(len); 638 BLASaxpy_(&dlen,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) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle); 650 651 vals[0]=vals[1]=xyt_handle->id; 652 giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 653 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); 654 PetscFunctionReturn(0); 655 } 656 657 /**************************************xyt.c***********************************/ 658 static PetscErrorCode det_separators(xyt_ADT xyt_handle) 659 { 660 PetscInt i, ct, id; 661 PetscInt mask, edge, *iptr; 662 PetscInt *dir, *used; 663 PetscInt sum[4], w[4]; 664 PetscScalar rsum[4], rw[4]; 665 PetscInt op[] = {GL_ADD,0}; 666 PetscScalar *lhs, *rhs; 667 PetscInt *nsep, *lnsep, *fo, nfo=0; 668 gs_ADT gs_handle=xyt_handle->mvi->gs_handle; 669 PetscInt *local2global=xyt_handle->mvi->local2global; 670 PetscInt n=xyt_handle->mvi->n; 671 PetscInt m=xyt_handle->mvi->m; 672 PetscInt level=xyt_handle->level; 673 PetscInt shared=FALSE; 674 PetscErrorCode ierr; 675 676 PetscFunctionBegin; 677 dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 678 nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 679 lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 680 fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 681 used = (PetscInt*)malloc(sizeof(PetscInt)*n); 682 683 ivec_zero(dir ,level+1); 684 ivec_zero(nsep ,level+1); 685 ivec_zero(lnsep,level+1); 686 ivec_set (fo ,-1,n+1); 687 ivec_zero(used,n); 688 689 lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 690 rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 691 692 /* determine the # of unique dof */ 693 rvec_zero(lhs,m); 694 rvec_set(lhs,1.0,n); 695 gs_gop_hc(gs_handle,lhs,"+\0",level); 696 ierr = PetscInfo(0,"done first gs_gop_hc\n");CHKERRQ(ierr); 697 rvec_zero(rsum,2); 698 for (ct=i=0;i<n;i++) 699 { 700 if (lhs[i]!=0.0) 701 {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];} 702 703 if (lhs[i]!=1.0) 704 { 705 shared=TRUE; 706 } 707 } 708 709 grop_hc(rsum,rw,2,op,level); 710 rsum[0]+=0.1; 711 rsum[1]+=0.1; 712 713 xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0]; 714 xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0]; 715 716 /* determine separator sets top down */ 717 if (shared) 718 { 719 /* solution is to do as in the symmetric shared case but then */ 720 /* pick the sub-hc with the most free dofs and do a mat-vec */ 721 /* and pick up the responses on the other sub-hc from the */ 722 /* initial separator set obtained from the symm. shared case */ 723 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n"); 724 for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) 725 { 726 /* set rsh of hc, fire, and collect lhs responses */ 727 (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m); 728 gs_gop_hc(gs_handle,lhs,"+\0",edge); 729 730 /* set lsh of hc, fire, and collect rhs responses */ 731 (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m); 732 gs_gop_hc(gs_handle,rhs,"+\0",edge); 733 734 for (i=0;i<n;i++) 735 { 736 if (id< mask) 737 { 738 if (lhs[i]!=0.0) 739 {lhs[i]=1.0;} 740 } 741 if (id>=mask) 742 { 743 if (rhs[i]!=0.0) 744 {rhs[i]=1.0;} 745 } 746 } 747 748 if (id< mask) 749 {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);} 750 else 751 {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);} 752 753 /* count number of dofs I own that have signal and not in sep set */ 754 rvec_zero(rsum,4); 755 for (ivec_zero(sum,4),ct=i=0;i<n;i++) 756 { 757 if (!used[i]) 758 { 759 /* number of unmarked dofs on node */ 760 ct++; 761 /* number of dofs to be marked on lhs hc */ 762 if (id< mask) 763 { 764 if (lhs[i]!=0.0) 765 {sum[0]++; rsum[0]+=1.0/lhs[i];} 766 } 767 /* number of dofs to be marked on rhs hc */ 768 if (id>=mask) 769 { 770 if (rhs[i]!=0.0) 771 {sum[1]++; rsum[1]+=1.0/rhs[i];} 772 } 773 } 774 } 775 776 /* go for load balance - choose half with most unmarked dofs, bias LHS */ 777 (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 778 (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct); 779 giop_hc(sum,w,4,op,edge); 780 grop_hc(rsum,rw,4,op,edge); 781 rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1; 782 783 if (id<mask) 784 { 785 /* mark dofs I own that have signal and not in sep set */ 786 for (ct=i=0;i<n;i++) 787 { 788 if ((!used[i])&&(lhs[i]!=0.0)) 789 { 790 ct++; nfo++; 791 792 if (nfo>n) 793 {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");} 794 795 *--iptr = local2global[i]; 796 used[i]=edge; 797 } 798 } 799 if (ct>1) {ivec_sort(iptr,ct);} 800 801 lnsep[edge]=ct; 802 nsep[edge]=(PetscInt) rsum[0]; 803 dir [edge]=LEFT; 804 } 805 806 if (id>=mask) 807 { 808 /* mark dofs I own that have signal and not in sep set */ 809 for (ct=i=0;i<n;i++) 810 { 811 if ((!used[i])&&(rhs[i]!=0.0)) 812 { 813 ct++; nfo++; 814 815 if (nfo>n) 816 {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");} 817 818 *--iptr = local2global[i]; 819 used[i]=edge; 820 } 821 } 822 if (ct>1) {ivec_sort(iptr,ct);} 823 824 lnsep[edge]=ct; 825 nsep[edge]= (PetscInt) rsum[1]; 826 dir [edge]=RIGHT; 827 } 828 829 /* LATER or we can recur on these to order seps at this level */ 830 /* do we need full set of separators for this? */ 831 832 /* fold rhs hc into lower */ 833 if (id>=mask) 834 {id-=mask;} 835 } 836 } 837 else 838 { 839 for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) 840 { 841 /* set rsh of hc, fire, and collect lhs responses */ 842 (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m); 843 gs_gop_hc(gs_handle,lhs,"+\0",edge); 844 845 /* set lsh of hc, fire, and collect rhs responses */ 846 (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m); 847 gs_gop_hc(gs_handle,rhs,"+\0",edge); 848 849 /* count number of dofs I own that have signal and not in sep set */ 850 for (ivec_zero(sum,4),ct=i=0;i<n;i++) 851 { 852 if (!used[i]) 853 { 854 /* number of unmarked dofs on node */ 855 ct++; 856 /* number of dofs to be marked on lhs hc */ 857 if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;} 858 /* number of dofs to be marked on rhs hc */ 859 if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;} 860 } 861 } 862 863 /* for the non-symmetric case we need separators of width 2 */ 864 /* so take both sides */ 865 (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 866 giop_hc(sum,w,4,op,edge); 867 868 ct=0; 869 if (id<mask) 870 { 871 /* mark dofs I own that have signal and not in sep set */ 872 for (i=0;i<n;i++) 873 { 874 if ((!used[i])&&(lhs[i]!=0.0)) 875 { 876 ct++; nfo++; 877 *--iptr = local2global[i]; 878 used[i]=edge; 879 } 880 } 881 /* LSH hc summation of ct should be sum[0] */ 882 } 883 else 884 { 885 /* mark dofs I own that have signal and not in sep set */ 886 for (i=0;i<n;i++) 887 { 888 if ((!used[i])&&(rhs[i]!=0.0)) 889 { 890 ct++; nfo++; 891 *--iptr = local2global[i]; 892 used[i]=edge; 893 } 894 } 895 /* RSH hc summation of ct should be sum[1] */ 896 } 897 898 if (ct>1) {ivec_sort(iptr,ct);} 899 lnsep[edge]=ct; 900 nsep[edge]=sum[0]+sum[1]; 901 dir [edge]=BOTH; 902 903 /* LATER or we can recur on these to order seps at this level */ 904 /* do we need full set of separators for this? */ 905 906 /* fold rhs hc into lower */ 907 if (id>=mask) 908 {id-=mask;} 909 } 910 } 911 912 /* level 0 is on processor case - so mark the remainder */ 913 for (ct=i=0;i<n;i++) 914 { 915 if (!used[i]) 916 { 917 ct++; nfo++; 918 *--iptr = local2global[i]; 919 used[i]=edge; 920 } 921 } 922 if (ct>1) {ivec_sort(iptr,ct);} 923 lnsep[edge]=ct; 924 nsep [edge]=ct; 925 dir [edge]=BOTH; 926 927 xyt_handle->info->nsep=nsep; 928 xyt_handle->info->lnsep=lnsep; 929 xyt_handle->info->fo=fo; 930 xyt_handle->info->nfo=nfo; 931 932 free(dir); 933 free(lhs); 934 free(rhs); 935 free(used); 936 PetscFunctionReturn(0); 937 } 938 939 /**************************************xyt.c***********************************/ 940 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data) 941 { 942 mv_info *mvi; 943 944 945 mvi = (mv_info*)malloc(sizeof(mv_info)); 946 mvi->n=n; 947 mvi->m=m; 948 mvi->n_global=-1; 949 mvi->m_global=-1; 950 mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt)); 951 ivec_copy(mvi->local2global,local2global,m); 952 mvi->local2global[m] = INT_MAX; 953 mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec; 954 mvi->grid_data=grid_data; 955 956 /* set xyt communication handle to perform restricted matvec */ 957 mvi->gs_handle = gs_init(local2global, m, num_nodes); 958 959 return(mvi); 960 } 961 962 /**************************************xyt.c***********************************/ 963 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 964 { 965 PetscFunctionBegin; 966 A->matvec((mv_info*)A->grid_data,v,u); 967 PetscFunctionReturn(0); 968 } 969 970 971 972