1 2 /*************************************xyt.c************************************ 3 Module Name: xyt 4 Module Info: 5 6 author: Henry M. Tufo III 7 e-mail: hmt@asci.uchicago.edu 8 contact: 9 +--------------------------------+--------------------------------+ 10 |MCS Division - Building 221 |Department of Computer Science | 11 |Argonne National Laboratory |Ryerson 152 | 12 |9700 S. Cass Avenue |The University of Chicago | 13 |Argonne, IL 60439 |Chicago, IL 60637 | 14 |(630) 252-5354/5986 ph/fx |(773) 702-6019/8487 ph/fx | 15 +--------------------------------+--------------------------------+ 16 17 Last Modification: 3.20.01 18 **************************************xyt.c***********************************/ 19 #include <../src/ksp/pc/impls/tfs/tfs.h> 20 21 #define LEFT -1 22 #define RIGHT 1 23 #define BOTH 0 24 25 typedef struct xyt_solver_info { 26 PetscInt n, m, n_global, m_global; 27 PetscInt nnz, max_nnz, msg_buf_sz; 28 PetscInt *nsep, *lnsep, *fo, nfo, *stages; 29 PetscInt *xcol_sz, *xcol_indices; 30 PetscScalar **xcol_vals, *x, *solve_uu, *solve_w; 31 PetscInt *ycol_sz, *ycol_indices; 32 PetscScalar **ycol_vals, *y; 33 PetscInt nsolves; 34 PetscScalar tot_solve_time; 35 } xyt_info; 36 37 typedef struct matvec_info { 38 PetscInt n, m, n_global, m_global; 39 PetscInt *local2global; 40 PCTFS_gs_ADT PCTFS_gs_handle; 41 PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*); 42 void *grid_data; 43 } mv_info; 44 45 struct xyt_CDT { 46 PetscInt id; 47 PetscInt ns; 48 PetscInt level; 49 xyt_info *info; 50 mv_info *mvi; 51 }; 52 53 static PetscInt n_xyt =0; 54 static PetscInt n_xyt_handles=0; 55 56 /* prototypes */ 57 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs); 58 static PetscErrorCode check_handle(xyt_ADT xyt_handle); 59 static PetscErrorCode det_separators(xyt_ADT xyt_handle); 60 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 61 static PetscErrorCode xyt_generate(xyt_ADT xyt_handle); 62 static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle); 63 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data); 64 65 /**************************************xyt.c***********************************/ 66 xyt_ADT XYT_new(void) 67 { 68 xyt_ADT xyt_handle; 69 70 /* rolling count on n_xyt ... pot. problem here */ 71 n_xyt_handles++; 72 xyt_handle = (xyt_ADT)malloc(sizeof(struct xyt_CDT)); 73 xyt_handle->id = ++n_xyt; 74 xyt_handle->info = NULL; 75 xyt_handle->mvi = NULL; 76 77 return(xyt_handle); 78 } 79 80 /**************************************xyt.c***********************************/ 81 PetscErrorCode XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt handle */ 82 PetscInt *local2global, /* global column mapping */ 83 PetscInt n, /* local num rows */ 84 PetscInt m, /* local num cols */ 85 PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */ 86 void *grid_data) /* grid data for matvec */ 87 { 88 89 PCTFS_comm_init(); 90 check_handle(xyt_handle); 91 92 /* only 2^k for now and all nodes participating */ 93 PetscCheck((1<<(xyt_handle->level=PCTFS_i_log2_num_nodes))==PCTFS_num_nodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %d != %d",1<<PCTFS_i_log2_num_nodes,PCTFS_num_nodes); 94 95 /* space for X info */ 96 xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info)); 97 98 /* set up matvec handles */ 99 xyt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data); 100 101 /* matrix is assumed to be of full rank */ 102 /* LATER we can reset to indicate rank def. */ 103 xyt_handle->ns=0; 104 105 /* determine separators and generate firing order - NB xyt info set here */ 106 det_separators(xyt_handle); 107 108 return(do_xyt_factor(xyt_handle)); 109 } 110 111 /**************************************xyt.c***********************************/ 112 PetscErrorCode XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b) 113 { 114 PCTFS_comm_init(); 115 check_handle(xyt_handle); 116 117 /* need to copy b into x? */ 118 if (b) PCTFS_rvec_copy(x,b,xyt_handle->mvi->n); 119 return do_xyt_solve(xyt_handle,x); 120 } 121 122 /**************************************xyt.c***********************************/ 123 PetscErrorCode XYT_free(xyt_ADT xyt_handle) 124 { 125 PCTFS_comm_init(); 126 check_handle(xyt_handle); 127 n_xyt_handles--; 128 129 free(xyt_handle->info->nsep); 130 free(xyt_handle->info->lnsep); 131 free(xyt_handle->info->fo); 132 free(xyt_handle->info->stages); 133 free(xyt_handle->info->solve_uu); 134 free(xyt_handle->info->solve_w); 135 free(xyt_handle->info->x); 136 free(xyt_handle->info->xcol_vals); 137 free(xyt_handle->info->xcol_sz); 138 free(xyt_handle->info->xcol_indices); 139 free(xyt_handle->info->y); 140 free(xyt_handle->info->ycol_vals); 141 free(xyt_handle->info->ycol_sz); 142 free(xyt_handle->info->ycol_indices); 143 free(xyt_handle->info); 144 free(xyt_handle->mvi->local2global); 145 PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle); 146 free(xyt_handle->mvi); 147 free(xyt_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 /**************************************xyt.c***********************************/ 156 PetscErrorCode XYT_stats(xyt_ADT xyt_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 PetscFunctionBegin; 164 PetscCall(PCTFS_comm_init()); 165 PetscCall(check_handle(xyt_handle)); 166 167 /* if factorization not done there are no stats */ 168 if (!xyt_handle->info||!xyt_handle->mvi) { 169 if (!PCTFS_my_id) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n")); 170 PetscFunctionReturn(1); 171 } 172 173 vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz; 174 vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n; 175 vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz; 176 PetscCall(PCTFS_giop(vals,work,PETSC_STATIC_ARRAY_LENGTH(op)-1,op)); 177 178 fvals[0]=fvals[1]=fvals[2]=xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++; 179 PetscCall(PCTFS_grop(fvals,fwork,PETSC_STATIC_ARRAY_LENGTH(fop)-1,fop)); 180 181 if (!PCTFS_my_id) { 182 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: min xyt_nnz=%" PetscInt_FMT "\n",PCTFS_my_id,vals[0])); 183 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: max xyt_nnz=%" PetscInt_FMT "\n",PCTFS_my_id,vals[1])); 184 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: avg xyt_nnz=%g\n",PCTFS_my_id,(double)(1.0*vals[2]/PCTFS_num_nodes))); 185 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: tot xyt_nnz=%" PetscInt_FMT "\n",PCTFS_my_id,vals[2])); 186 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: xyt C(2d) =%g\n",PCTFS_my_id,(double)(vals[2]/(PetscPowReal(1.0*vals[5],1.5))))); 187 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: xyt C(3d) =%g\n",PCTFS_my_id,(double)(vals[2]/(PetscPowReal(1.0*vals[5],1.6667))))); 188 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: min xyt_n =%" PetscInt_FMT "\n",PCTFS_my_id,vals[3])); 189 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: max xyt_n =%" PetscInt_FMT "\n",PCTFS_my_id,vals[4])); 190 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: avg xyt_n =%g\n",PCTFS_my_id,(double)(1.0*vals[5]/PCTFS_num_nodes))); 191 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: tot xyt_n =%" PetscInt_FMT "\n",PCTFS_my_id,vals[5])); 192 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: min xyt_buf=%" PetscInt_FMT "\n",PCTFS_my_id,vals[6])); 193 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: max xyt_buf=%" PetscInt_FMT "\n",PCTFS_my_id,vals[7])); 194 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: avg xyt_buf=%g\n",PCTFS_my_id,(double)(1.0*vals[8]/PCTFS_num_nodes))); 195 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: min xyt_slv=%g\n",PCTFS_my_id,(double)PetscRealPart(fvals[0]))); 196 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: max xyt_slv=%g\n",PCTFS_my_id,(double)PetscRealPart(fvals[1]))); 197 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%d :: avg xyt_slv=%g\n",PCTFS_my_id,(double)PetscRealPart(fvals[2]/PCTFS_num_nodes))); 198 } 199 PetscFunctionReturn(0); 200 } 201 202 /*************************************xyt.c************************************ 203 204 Description: get A_local, local portion of global coarse matrix which 205 is a row dist. nxm matrix w/ n<m. 206 o my_ml holds address of ML struct associated w/A_local and coarse grid 207 o local2global holds global number of column i (i=0,...,m-1) 208 o local2global holds global number of row i (i=0,...,n-1) 209 o mylocmatvec performs A_local . vec_local (note that gs is performed using 210 PCTFS_gs_init/gop). 211 212 mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 213 mylocmatvec (void :: void *data, double *in, double *out) 214 **************************************xyt.c***********************************/ 215 static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle) 216 { 217 return xyt_generate(xyt_handle); 218 } 219 220 /**************************************xyt.c***********************************/ 221 static PetscErrorCode xyt_generate(xyt_ADT xyt_handle) 222 { 223 PetscInt i,j,k,idx; 224 PetscInt dim, col; 225 PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 226 PetscInt *segs; 227 PetscInt op[] = {GL_ADD,0}; 228 PetscInt off, len; 229 PetscScalar *x_ptr, *y_ptr; 230 PetscInt *iptr, flag; 231 PetscInt start =0, end, work; 232 PetscInt op2[] = {GL_MIN,0}; 233 PCTFS_gs_ADT PCTFS_gs_handle; 234 PetscInt *nsep, *lnsep, *fo; 235 PetscInt a_n =xyt_handle->mvi->n; 236 PetscInt a_m =xyt_handle->mvi->m; 237 PetscInt *a_local2global=xyt_handle->mvi->local2global; 238 PetscInt level; 239 PetscInt n, m; 240 PetscInt *xcol_sz, *xcol_indices, *stages; 241 PetscScalar **xcol_vals, *x; 242 PetscInt *ycol_sz, *ycol_indices; 243 PetscScalar **ycol_vals, *y; 244 PetscInt n_global; 245 PetscInt xt_nnz =0, xt_max_nnz=0; 246 PetscInt yt_nnz =0, yt_max_nnz=0; 247 PetscInt xt_zero_nnz =0; 248 PetscInt xt_zero_nnz_0=0; 249 PetscInt yt_zero_nnz =0; 250 PetscInt yt_zero_nnz_0=0; 251 PetscBLASInt i1 = 1,dlen; 252 PetscScalar dm1 = -1.0; 253 254 n =xyt_handle->mvi->n; 255 nsep =xyt_handle->info->nsep; 256 lnsep =xyt_handle->info->lnsep; 257 fo =xyt_handle->info->fo; 258 end =lnsep[0]; 259 level =xyt_handle->level; 260 PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 261 262 /* is there a null space? */ 263 /* LATER add in ability to detect null space by checking alpha */ 264 for (i=0, j=0; i<=level; i++) j+=nsep[i]; 265 266 m = j-xyt_handle->ns; 267 if (m!=j) { 268 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",m,j,xyt_handle->ns)); 269 } 270 271 PetscCall(PetscInfo(0,"xyt_generate() :: X(%" PetscInt_FMT ",%" PetscInt_FMT ")\n",n,m)); 272 273 /* get and initialize storage for x local */ 274 /* note that x local is nxm and stored by columns */ 275 xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 276 xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 277 xcol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 278 for (i=j=0; i<m; i++, j+=2) { 279 xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1; 280 xcol_vals[i] = NULL; 281 } 282 xcol_indices[j]=-1; 283 284 /* get and initialize storage for y local */ 285 /* note that y local is nxm and stored by columns */ 286 ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 287 ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 288 ycol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 289 for (i=j=0; i<m; i++, j+=2) { 290 ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1; 291 ycol_vals[i] = NULL; 292 } 293 ycol_indices[j]=-1; 294 295 /* size of separators for each sub-hc working from bottom of tree to top */ 296 /* this looks like nsep[]=segments */ 297 stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 298 segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 299 PCTFS_ivec_zero(stages,level+1); 300 PCTFS_ivec_copy(segs,nsep,level+1); 301 for (i=0; i<level; i++) segs[i+1] += segs[i]; 302 stages[0] = segs[0]; 303 304 /* temporary vectors */ 305 u = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 306 z = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 307 v = (PetscScalar*) malloc(a_m*sizeof(PetscScalar)); 308 uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 309 w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 310 311 /* extra nnz due to replication of vertices across separators */ 312 for (i=1, j=0; i<=level; i++) j+=nsep[i]; 313 314 /* storage for sparse x values */ 315 n_global = xyt_handle->info->n_global; 316 xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 317 x = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 318 y = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 319 320 /* LATER - can embed next sep to fire in gs */ 321 /* time to make the donuts - generate X factor */ 322 for (dim=i=j=0; i<m; i++) { 323 /* time to move to the next level? */ 324 while (i==segs[dim]) { 325 PetscCheck(dim!=level,PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level"); 326 stages[dim++]=i; 327 end +=lnsep[dim]; 328 } 329 stages[dim]=i; 330 331 /* which column are we firing? */ 332 /* i.e. set v_l */ 333 /* use new seps and do global min across hc to determine which one to fire */ 334 (start<end) ? (col=fo[start]) : (col=INT_MAX); 335 PCTFS_giop_hc(&col,&work,1,op2,dim); 336 337 /* shouldn't need this */ 338 if (col==INT_MAX) { 339 PetscCall(PetscInfo(0,"hey ... col==INT_MAX??\n")); 340 continue; 341 } 342 343 /* do I own it? I should */ 344 PCTFS_rvec_zero(v,a_m); 345 if (col==fo[start]) { 346 start++; 347 idx=PCTFS_ivec_linear_search(col, a_local2global, a_n); 348 if (idx!=-1) { 349 v[idx] = 1.0; 350 j++; 351 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!"); 352 } else { 353 idx=PCTFS_ivec_linear_search(col, a_local2global, a_m); 354 if (idx!=-1) v[idx] = 1.0; 355 } 356 357 /* perform u = A.v_l */ 358 PCTFS_rvec_zero(u,n); 359 do_matvec(xyt_handle->mvi,v,u); 360 361 /* uu = X^T.u_l (local portion) */ 362 /* technically only need to zero out first i entries */ 363 /* later turn this into an XYT_solve call ? */ 364 PCTFS_rvec_zero(uu,m); 365 y_ptr=y; 366 iptr = ycol_indices; 367 for (k=0; k<i; k++) { 368 off = *iptr++; 369 len = *iptr++; 370 PetscCall(PetscBLASIntCast(len,&dlen)); 371 PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1)); 372 y_ptr+=len; 373 } 374 375 /* uu = X^T.u_l (comm portion) */ 376 PetscCall(PCTFS_ssgl_radd (uu, w, dim, stages)); 377 378 /* z = X.uu */ 379 PCTFS_rvec_zero(z,n); 380 x_ptr=x; 381 iptr = xcol_indices; 382 for (k=0; k<i; k++) { 383 off = *iptr++; 384 len = *iptr++; 385 PetscCall(PetscBLASIntCast(len,&dlen)); 386 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1)); 387 x_ptr+=len; 388 } 389 390 /* compute v_l = v_l - z */ 391 PCTFS_rvec_zero(v+a_n,a_m-a_n); 392 PetscCall(PetscBLASIntCast(n,&dlen)); 393 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1)); 394 395 /* compute u_l = A.v_l */ 396 if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); 397 PCTFS_rvec_zero(u,n); 398 do_matvec(xyt_handle->mvi,v,u); 399 400 /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 401 PetscCall(PetscBLASIntCast(n,&dlen)); 402 PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,u,&i1)); 403 /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 404 PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 405 406 alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 407 408 /* check for small alpha */ 409 /* LATER use this to detect and determine null space */ 410 PetscCheck(PetscAbsScalar(alpha)>=1.0e-14,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g",(double)PetscAbsScalar(alpha)); 411 412 /* compute v_l = v_l/sqrt(alpha) */ 413 PCTFS_rvec_scale(v,1.0/alpha,n); 414 PCTFS_rvec_scale(u,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 if (v[k]!=0.0) { 421 len=k; 422 if (flag) {off=k; flag=0;} 423 } 424 } 425 426 len -= (off-1); 427 428 if (len>0) { 429 if ((xt_nnz+len)>xt_max_nnz) { 430 PetscCall(PetscInfo(0,"increasing space for X by 2x!\n")); 431 xt_max_nnz *= 2; 432 x_ptr = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 433 PCTFS_rvec_copy(x_ptr,x,xt_nnz); 434 free(x); 435 x = x_ptr; 436 x_ptr+=xt_nnz; 437 } 438 xt_nnz += len; 439 PCTFS_rvec_copy(x_ptr,v+off,len); 440 441 /* keep track of number of zeros */ 442 if (dim) { 443 for (k=0; k<len; k++) { 444 if (x_ptr[k]==0.0) xt_zero_nnz++; 445 } 446 } else { 447 for (k=0; k<len; k++) { 448 if (x_ptr[k]==0.0) xt_zero_nnz_0++; 449 } 450 } 451 xcol_indices[2*i] = off; 452 xcol_sz[i] = xcol_indices[2*i+1] = len; 453 xcol_vals[i] = x_ptr; 454 } else { 455 xcol_indices[2*i] = 0; 456 xcol_sz[i] = xcol_indices[2*i+1] = 0; 457 xcol_vals[i] = x_ptr; 458 } 459 460 /* add newly generated column, u_l, to Y */ 461 flag = 1; 462 off =len=0; 463 for (k=0; k<n; k++) { 464 if (u[k]!=0.0) { 465 len=k; 466 if (flag) { off=k; flag=0; } 467 } 468 } 469 470 len -= (off-1); 471 472 if (len>0) { 473 if ((yt_nnz+len)>yt_max_nnz) { 474 PetscCall(PetscInfo(0,"increasing space for Y by 2x!\n")); 475 yt_max_nnz *= 2; 476 y_ptr = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 477 PCTFS_rvec_copy(y_ptr,y,yt_nnz); 478 free(y); 479 y = y_ptr; 480 y_ptr+=yt_nnz; 481 } 482 yt_nnz += len; 483 PCTFS_rvec_copy(y_ptr,u+off,len); 484 485 /* keep track of number of zeros */ 486 if (dim) { 487 for (k=0; k<len; k++) { 488 if (y_ptr[k]==0.0) yt_zero_nnz++; 489 } 490 } else { 491 for (k=0; k<len; k++) { 492 if (y_ptr[k]==0.0) yt_zero_nnz_0++; 493 } 494 } 495 ycol_indices[2*i] = off; 496 ycol_sz[i] = ycol_indices[2*i+1] = len; 497 ycol_vals[i] = y_ptr; 498 } else { 499 ycol_indices[2*i] = 0; 500 ycol_sz[i] = ycol_indices[2*i+1] = 0; 501 ycol_vals[i] = y_ptr; 502 } 503 } 504 505 /* close off stages for execution phase */ 506 while (dim!=level) { 507 stages[dim++]=i; 508 PetscCall(PetscInfo(0,"disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n",dim,level)); 509 } 510 stages[dim]=i; 511 512 xyt_handle->info->n =xyt_handle->mvi->n; 513 xyt_handle->info->m =m; 514 xyt_handle->info->nnz =xt_nnz + yt_nnz; 515 xyt_handle->info->max_nnz =xt_max_nnz + yt_max_nnz; 516 xyt_handle->info->msg_buf_sz =stages[level]-stages[0]; 517 xyt_handle->info->solve_uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 518 xyt_handle->info->solve_w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 519 xyt_handle->info->x =x; 520 xyt_handle->info->xcol_vals =xcol_vals; 521 xyt_handle->info->xcol_sz =xcol_sz; 522 xyt_handle->info->xcol_indices=xcol_indices; 523 xyt_handle->info->stages =stages; 524 xyt_handle->info->y =y; 525 xyt_handle->info->ycol_vals =ycol_vals; 526 xyt_handle->info->ycol_sz =ycol_sz; 527 xyt_handle->info->ycol_indices=ycol_indices; 528 529 free(segs); 530 free(u); 531 free(v); 532 free(uu); 533 free(z); 534 free(w); 535 536 return(0); 537 } 538 539 /**************************************xyt.c***********************************/ 540 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 541 { 542 PetscInt off, len, *iptr; 543 PetscInt level =xyt_handle->level; 544 PetscInt n =xyt_handle->info->n; 545 PetscInt m =xyt_handle->info->m; 546 PetscInt *stages =xyt_handle->info->stages; 547 PetscInt *xcol_indices=xyt_handle->info->xcol_indices; 548 PetscInt *ycol_indices=xyt_handle->info->ycol_indices; 549 PetscScalar *x_ptr, *y_ptr, *uu_ptr; 550 PetscScalar *solve_uu=xyt_handle->info->solve_uu; 551 PetscScalar *solve_w =xyt_handle->info->solve_w; 552 PetscScalar *x =xyt_handle->info->x; 553 PetscScalar *y =xyt_handle->info->y; 554 PetscBLASInt i1 = 1,dlen; 555 556 PetscFunctionBegin; 557 uu_ptr=solve_uu; 558 PCTFS_rvec_zero(uu_ptr,m); 559 560 /* x = X.Y^T.b */ 561 /* uu = Y^T.b */ 562 for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) { 563 off =*iptr++; 564 len =*iptr++; 565 PetscCall(PetscBLASIntCast(len,&dlen)); 566 PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1)); 567 } 568 569 /* comunication of beta */ 570 uu_ptr=solve_uu; 571 if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages)); 572 PCTFS_rvec_zero(uc,n); 573 574 /* x = X.uu */ 575 for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) { 576 off =*iptr++; 577 len =*iptr++; 578 PetscCall(PetscBLASIntCast(len,&dlen)); 579 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1)); 580 } 581 PetscFunctionReturn(0); 582 } 583 584 /**************************************xyt.c***********************************/ 585 static PetscErrorCode check_handle(xyt_ADT xyt_handle) 586 { 587 PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 588 589 PetscFunctionBegin; 590 PetscCheck(xyt_handle,PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %p",xyt_handle); 591 592 vals[0]=vals[1]=xyt_handle->id; 593 PCTFS_giop(vals,work,PETSC_STATIC_ARRAY_LENGTH(op)-1,op); 594 PetscCheck(!(vals[0]!=vals[1])&&!(xyt_handle->id<=0),PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %" PetscInt_FMT "/%" PetscInt_FMT " %" PetscInt_FMT, vals[0],vals[1], xyt_handle->id); 595 PetscFunctionReturn(0); 596 } 597 598 /**************************************xyt.c***********************************/ 599 static PetscErrorCode det_separators(xyt_ADT xyt_handle) 600 { 601 PetscInt i, ct, id; 602 PetscInt mask, edge, *iptr; 603 PetscInt *dir, *used; 604 PetscInt sum[4], w[4]; 605 PetscScalar rsum[4], rw[4]; 606 PetscInt op[] = {GL_ADD,0}; 607 PetscScalar *lhs, *rhs; 608 PetscInt *nsep, *lnsep, *fo, nfo=0; 609 PCTFS_gs_ADT PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 610 PetscInt *local2global =xyt_handle->mvi->local2global; 611 PetscInt n =xyt_handle->mvi->n; 612 PetscInt m =xyt_handle->mvi->m; 613 PetscInt level =xyt_handle->level; 614 PetscInt shared =0; 615 616 PetscFunctionBegin; 617 dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 618 nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 619 lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 620 fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 621 used = (PetscInt*)malloc(sizeof(PetscInt)*n); 622 623 PCTFS_ivec_zero(dir,level+1); 624 PCTFS_ivec_zero(nsep,level+1); 625 PCTFS_ivec_zero(lnsep,level+1); 626 PCTFS_ivec_set (fo,-1,n+1); 627 PCTFS_ivec_zero(used,n); 628 629 lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 630 rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 631 632 /* determine the # of unique dof */ 633 PCTFS_rvec_zero(lhs,m); 634 PCTFS_rvec_set(lhs,1.0,n); 635 PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 636 PetscCall(PetscInfo(0,"done first PCTFS_gs_gop_hc\n")); 637 PCTFS_rvec_zero(rsum,2); 638 for (i=0; i<n; i++) { 639 if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; } 640 if (lhs[i]!=1.0) shared=1; 641 } 642 643 PCTFS_grop_hc(rsum,rw,2,op,level); 644 rsum[0]+=0.1; 645 rsum[1]+=0.1; 646 647 xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0]; 648 xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0]; 649 650 /* determine separator sets top down */ 651 if (shared) { 652 /* solution is to do as in the symmetric shared case but then */ 653 /* pick the sub-hc with the most free dofs and do a mat-vec */ 654 /* and pick up the responses on the other sub-hc from the */ 655 /* initial separator set obtained from the symm. shared case */ 656 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!"); 657 /* [dead code deleted since it is unlikely to be completed] */ 658 } else { 659 for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 660 /* set rsh of hc, fire, and collect lhs responses */ 661 (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 662 PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 663 664 /* set lsh of hc, fire, and collect rhs responses */ 665 (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 666 PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 667 668 /* count number of dofs I own that have signal and not in sep set */ 669 for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 670 if (!used[i]) { 671 /* number of unmarked dofs on node */ 672 ct++; 673 /* number of dofs to be marked on lhs hc */ 674 if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++; 675 /* number of dofs to be marked on rhs hc */ 676 if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++; 677 } 678 } 679 680 /* for the non-symmetric case we need separators of width 2 */ 681 /* so take both sides */ 682 (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 683 PCTFS_giop_hc(sum,w,4,op,edge); 684 685 ct=0; 686 if (id<mask) { 687 /* mark dofs I own that have signal and not in sep set */ 688 for (i=0;i<n;i++) { 689 if ((!used[i])&&(lhs[i]!=0.0)) { 690 ct++; nfo++; 691 *--iptr = local2global[i]; 692 used[i] =edge; 693 } 694 } 695 /* LSH hc summation of ct should be sum[0] */ 696 } else { 697 /* mark dofs I own that have signal and not in sep set */ 698 for (i=0; i<n; i++) { 699 if ((!used[i])&&(rhs[i]!=0.0)) { 700 ct++; nfo++; 701 *--iptr = local2global[i]; 702 used[i] = edge; 703 } 704 } 705 /* RSH hc summation of ct should be sum[1] */ 706 } 707 708 if (ct>1) PCTFS_ivec_sort(iptr,ct); 709 lnsep[edge]=ct; 710 nsep[edge] =sum[0]+sum[1]; 711 dir [edge] =BOTH; 712 713 /* LATER or we can recur on these to order seps at this level */ 714 /* do we need full set of separators for this? */ 715 716 /* fold rhs hc into lower */ 717 if (id>=mask) id-=mask; 718 } 719 } 720 721 /* level 0 is on processor case - so mark the remainder */ 722 for (ct=i=0;i<n;i++) { 723 if (!used[i]) { 724 ct++; nfo++; 725 *--iptr = local2global[i]; 726 used[i] = edge; 727 } 728 } 729 if (ct>1) PCTFS_ivec_sort(iptr,ct); 730 lnsep[edge]=ct; 731 nsep [edge]=ct; 732 dir [edge]=BOTH; 733 734 xyt_handle->info->nsep = nsep; 735 xyt_handle->info->lnsep = lnsep; 736 xyt_handle->info->fo = fo; 737 xyt_handle->info->nfo = nfo; 738 739 free(dir); 740 free(lhs); 741 free(rhs); 742 free(used); 743 PetscFunctionReturn(0); 744 } 745 746 /**************************************xyt.c***********************************/ 747 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data) 748 { 749 mv_info *mvi; 750 751 mvi = (mv_info*)malloc(sizeof(mv_info)); 752 mvi->n = n; 753 mvi->m = m; 754 mvi->n_global = -1; 755 mvi->m_global = -1; 756 mvi->local2global= (PetscInt*)malloc((m+1)*sizeof(PetscInt)); 757 758 PCTFS_ivec_copy(mvi->local2global,local2global,m); 759 mvi->local2global[m] = INT_MAX; 760 mvi->matvec = matvec; 761 mvi->grid_data = grid_data; 762 763 /* set xyt communication handle to perform restricted matvec */ 764 mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 765 766 return(mvi); 767 } 768 769 /**************************************xyt.c***********************************/ 770 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 771 { 772 PetscFunctionBegin; 773 A->matvec((mv_info*)A->grid_data,v,u); 774 PetscFunctionReturn(0); 775 } 776