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