1827bd09bSSatish Balay 211cc89d2SBarry Smith /* 3827bd09bSSatish Balay Module Name: xxt 4827bd09bSSatish Balay Module Info: 5827bd09bSSatish Balay 6827bd09bSSatish Balay author: Henry M. Tufo III 7827bd09bSSatish Balay e-mail: hmt@asci.uchicago.edu 8827bd09bSSatish Balay contact: 9827bd09bSSatish Balay +--------------------------------+--------------------------------+ 10827bd09bSSatish Balay |MCS Division - Building 221 |Department of Computer Science | 11827bd09bSSatish Balay |Argonne National Laboratory |Ryerson 152 | 12827bd09bSSatish Balay |9700 S. Cass Avenue |The University of Chicago | 13827bd09bSSatish Balay |Argonne, IL 60439 |Chicago, IL 60637 | 14827bd09bSSatish Balay |(630) 252-5354/5986 ph/fx |(773) 702-6019/8487 ph/fx | 15827bd09bSSatish Balay +--------------------------------+--------------------------------+ 16827bd09bSSatish Balay 17827bd09bSSatish Balay Last Modification: 3.20.01 1811cc89d2SBarry Smith */ 19c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h> 20827bd09bSSatish Balay 21827bd09bSSatish Balay #define LEFT -1 22827bd09bSSatish Balay #define RIGHT 1 23827bd09bSSatish Balay #define BOTH 0 24827bd09bSSatish Balay 25827bd09bSSatish Balay typedef struct xxt_solver_info { 2652f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 2752f87cdaSBarry Smith PetscInt nnz, max_nnz, msg_buf_sz; 2852f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo, *stages; 2952f87cdaSBarry Smith PetscInt *col_sz, *col_indices; 30a501084fSBarry Smith PetscScalar **col_vals, *x, *solve_uu, *solve_w; 3152f87cdaSBarry Smith PetscInt nsolves; 32a501084fSBarry Smith PetscScalar tot_solve_time; 33827bd09bSSatish Balay } xxt_info; 34827bd09bSSatish Balay 35827bd09bSSatish Balay typedef struct matvec_info { 3652f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 3752f87cdaSBarry Smith PetscInt *local2global; 38ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 39a501084fSBarry Smith PetscErrorCode (*matvec)(struct matvec_info *, PetscScalar *, PetscScalar *); 40827bd09bSSatish Balay void *grid_data; 41827bd09bSSatish Balay } mv_info; 42827bd09bSSatish Balay 43827bd09bSSatish Balay struct xxt_CDT { 4452f87cdaSBarry Smith PetscInt id; 4552f87cdaSBarry Smith PetscInt ns; 4652f87cdaSBarry Smith PetscInt level; 47827bd09bSSatish Balay xxt_info *info; 48827bd09bSSatish Balay mv_info *mvi; 49827bd09bSSatish Balay }; 50827bd09bSSatish Balay 5152f87cdaSBarry Smith static PetscInt n_xxt = 0; 5252f87cdaSBarry Smith static PetscInt n_xxt_handles = 0; 53827bd09bSSatish Balay 54827bd09bSSatish Balay /* prototypes */ 553fdc5746SBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *rhs); 563fdc5746SBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle); 573fdc5746SBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle); 583fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 5938dcbb09SBarry Smith static PetscErrorCode xxt_generate(xxt_ADT xxt_handle); 6038dcbb09SBarry Smith static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle); 615c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data); 62a501084fSBarry Smith 63d71ae5a4SJacob Faibussowitsch xxt_ADT XXT_new(void) 64d71ae5a4SJacob Faibussowitsch { 65827bd09bSSatish Balay xxt_ADT xxt_handle; 66827bd09bSSatish Balay 67827bd09bSSatish Balay /* rolling count on n_xxt ... pot. problem here */ 68827bd09bSSatish Balay n_xxt_handles++; 69a501084fSBarry Smith xxt_handle = (xxt_ADT)malloc(sizeof(struct xxt_CDT)); 70827bd09bSSatish Balay xxt_handle->id = ++n_xxt; 719371c9d4SSatish Balay xxt_handle->info = NULL; 729371c9d4SSatish Balay xxt_handle->mvi = NULL; 73827bd09bSSatish Balay 74827bd09bSSatish Balay return (xxt_handle); 75827bd09bSSatish Balay } 76827bd09bSSatish Balay 7738dcbb09SBarry Smith PetscErrorCode XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt handle */ 7852f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */ 7952f87cdaSBarry Smith PetscInt n, /* local num rows */ 8052f87cdaSBarry Smith PetscInt m, /* local num cols */ 815c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc */ 821147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */ 83827bd09bSSatish Balay { 843ba16761SJacob Faibussowitsch PetscFunctionBegin; 853ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 863ba16761SJacob Faibussowitsch PetscCall(check_handle(xxt_handle)); 87827bd09bSSatish Balay 88827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 8963a3b9bcSJacob Faibussowitsch PetscCheck((1 << (xxt_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); 90827bd09bSSatish Balay 91827bd09bSSatish Balay /* space for X info */ 92a501084fSBarry Smith xxt_handle->info = (xxt_info *)malloc(sizeof(xxt_info)); 93827bd09bSSatish Balay 94827bd09bSSatish Balay /* set up matvec handles */ 955c8f6a95SKarl Rupp xxt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode(*)(mv_info *, PetscScalar *, PetscScalar *))matvec, grid_data); 96827bd09bSSatish Balay 97827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 98827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 99827bd09bSSatish Balay xxt_handle->ns = 0; 100827bd09bSSatish Balay 101827bd09bSSatish Balay /* determine separators and generate firing order - NB xxt info set here */ 1023ba16761SJacob Faibussowitsch PetscCall(det_separators(xxt_handle)); 103827bd09bSSatish Balay 1043ba16761SJacob Faibussowitsch PetscCall(do_xxt_factor(xxt_handle)); 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106827bd09bSSatish Balay } 107827bd09bSSatish Balay 108d71ae5a4SJacob Faibussowitsch PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b) 109d71ae5a4SJacob Faibussowitsch { 1103ba16761SJacob Faibussowitsch PetscFunctionBegin; 1113ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 1123ba16761SJacob Faibussowitsch PetscCall(check_handle(xxt_handle)); 113827bd09bSSatish Balay 114827bd09bSSatish Balay /* need to copy b into x? */ 1153ba16761SJacob Faibussowitsch if (b) PetscCall(PCTFS_rvec_copy(x, b, xxt_handle->mvi->n)); 1163ba16761SJacob Faibussowitsch PetscCall(do_xxt_solve(xxt_handle, x)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118827bd09bSSatish Balay } 119827bd09bSSatish Balay 1203ba16761SJacob Faibussowitsch PetscErrorCode XXT_free(xxt_ADT xxt_handle) 121d71ae5a4SJacob Faibussowitsch { 1223ba16761SJacob Faibussowitsch PetscFunctionBegin; 1233ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 1243ba16761SJacob Faibussowitsch PetscCall(check_handle(xxt_handle)); 125827bd09bSSatish Balay n_xxt_handles--; 126827bd09bSSatish Balay 127a501084fSBarry Smith free(xxt_handle->info->nsep); 128a501084fSBarry Smith free(xxt_handle->info->lnsep); 129a501084fSBarry Smith free(xxt_handle->info->fo); 130a501084fSBarry Smith free(xxt_handle->info->stages); 131a501084fSBarry Smith free(xxt_handle->info->solve_uu); 132a501084fSBarry Smith free(xxt_handle->info->solve_w); 133a501084fSBarry Smith free(xxt_handle->info->x); 134a501084fSBarry Smith free(xxt_handle->info->col_vals); 135a501084fSBarry Smith free(xxt_handle->info->col_sz); 136a501084fSBarry Smith free(xxt_handle->info->col_indices); 137a501084fSBarry Smith free(xxt_handle->info); 138a501084fSBarry Smith free(xxt_handle->mvi->local2global); 1393ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle)); 140a501084fSBarry Smith free(xxt_handle->mvi); 141a501084fSBarry Smith free(xxt_handle); 142827bd09bSSatish Balay 143827bd09bSSatish Balay /* if the check fails we nuke */ 144a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 145827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147827bd09bSSatish Balay } 148827bd09bSSatish Balay 14911cc89d2SBarry Smith /* 150827bd09bSSatish Balay 151827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 152827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 153827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 154827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 155827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 156827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 157ca8e9878SJed Brown PCTFS_gs_init/gop). 158827bd09bSSatish Balay 159827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 160827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 16111cc89d2SBarry Smith */ 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle) 163d71ae5a4SJacob Faibussowitsch { 1647b1ae94cSBarry Smith return xxt_generate(xxt_handle); 165827bd09bSSatish Balay } 166827bd09bSSatish Balay 167d71ae5a4SJacob Faibussowitsch static PetscErrorCode xxt_generate(xxt_ADT xxt_handle) 168d71ae5a4SJacob Faibussowitsch { 16952f87cdaSBarry Smith PetscInt i, j, k, idex; 17052f87cdaSBarry Smith PetscInt dim, col; 171a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 17252f87cdaSBarry Smith PetscInt *segs; 17352f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 17452f87cdaSBarry Smith PetscInt off, len; 175a501084fSBarry Smith PetscScalar *x_ptr; 17652f87cdaSBarry Smith PetscInt *iptr, flag; 17752f87cdaSBarry Smith PetscInt start = 0, end, work; 17852f87cdaSBarry Smith PetscInt op2[] = {GL_MIN, 0}; 179ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 18052f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 18152f87cdaSBarry Smith PetscInt a_n = xxt_handle->mvi->n; 18252f87cdaSBarry Smith PetscInt a_m = xxt_handle->mvi->m; 18352f87cdaSBarry Smith PetscInt *a_local2global = xxt_handle->mvi->local2global; 18452f87cdaSBarry Smith PetscInt level; 18552f87cdaSBarry Smith PetscInt xxt_nnz = 0, xxt_max_nnz = 0; 18652f87cdaSBarry Smith PetscInt n, m; 18752f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages; 188a501084fSBarry Smith PetscScalar **col_vals, *x; 18952f87cdaSBarry Smith PetscInt n_global; 19071a0148aSBarry Smith PetscBLASInt i1 = 1, dlen; 191a501084fSBarry Smith PetscScalar dm1 = -1.0; 192827bd09bSSatish Balay 193827bd09bSSatish Balay n = xxt_handle->mvi->n; 194827bd09bSSatish Balay nsep = xxt_handle->info->nsep; 195827bd09bSSatish Balay lnsep = xxt_handle->info->lnsep; 196827bd09bSSatish Balay fo = xxt_handle->info->fo; 197827bd09bSSatish Balay end = lnsep[0]; 198827bd09bSSatish Balay level = xxt_handle->level; 199ca8e9878SJed Brown PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 200827bd09bSSatish Balay 201827bd09bSSatish Balay /* is there a null space? */ 202827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2032fa5cd67SKarl Rupp for (i = 0, j = 0; i <= level; i++) j += nsep[i]; 204827bd09bSSatish Balay 205827bd09bSSatish Balay m = j - xxt_handle->ns; 20648a46eb9SPierre Jolivet if (m != j) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "xxt_generate() :: null space exists %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, j, xxt_handle->ns)); 207827bd09bSSatish Balay 208827bd09bSSatish Balay /* get and initialize storage for x local */ 209827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 21052f87cdaSBarry Smith col_sz = (PetscInt *)malloc(m * sizeof(PetscInt)); 21152f87cdaSBarry Smith col_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt)); 212a501084fSBarry Smith col_vals = (PetscScalar **)malloc(m * sizeof(PetscScalar *)); 213db4deed7SKarl Rupp for (i = j = 0; i < m; i++, j += 2) { 214827bd09bSSatish Balay col_indices[j] = col_indices[j + 1] = col_sz[i] = -1; 215827bd09bSSatish Balay col_vals[i] = NULL; 216827bd09bSSatish Balay } 217827bd09bSSatish Balay col_indices[j] = -1; 218827bd09bSSatish Balay 219827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 220827bd09bSSatish Balay /* this looks like nsep[]=segments */ 22152f87cdaSBarry Smith stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 22252f87cdaSBarry Smith segs = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 2233ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(stages, level + 1)); 224ca8e9878SJed Brown PCTFS_ivec_copy(segs, nsep, level + 1); 2252fa5cd67SKarl Rupp for (i = 0; i < level; i++) segs[i + 1] += segs[i]; 226827bd09bSSatish Balay stages[0] = segs[0]; 227827bd09bSSatish Balay 228827bd09bSSatish Balay /* temporary vectors */ 229a501084fSBarry Smith u = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 230a501084fSBarry Smith z = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 231a501084fSBarry Smith v = (PetscScalar *)malloc(a_m * sizeof(PetscScalar)); 232a501084fSBarry Smith uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 233a501084fSBarry Smith w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 234827bd09bSSatish Balay 235827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 2362fa5cd67SKarl Rupp for (i = 1, j = 0; i <= level; i++) j += nsep[i]; 237827bd09bSSatish Balay 238827bd09bSSatish Balay /* storage for sparse x values */ 239827bd09bSSatish Balay n_global = xxt_handle->info->n_global; 24085ec1a3cSBarry Smith xxt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes; 241a501084fSBarry Smith x = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar)); 242827bd09bSSatish Balay xxt_nnz = 0; 243827bd09bSSatish Balay 244827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 245827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 2462fa5cd67SKarl Rupp for (dim = i = j = 0; i < m; i++) { 247827bd09bSSatish Balay /* time to move to the next level? */ 248d4af0d30SBarry Smith while (i == segs[dim]) { 24908401ef6SPierre Jolivet PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level"); 250827bd09bSSatish Balay stages[dim++] = i; 251827bd09bSSatish Balay end += lnsep[dim]; 252827bd09bSSatish Balay } 253827bd09bSSatish Balay stages[dim] = i; 254827bd09bSSatish Balay 255827bd09bSSatish Balay /* which column are we firing? */ 256827bd09bSSatish Balay /* i.e. set v_l */ 257827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 258827bd09bSSatish Balay (start < end) ? (col = fo[start]) : (col = INT_MAX); 2593ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(&col, &work, 1, op2, dim)); 260827bd09bSSatish Balay 261827bd09bSSatish Balay /* shouldn't need this */ 262db4deed7SKarl Rupp if (col == INT_MAX) { 2639566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n")); 264827bd09bSSatish Balay continue; 265827bd09bSSatish Balay } 266827bd09bSSatish Balay 267827bd09bSSatish Balay /* do I own it? I should */ 2683ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(v, a_m)); 269db4deed7SKarl Rupp if (col == fo[start]) { 270827bd09bSSatish Balay start++; 271ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_n); 2720fdf79fbSJacob Faibussowitsch PetscCheck(idex != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!"); 2739371c9d4SSatish Balay v[idex] = 1.0; 2749371c9d4SSatish Balay j++; 275db4deed7SKarl Rupp } else { 276ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_m); 2772fa5cd67SKarl Rupp if (idex != -1) v[idex] = 1.0; 278827bd09bSSatish Balay } 279827bd09bSSatish Balay 280827bd09bSSatish Balay /* perform u = A.v_l */ 2813ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(u, n)); 2823ba16761SJacob Faibussowitsch PetscCall(do_matvec(xxt_handle->mvi, v, u)); 283827bd09bSSatish Balay 284827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 285827bd09bSSatish Balay /* technically only need to zero out first i entries */ 286827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */ 2873ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uu, m)); 288827bd09bSSatish Balay x_ptr = x; 289827bd09bSSatish Balay iptr = col_indices; 290db4deed7SKarl Rupp for (k = 0; k < i; k++) { 291827bd09bSSatish Balay off = *iptr++; 292827bd09bSSatish Balay len = *iptr++; 2939566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 294792fecdfSBarry Smith PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, x_ptr, &i1)); 295827bd09bSSatish Balay x_ptr += len; 296827bd09bSSatish Balay } 297827bd09bSSatish Balay 298827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 2999566063dSJacob Faibussowitsch PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages)); 300827bd09bSSatish Balay 301827bd09bSSatish Balay /* z = X.uu */ 3023ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(z, n)); 303827bd09bSSatish Balay x_ptr = x; 304827bd09bSSatish Balay iptr = col_indices; 305db4deed7SKarl Rupp for (k = 0; k < i; k++) { 306827bd09bSSatish Balay off = *iptr++; 307827bd09bSSatish Balay len = *iptr++; 3089566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 309792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1)); 310827bd09bSSatish Balay x_ptr += len; 311827bd09bSSatish Balay } 312827bd09bSSatish Balay 313827bd09bSSatish Balay /* compute v_l = v_l - z */ 3143ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(v + a_n, a_m - a_n)); 3159566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 316792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1)); 317827bd09bSSatish Balay 318827bd09bSSatish Balay /* compute u_l = A.v_l */ 3193ba16761SJacob Faibussowitsch if (a_n != a_m) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim)); 3203ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(u, n)); 3213ba16761SJacob Faibussowitsch PetscCall(do_matvec(xxt_handle->mvi, v, u)); 322827bd09bSSatish Balay 323827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */ 3249566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 325792fecdfSBarry Smith PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, v, &i1)); 326827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */ 3273ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim)); 328827bd09bSSatish Balay 3298f1a2a5eSBarry Smith alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha); 330827bd09bSSatish Balay 331827bd09bSSatish Balay /* check for small alpha */ 332827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 33363a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha)); 334827bd09bSSatish Balay 335827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 3363ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_scale(v, 1.0 / alpha, n)); 337827bd09bSSatish Balay 338827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 339827bd09bSSatish Balay flag = 1; 340827bd09bSSatish Balay off = len = 0; 341db4deed7SKarl Rupp for (k = 0; k < n; k++) { 342db4deed7SKarl Rupp if (v[k] != 0.0) { 343827bd09bSSatish Balay len = k; 3449371c9d4SSatish Balay if (flag) { 3459371c9d4SSatish Balay off = k; 3469371c9d4SSatish Balay flag = 0; 3479371c9d4SSatish Balay } 348827bd09bSSatish Balay } 349827bd09bSSatish Balay } 350827bd09bSSatish Balay 351827bd09bSSatish Balay len -= (off - 1); 352827bd09bSSatish Balay 353db4deed7SKarl Rupp if (len > 0) { 354db4deed7SKarl Rupp if ((xxt_nnz + len) > xxt_max_nnz) { 3559566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "increasing space for X by 2x!\n")); 356827bd09bSSatish Balay xxt_max_nnz *= 2; 357a501084fSBarry Smith x_ptr = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar)); 3583ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_copy(x_ptr, x, xxt_nnz)); 359a501084fSBarry Smith free(x); 360827bd09bSSatish Balay x = x_ptr; 361827bd09bSSatish Balay x_ptr += xxt_nnz; 362827bd09bSSatish Balay } 363827bd09bSSatish Balay xxt_nnz += len; 3643ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_copy(x_ptr, v + off, len)); 365827bd09bSSatish Balay 366827bd09bSSatish Balay col_indices[2 * i] = off; 367827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = len; 368827bd09bSSatish Balay col_vals[i] = x_ptr; 3699371c9d4SSatish Balay } else { 370827bd09bSSatish Balay col_indices[2 * i] = 0; 371827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = 0; 372827bd09bSSatish Balay col_vals[i] = x_ptr; 373827bd09bSSatish Balay } 374827bd09bSSatish Balay } 375827bd09bSSatish Balay 376827bd09bSSatish Balay /* close off stages for execution phase */ 377db4deed7SKarl Rupp while (dim != level) { 378827bd09bSSatish Balay stages[dim++] = i; 37963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level)); 380827bd09bSSatish Balay } 381827bd09bSSatish Balay stages[dim] = i; 382827bd09bSSatish Balay 383827bd09bSSatish Balay xxt_handle->info->n = xxt_handle->mvi->n; 384827bd09bSSatish Balay xxt_handle->info->m = m; 385827bd09bSSatish Balay xxt_handle->info->nnz = xxt_nnz; 386827bd09bSSatish Balay xxt_handle->info->max_nnz = xxt_max_nnz; 387827bd09bSSatish Balay xxt_handle->info->msg_buf_sz = stages[level] - stages[0]; 388a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 389a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 390827bd09bSSatish Balay xxt_handle->info->x = x; 391827bd09bSSatish Balay xxt_handle->info->col_vals = col_vals; 392827bd09bSSatish Balay xxt_handle->info->col_sz = col_sz; 393827bd09bSSatish Balay xxt_handle->info->col_indices = col_indices; 394827bd09bSSatish Balay xxt_handle->info->stages = stages; 395827bd09bSSatish Balay xxt_handle->info->nsolves = 0; 396827bd09bSSatish Balay xxt_handle->info->tot_solve_time = 0.0; 397827bd09bSSatish Balay 398a501084fSBarry Smith free(segs); 399a501084fSBarry Smith free(u); 400a501084fSBarry Smith free(v); 401a501084fSBarry Smith free(uu); 402a501084fSBarry Smith free(z); 403a501084fSBarry Smith free(w); 404827bd09bSSatish Balay 4053ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 406827bd09bSSatish Balay } 407827bd09bSSatish Balay 408d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc) 409d71ae5a4SJacob Faibussowitsch { 41052f87cdaSBarry Smith PetscInt off, len, *iptr; 41152f87cdaSBarry Smith PetscInt level = xxt_handle->level; 41252f87cdaSBarry Smith PetscInt n = xxt_handle->info->n; 41352f87cdaSBarry Smith PetscInt m = xxt_handle->info->m; 41452f87cdaSBarry Smith PetscInt *stages = xxt_handle->info->stages; 41552f87cdaSBarry Smith PetscInt *col_indices = xxt_handle->info->col_indices; 416a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr; 417a501084fSBarry Smith PetscScalar *solve_uu = xxt_handle->info->solve_uu; 418a501084fSBarry Smith PetscScalar *solve_w = xxt_handle->info->solve_w; 419a501084fSBarry Smith PetscScalar *x = xxt_handle->info->x; 42071a0148aSBarry Smith PetscBLASInt i1 = 1, dlen; 421827bd09bSSatish Balay 4220924e98cSBarry Smith PetscFunctionBegin; 423827bd09bSSatish Balay uu_ptr = solve_uu; 4243ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uu_ptr, m)); 425827bd09bSSatish Balay 426827bd09bSSatish Balay /* x = X.Y^T.b */ 427827bd09bSSatish Balay /* uu = Y^T.b */ 428db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) { 429c5df96a5SBarry Smith off = *iptr++; 430c5df96a5SBarry Smith len = *iptr++; 4319566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 432792fecdfSBarry Smith PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, x_ptr, &i1)); 433827bd09bSSatish Balay } 434827bd09bSSatish Balay 435d5b43468SJose E. Roman /* communication of beta */ 436827bd09bSSatish Balay uu_ptr = solve_uu; 4379566063dSJacob Faibussowitsch if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages)); 438827bd09bSSatish Balay 4393ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uc, n)); 440827bd09bSSatish Balay 441827bd09bSSatish Balay /* x = X.uu */ 442db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) { 443c5df96a5SBarry Smith off = *iptr++; 444c5df96a5SBarry Smith len = *iptr++; 4459566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 446792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1)); 447827bd09bSSatish Balay } 4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 449827bd09bSSatish Balay } 450827bd09bSSatish Balay 451d71ae5a4SJacob Faibussowitsch static PetscErrorCode check_handle(xxt_ADT xxt_handle) 452d71ae5a4SJacob Faibussowitsch { 45352f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX}; 454827bd09bSSatish Balay 4550924e98cSBarry Smith PetscFunctionBegin; 4567ee088dcSPierre Jolivet PetscCheck(xxt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", (void *)xxt_handle); 457827bd09bSSatish Balay 458827bd09bSSatish Balay vals[0] = vals[1] = xxt_handle->id; 4593ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op)); 46063a3b9bcSJacob Faibussowitsch PetscCheck(!(vals[0] != vals[1]) && !(xxt_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], xxt_handle->id); 4613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 462827bd09bSSatish Balay } 463827bd09bSSatish Balay 464d71ae5a4SJacob Faibussowitsch static PetscErrorCode det_separators(xxt_ADT xxt_handle) 465d71ae5a4SJacob Faibussowitsch { 46652f87cdaSBarry Smith PetscInt i, ct, id; 46752f87cdaSBarry Smith PetscInt mask, edge, *iptr; 46852f87cdaSBarry Smith PetscInt *dir, *used; 46952f87cdaSBarry Smith PetscInt sum[4], w[4]; 470a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 47152f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 472a501084fSBarry Smith PetscScalar *lhs, *rhs; 47352f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo = 0; 474ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 47552f87cdaSBarry Smith PetscInt *local2global = xxt_handle->mvi->local2global; 47652f87cdaSBarry Smith PetscInt n = xxt_handle->mvi->n; 47752f87cdaSBarry Smith PetscInt m = xxt_handle->mvi->m; 47852f87cdaSBarry Smith PetscInt level = xxt_handle->level; 479ab824b78SBarry Smith PetscInt shared = 0; 480827bd09bSSatish Balay 4810924e98cSBarry Smith PetscFunctionBegin; 48252f87cdaSBarry Smith dir = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 48352f87cdaSBarry Smith nsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 48452f87cdaSBarry Smith lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 48552f87cdaSBarry Smith fo = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1)); 48652f87cdaSBarry Smith used = (PetscInt *)malloc(sizeof(PetscInt) * n); 487827bd09bSSatish Balay 4883ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(dir, level + 1)); 4893ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(nsep, level + 1)); 4903ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(lnsep, level + 1)); 4913ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_set(fo, -1, n + 1)); 4923ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(used, n)); 493827bd09bSSatish Balay 4948cda6cd7SBarry Smith lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 4958cda6cd7SBarry Smith rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 496827bd09bSSatish Balay 497827bd09bSSatish Balay /* determine the # of unique dof */ 4983ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(lhs, m)); 4993ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_set(lhs, 1.0, n)); 5003ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level)); 5013ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(rsum, 2)); 5023d3eaba7SBarry Smith for (i = 0; i < n; i++) { 5032fa5cd67SKarl Rupp if (lhs[i] != 0.0) { 5049371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 5059371c9d4SSatish Balay rsum[1] += lhs[i]; 5062fa5cd67SKarl Rupp } 507827bd09bSSatish Balay } 5083ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(rsum, rw, 2, op, level)); 509827bd09bSSatish Balay rsum[0] += 0.1; 510827bd09bSSatish Balay rsum[1] += 0.1; 511827bd09bSSatish Balay 51277b4d14cSPeter Brune if (PetscAbsScalar(rsum[0] - rsum[1]) > EPS) shared = 1; 513827bd09bSSatish Balay 51452f87cdaSBarry Smith xxt_handle->info->n_global = xxt_handle->info->m_global = (PetscInt)rsum[0]; 51552f87cdaSBarry Smith xxt_handle->mvi->n_global = xxt_handle->mvi->m_global = (PetscInt)rsum[0]; 516827bd09bSSatish Balay 517827bd09bSSatish Balay /* determine separator sets top down */ 5182fa5cd67SKarl Rupp if (shared) { 519db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 520827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 5213ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m)); 5223ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge)); 523827bd09bSSatish Balay 524827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 5253ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m)); 5263ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge)); 527827bd09bSSatish Balay 528db4deed7SKarl Rupp for (i = 0; i < n; i++) { 529db4deed7SKarl Rupp if (id < mask) { 5302fa5cd67SKarl Rupp if (lhs[i] != 0.0) lhs[i] = 1.0; 531827bd09bSSatish Balay } 532db4deed7SKarl Rupp if (id >= mask) { 5332fa5cd67SKarl Rupp if (rhs[i] != 0.0) rhs[i] = 1.0; 534827bd09bSSatish Balay } 535827bd09bSSatish Balay } 536827bd09bSSatish Balay 5373ba16761SJacob Faibussowitsch if (id < mask) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge - 1)); 5383ba16761SJacob Faibussowitsch else PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge - 1)); 539827bd09bSSatish Balay 540827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 5413ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(rsum, 4)); 5423ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(sum, 4)); 5433ba16761SJacob Faibussowitsch for (ct = i = 0; i < n; i++) { 544db4deed7SKarl Rupp if (!used[i]) { 545827bd09bSSatish Balay /* number of unmarked dofs on node */ 546827bd09bSSatish Balay ct++; 547827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 548db4deed7SKarl Rupp if (id < mask) { 5499371c9d4SSatish Balay if (lhs[i] != 0.0) { 5509371c9d4SSatish Balay sum[0]++; 5519371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 5529371c9d4SSatish Balay } 553827bd09bSSatish Balay } 554827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 555db4deed7SKarl Rupp if (id >= mask) { 5569371c9d4SSatish Balay if (rhs[i] != 0.0) { 5579371c9d4SSatish Balay sum[1]++; 5589371c9d4SSatish Balay rsum[1] += 1.0 / rhs[i]; 5599371c9d4SSatish Balay } 560827bd09bSSatish Balay } 561827bd09bSSatish Balay } 562827bd09bSSatish Balay } 563827bd09bSSatish Balay 564827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 565827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 566827bd09bSSatish Balay (id < mask) ? (rsum[2] = ct) : (rsum[3] = ct); 5673ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge)); 5683ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(rsum, rw, 4, op, edge)); 5699371c9d4SSatish Balay rsum[0] += 0.1; 5709371c9d4SSatish Balay rsum[1] += 0.1; 5719371c9d4SSatish Balay rsum[2] += 0.1; 5729371c9d4SSatish Balay rsum[3] += 0.1; 573827bd09bSSatish Balay 574db4deed7SKarl Rupp if (id < mask) { 575827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 576db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 577db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 5789371c9d4SSatish Balay ct++; 5799371c9d4SSatish Balay nfo++; 580827bd09bSSatish Balay 58108401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n"); 582827bd09bSSatish Balay 583827bd09bSSatish Balay *--iptr = local2global[i]; 584827bd09bSSatish Balay used[i] = edge; 585827bd09bSSatish Balay } 586827bd09bSSatish Balay } 5873ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 588827bd09bSSatish Balay 589827bd09bSSatish Balay lnsep[edge] = ct; 59052f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[0]; 591827bd09bSSatish Balay dir[edge] = LEFT; 592827bd09bSSatish Balay } 593827bd09bSSatish Balay 594db4deed7SKarl Rupp if (id >= mask) { 595827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 596db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 597db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 5989371c9d4SSatish Balay ct++; 5999371c9d4SSatish Balay nfo++; 600827bd09bSSatish Balay 60108401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n"); 602827bd09bSSatish Balay 603827bd09bSSatish Balay *--iptr = local2global[i]; 604827bd09bSSatish Balay used[i] = edge; 605827bd09bSSatish Balay } 606827bd09bSSatish Balay } 6073ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 608827bd09bSSatish Balay 609827bd09bSSatish Balay lnsep[edge] = ct; 61052f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[1]; 611827bd09bSSatish Balay dir[edge] = RIGHT; 612827bd09bSSatish Balay } 613827bd09bSSatish Balay 614827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 615827bd09bSSatish Balay /* do we need full set of separators for this? */ 616827bd09bSSatish Balay 617827bd09bSSatish Balay /* fold rhs hc into lower */ 6182fa5cd67SKarl Rupp if (id >= mask) id -= mask; 619827bd09bSSatish Balay } 620db4deed7SKarl Rupp } else { 621db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 622827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 6233ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m)); 6243ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge)); 625827bd09bSSatish Balay 626827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 6273ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m)); 6283ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge)); 629827bd09bSSatish Balay 630827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 6313ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(sum, 4)); 6323ba16761SJacob Faibussowitsch for (ct = i = 0; i < n; i++) { 633db4deed7SKarl Rupp if (!used[i]) { 634827bd09bSSatish Balay /* number of unmarked dofs on node */ 635827bd09bSSatish Balay ct++; 636827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6372fa5cd67SKarl Rupp if ((id < mask) && (lhs[i] != 0.0)) sum[0]++; 638827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6392fa5cd67SKarl Rupp if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++; 640827bd09bSSatish Balay } 641827bd09bSSatish Balay } 642827bd09bSSatish Balay 643827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 644827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 6453ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge)); 646827bd09bSSatish Balay 647827bd09bSSatish Balay /* lhs hc wins */ 648db4deed7SKarl Rupp if (sum[2] >= sum[3]) { 649db4deed7SKarl Rupp if (id < mask) { 650827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 651db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 652db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 6539371c9d4SSatish Balay ct++; 6549371c9d4SSatish Balay nfo++; 655827bd09bSSatish Balay *--iptr = local2global[i]; 656827bd09bSSatish Balay used[i] = edge; 657827bd09bSSatish Balay } 658827bd09bSSatish Balay } 6593ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 660827bd09bSSatish Balay lnsep[edge] = ct; 661827bd09bSSatish Balay } 662827bd09bSSatish Balay nsep[edge] = sum[0]; 663827bd09bSSatish Balay dir[edge] = LEFT; 664db4deed7SKarl Rupp } else { /* rhs hc wins */ 665db4deed7SKarl Rupp if (id >= mask) { 666827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 667db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 668db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 6699371c9d4SSatish Balay ct++; 6709371c9d4SSatish Balay nfo++; 671827bd09bSSatish Balay *--iptr = local2global[i]; 672827bd09bSSatish Balay used[i] = edge; 673827bd09bSSatish Balay } 674827bd09bSSatish Balay } 6753ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 676827bd09bSSatish Balay lnsep[edge] = ct; 677827bd09bSSatish Balay } 678827bd09bSSatish Balay nsep[edge] = sum[1]; 679827bd09bSSatish Balay dir[edge] = RIGHT; 680827bd09bSSatish Balay } 681827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 682827bd09bSSatish Balay /* do we need full set of separators for this? */ 683827bd09bSSatish Balay 684827bd09bSSatish Balay /* fold rhs hc into lower */ 6852fa5cd67SKarl Rupp if (id >= mask) id -= mask; 686827bd09bSSatish Balay } 687827bd09bSSatish Balay } 688827bd09bSSatish Balay 689827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 690db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 691db4deed7SKarl Rupp if (!used[i]) { 6929371c9d4SSatish Balay ct++; 6939371c9d4SSatish Balay nfo++; 694827bd09bSSatish Balay *--iptr = local2global[i]; 695827bd09bSSatish Balay used[i] = edge; 696827bd09bSSatish Balay } 697827bd09bSSatish Balay } 6983ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 699827bd09bSSatish Balay lnsep[edge] = ct; 700827bd09bSSatish Balay nsep[edge] = ct; 701827bd09bSSatish Balay dir[edge] = LEFT; 702827bd09bSSatish Balay 703827bd09bSSatish Balay xxt_handle->info->nsep = nsep; 704827bd09bSSatish Balay xxt_handle->info->lnsep = lnsep; 705827bd09bSSatish Balay xxt_handle->info->fo = fo; 706827bd09bSSatish Balay xxt_handle->info->nfo = nfo; 707827bd09bSSatish Balay 708a501084fSBarry Smith free(dir); 709a501084fSBarry Smith free(lhs); 710a501084fSBarry Smith free(rhs); 711a501084fSBarry Smith free(used); 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 713827bd09bSSatish Balay } 714827bd09bSSatish Balay 715d71ae5a4SJacob Faibussowitsch static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data) 716d71ae5a4SJacob Faibussowitsch { 717*ff6a9541SJacob Faibussowitsch mv_info *mvi = (mv_info *)malloc(sizeof(mv_info)); 718827bd09bSSatish Balay 719827bd09bSSatish Balay mvi->n = n; 720827bd09bSSatish Balay mvi->m = m; 721827bd09bSSatish Balay mvi->n_global = -1; 722827bd09bSSatish Balay mvi->m_global = -1; 72352f87cdaSBarry Smith mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt)); 724ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global, local2global, m); 725827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7265c8f6a95SKarl Rupp mvi->matvec = matvec; 727827bd09bSSatish Balay mvi->grid_data = grid_data; 728827bd09bSSatish Balay 729827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */ 730ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 731827bd09bSSatish Balay 732827bd09bSSatish Balay return (mvi); 733827bd09bSSatish Balay } 734827bd09bSSatish Balay 735d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 736d71ae5a4SJacob Faibussowitsch { 7370924e98cSBarry Smith PetscFunctionBegin; 7383ba16761SJacob Faibussowitsch PetscCall(A->matvec((mv_info *)A->grid_data, v, u)); 7393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 740827bd09bSSatish Balay } 741