1827bd09bSSatish Balay 2*11cc89d2SBarry 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 18*11cc89d2SBarry 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 639371c9d4SSatish Balay xxt_ADT XXT_new(void) { 64827bd09bSSatish Balay xxt_ADT xxt_handle; 65827bd09bSSatish Balay 66827bd09bSSatish Balay /* rolling count on n_xxt ... pot. problem here */ 67827bd09bSSatish Balay n_xxt_handles++; 68a501084fSBarry Smith xxt_handle = (xxt_ADT)malloc(sizeof(struct xxt_CDT)); 69827bd09bSSatish Balay xxt_handle->id = ++n_xxt; 709371c9d4SSatish Balay xxt_handle->info = NULL; 719371c9d4SSatish Balay xxt_handle->mvi = NULL; 72827bd09bSSatish Balay 73827bd09bSSatish Balay return (xxt_handle); 74827bd09bSSatish Balay } 75827bd09bSSatish Balay 7638dcbb09SBarry Smith PetscErrorCode XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt handle */ 7752f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */ 7852f87cdaSBarry Smith PetscInt n, /* local num rows */ 7952f87cdaSBarry Smith PetscInt m, /* local num cols */ 805c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc */ 811147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */ 82827bd09bSSatish Balay { 83b1c944f5SJed Brown PCTFS_comm_init(); 84827bd09bSSatish Balay check_handle(xxt_handle); 85827bd09bSSatish Balay 86827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 8763a3b9bcSJacob 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); 88827bd09bSSatish Balay 89827bd09bSSatish Balay /* space for X info */ 90a501084fSBarry Smith xxt_handle->info = (xxt_info *)malloc(sizeof(xxt_info)); 91827bd09bSSatish Balay 92827bd09bSSatish Balay /* set up matvec handles */ 935c8f6a95SKarl Rupp xxt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode(*)(mv_info *, PetscScalar *, PetscScalar *))matvec, grid_data); 94827bd09bSSatish Balay 95827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 96827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 97827bd09bSSatish Balay xxt_handle->ns = 0; 98827bd09bSSatish Balay 99827bd09bSSatish Balay /* determine separators and generate firing order - NB xxt info set here */ 100827bd09bSSatish Balay det_separators(xxt_handle); 101827bd09bSSatish Balay 102827bd09bSSatish Balay return (do_xxt_factor(xxt_handle)); 103827bd09bSSatish Balay } 104827bd09bSSatish Balay 1059371c9d4SSatish Balay PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b) { 106b1c944f5SJed Brown PCTFS_comm_init(); 107827bd09bSSatish Balay check_handle(xxt_handle); 108827bd09bSSatish Balay 109827bd09bSSatish Balay /* need to copy b into x? */ 1102fa5cd67SKarl Rupp if (b) PCTFS_rvec_copy(x, b, xxt_handle->mvi->n); 11138dcbb09SBarry Smith return do_xxt_solve(xxt_handle, x); 112827bd09bSSatish Balay } 113827bd09bSSatish Balay 1149371c9d4SSatish Balay PetscInt XXT_free(xxt_ADT xxt_handle) { 115b1c944f5SJed Brown PCTFS_comm_init(); 116827bd09bSSatish Balay check_handle(xxt_handle); 117827bd09bSSatish Balay n_xxt_handles--; 118827bd09bSSatish Balay 119a501084fSBarry Smith free(xxt_handle->info->nsep); 120a501084fSBarry Smith free(xxt_handle->info->lnsep); 121a501084fSBarry Smith free(xxt_handle->info->fo); 122a501084fSBarry Smith free(xxt_handle->info->stages); 123a501084fSBarry Smith free(xxt_handle->info->solve_uu); 124a501084fSBarry Smith free(xxt_handle->info->solve_w); 125a501084fSBarry Smith free(xxt_handle->info->x); 126a501084fSBarry Smith free(xxt_handle->info->col_vals); 127a501084fSBarry Smith free(xxt_handle->info->col_sz); 128a501084fSBarry Smith free(xxt_handle->info->col_indices); 129a501084fSBarry Smith free(xxt_handle->info); 130a501084fSBarry Smith free(xxt_handle->mvi->local2global); 131ca8e9878SJed Brown PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle); 132a501084fSBarry Smith free(xxt_handle->mvi); 133a501084fSBarry Smith free(xxt_handle); 134827bd09bSSatish Balay 135827bd09bSSatish Balay /* if the check fails we nuke */ 136a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 137827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 138827bd09bSSatish Balay return (0); 139827bd09bSSatish Balay } 140827bd09bSSatish Balay 141*11cc89d2SBarry Smith /* This function is currently unused */ 1429371c9d4SSatish Balay PetscErrorCode XXT_stats(xxt_ADT xxt_handle) { 14352f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD}; 14452f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD}; 14552f87cdaSBarry Smith PetscInt vals[9], work[9]; 146a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 147827bd09bSSatish Balay 14863a3b9bcSJacob Faibussowitsch PetscFunctionBegin; 14963a3b9bcSJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 15063a3b9bcSJacob Faibussowitsch PetscCall(check_handle(xxt_handle)); 151827bd09bSSatish Balay 152827bd09bSSatish Balay /* if factorization not done there are no stats */ 153db4deed7SKarl Rupp if (!xxt_handle->info || !xxt_handle->mvi) { 1549566063dSJacob Faibussowitsch if (!PCTFS_my_id) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "XXT_stats() :: no stats available!\n")); 155*11cc89d2SBarry Smith PetscFunctionReturn(0); 156827bd09bSSatish Balay } 157827bd09bSSatish Balay 158827bd09bSSatish Balay vals[0] = vals[1] = vals[2] = xxt_handle->info->nnz; 159827bd09bSSatish Balay vals[3] = vals[4] = vals[5] = xxt_handle->mvi->n; 160827bd09bSSatish Balay vals[6] = vals[7] = vals[8] = xxt_handle->info->msg_buf_sz; 161dd39110bSPierre Jolivet PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op); 162827bd09bSSatish Balay 1632fa5cd67SKarl Rupp fvals[0] = fvals[1] = fvals[2] = xxt_handle->info->tot_solve_time / xxt_handle->info->nsolves++; 164dd39110bSPierre Jolivet PCTFS_grop(fvals, fwork, PETSC_STATIC_ARRAY_LENGTH(fop) - 1, fop); 165827bd09bSSatish Balay 166db4deed7SKarl Rupp if (!PCTFS_my_id) { 16763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[0])); 16863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[1])); 16963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_nnz=%g\n", PCTFS_my_id, (double)(1.0 * vals[2] / PCTFS_num_nodes))); 17063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[2])); 17163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: xxt C(2d) =%g\n", PCTFS_my_id, (double)(vals[2] / (PetscPowReal(1.0 * vals[5], 1.5))))); 17263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: xxt C(3d) =%g\n", PCTFS_my_id, (double)(vals[2] / (PetscPowReal(1.0 * vals[5], 1.6667))))); 17363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[3])); 17463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[4])); 17563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_n =%g\n", PCTFS_my_id, (double)(1.0 * vals[5] / PCTFS_num_nodes))); 17663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot xxt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[5])); 17763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[6])); 17863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[7])); 17963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_buf=%g\n", PCTFS_my_id, (double)(1.0 * vals[8] / PCTFS_num_nodes))); 18063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[0]))); 18163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[1]))); 18263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[2] / PCTFS_num_nodes))); 183827bd09bSSatish Balay } 18463a3b9bcSJacob Faibussowitsch PetscFunctionReturn(0); 185827bd09bSSatish Balay } 186827bd09bSSatish Balay 187*11cc89d2SBarry Smith /* 188827bd09bSSatish Balay 189827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 190827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 191827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 192827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 193827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 194827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 195ca8e9878SJed Brown PCTFS_gs_init/gop). 196827bd09bSSatish Balay 197827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 198827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 199*11cc89d2SBarry Smith */ 2009371c9d4SSatish Balay static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle) { 2017b1ae94cSBarry Smith return xxt_generate(xxt_handle); 202827bd09bSSatish Balay } 203827bd09bSSatish Balay 2049371c9d4SSatish Balay static PetscErrorCode xxt_generate(xxt_ADT xxt_handle) { 20552f87cdaSBarry Smith PetscInt i, j, k, idex; 20652f87cdaSBarry Smith PetscInt dim, col; 207a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 20852f87cdaSBarry Smith PetscInt *segs; 20952f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 21052f87cdaSBarry Smith PetscInt off, len; 211a501084fSBarry Smith PetscScalar *x_ptr; 21252f87cdaSBarry Smith PetscInt *iptr, flag; 21352f87cdaSBarry Smith PetscInt start = 0, end, work; 21452f87cdaSBarry Smith PetscInt op2[] = {GL_MIN, 0}; 215ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 21652f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 21752f87cdaSBarry Smith PetscInt a_n = xxt_handle->mvi->n; 21852f87cdaSBarry Smith PetscInt a_m = xxt_handle->mvi->m; 21952f87cdaSBarry Smith PetscInt *a_local2global = xxt_handle->mvi->local2global; 22052f87cdaSBarry Smith PetscInt level; 22152f87cdaSBarry Smith PetscInt xxt_nnz = 0, xxt_max_nnz = 0; 22252f87cdaSBarry Smith PetscInt n, m; 22352f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages; 224a501084fSBarry Smith PetscScalar **col_vals, *x; 22552f87cdaSBarry Smith PetscInt n_global; 22671a0148aSBarry Smith PetscBLASInt i1 = 1, dlen; 227a501084fSBarry Smith PetscScalar dm1 = -1.0; 228827bd09bSSatish Balay 229827bd09bSSatish Balay n = xxt_handle->mvi->n; 230827bd09bSSatish Balay nsep = xxt_handle->info->nsep; 231827bd09bSSatish Balay lnsep = xxt_handle->info->lnsep; 232827bd09bSSatish Balay fo = xxt_handle->info->fo; 233827bd09bSSatish Balay end = lnsep[0]; 234827bd09bSSatish Balay level = xxt_handle->level; 235ca8e9878SJed Brown PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 236827bd09bSSatish Balay 237827bd09bSSatish Balay /* is there a null space? */ 238827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2392fa5cd67SKarl Rupp for (i = 0, j = 0; i <= level; i++) j += nsep[i]; 240827bd09bSSatish Balay 241827bd09bSSatish Balay m = j - xxt_handle->ns; 24248a46eb9SPierre 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)); 243827bd09bSSatish Balay 244827bd09bSSatish Balay /* get and initialize storage for x local */ 245827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 24652f87cdaSBarry Smith col_sz = (PetscInt *)malloc(m * sizeof(PetscInt)); 24752f87cdaSBarry Smith col_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt)); 248a501084fSBarry Smith col_vals = (PetscScalar **)malloc(m * sizeof(PetscScalar *)); 249db4deed7SKarl Rupp for (i = j = 0; i < m; i++, j += 2) { 250827bd09bSSatish Balay col_indices[j] = col_indices[j + 1] = col_sz[i] = -1; 251827bd09bSSatish Balay col_vals[i] = NULL; 252827bd09bSSatish Balay } 253827bd09bSSatish Balay col_indices[j] = -1; 254827bd09bSSatish Balay 255827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 256827bd09bSSatish Balay /* this looks like nsep[]=segments */ 25752f87cdaSBarry Smith stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 25852f87cdaSBarry Smith segs = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 259ca8e9878SJed Brown PCTFS_ivec_zero(stages, level + 1); 260ca8e9878SJed Brown PCTFS_ivec_copy(segs, nsep, level + 1); 2612fa5cd67SKarl Rupp for (i = 0; i < level; i++) segs[i + 1] += segs[i]; 262827bd09bSSatish Balay stages[0] = segs[0]; 263827bd09bSSatish Balay 264827bd09bSSatish Balay /* temporary vectors */ 265a501084fSBarry Smith u = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 266a501084fSBarry Smith z = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 267a501084fSBarry Smith v = (PetscScalar *)malloc(a_m * sizeof(PetscScalar)); 268a501084fSBarry Smith uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 269a501084fSBarry Smith w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 270827bd09bSSatish Balay 271827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 2722fa5cd67SKarl Rupp for (i = 1, j = 0; i <= level; i++) j += nsep[i]; 273827bd09bSSatish Balay 274827bd09bSSatish Balay /* storage for sparse x values */ 275827bd09bSSatish Balay n_global = xxt_handle->info->n_global; 27685ec1a3cSBarry Smith xxt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes; 277a501084fSBarry Smith x = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar)); 278827bd09bSSatish Balay xxt_nnz = 0; 279827bd09bSSatish Balay 280827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 281827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 2822fa5cd67SKarl Rupp for (dim = i = j = 0; i < m; i++) { 283827bd09bSSatish Balay /* time to move to the next level? */ 284d4af0d30SBarry Smith while (i == segs[dim]) { 28508401ef6SPierre Jolivet PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level"); 286827bd09bSSatish Balay stages[dim++] = i; 287827bd09bSSatish Balay end += lnsep[dim]; 288827bd09bSSatish Balay } 289827bd09bSSatish Balay stages[dim] = i; 290827bd09bSSatish Balay 291827bd09bSSatish Balay /* which column are we firing? */ 292827bd09bSSatish Balay /* i.e. set v_l */ 293827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 294827bd09bSSatish Balay (start < end) ? (col = fo[start]) : (col = INT_MAX); 295b1c944f5SJed Brown PCTFS_giop_hc(&col, &work, 1, op2, dim); 296827bd09bSSatish Balay 297827bd09bSSatish Balay /* shouldn't need this */ 298db4deed7SKarl Rupp if (col == INT_MAX) { 2999566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n")); 300827bd09bSSatish Balay continue; 301827bd09bSSatish Balay } 302827bd09bSSatish Balay 303827bd09bSSatish Balay /* do I own it? I should */ 304ca8e9878SJed Brown PCTFS_rvec_zero(v, a_m); 305db4deed7SKarl Rupp if (col == fo[start]) { 306827bd09bSSatish Balay start++; 307ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_n); 308e7e72b3dSBarry Smith if (idex != -1) { 3099371c9d4SSatish Balay v[idex] = 1.0; 3109371c9d4SSatish Balay j++; 311546078acSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!"); 312db4deed7SKarl Rupp } else { 313ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_m); 3142fa5cd67SKarl Rupp if (idex != -1) v[idex] = 1.0; 315827bd09bSSatish Balay } 316827bd09bSSatish Balay 317827bd09bSSatish Balay /* perform u = A.v_l */ 318ca8e9878SJed Brown PCTFS_rvec_zero(u, n); 319827bd09bSSatish Balay do_matvec(xxt_handle->mvi, v, u); 320827bd09bSSatish Balay 321827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 322827bd09bSSatish Balay /* technically only need to zero out first i entries */ 323827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */ 324ca8e9878SJed Brown PCTFS_rvec_zero(uu, m); 325827bd09bSSatish Balay x_ptr = x; 326827bd09bSSatish Balay iptr = col_indices; 327db4deed7SKarl Rupp for (k = 0; k < i; k++) { 328827bd09bSSatish Balay off = *iptr++; 329827bd09bSSatish Balay len = *iptr++; 3309566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 331792fecdfSBarry Smith PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, x_ptr, &i1)); 332827bd09bSSatish Balay x_ptr += len; 333827bd09bSSatish Balay } 334827bd09bSSatish Balay 335827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 3369566063dSJacob Faibussowitsch PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages)); 337827bd09bSSatish Balay 338827bd09bSSatish Balay /* z = X.uu */ 339ca8e9878SJed Brown PCTFS_rvec_zero(z, n); 340827bd09bSSatish Balay x_ptr = x; 341827bd09bSSatish Balay iptr = col_indices; 342db4deed7SKarl Rupp for (k = 0; k < i; k++) { 343827bd09bSSatish Balay off = *iptr++; 344827bd09bSSatish Balay len = *iptr++; 3459566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 346792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1)); 347827bd09bSSatish Balay x_ptr += len; 348827bd09bSSatish Balay } 349827bd09bSSatish Balay 350827bd09bSSatish Balay /* compute v_l = v_l - z */ 351ca8e9878SJed Brown PCTFS_rvec_zero(v + a_n, a_m - a_n); 3529566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 353792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1)); 354827bd09bSSatish Balay 355827bd09bSSatish Balay /* compute u_l = A.v_l */ 3562fa5cd67SKarl Rupp if (a_n != a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim); 357ca8e9878SJed Brown PCTFS_rvec_zero(u, n); 358827bd09bSSatish Balay do_matvec(xxt_handle->mvi, v, u); 359827bd09bSSatish Balay 360827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */ 3619566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 362792fecdfSBarry Smith PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, v, &i1)); 363827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */ 364b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 365827bd09bSSatish Balay 3668f1a2a5eSBarry Smith alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha); 367827bd09bSSatish Balay 368827bd09bSSatish Balay /* check for small alpha */ 369827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 37063a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha)); 371827bd09bSSatish Balay 372827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 373ca8e9878SJed Brown PCTFS_rvec_scale(v, 1.0 / alpha, n); 374827bd09bSSatish Balay 375827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 376827bd09bSSatish Balay flag = 1; 377827bd09bSSatish Balay off = len = 0; 378db4deed7SKarl Rupp for (k = 0; k < n; k++) { 379db4deed7SKarl Rupp if (v[k] != 0.0) { 380827bd09bSSatish Balay len = k; 3819371c9d4SSatish Balay if (flag) { 3829371c9d4SSatish Balay off = k; 3839371c9d4SSatish Balay flag = 0; 3849371c9d4SSatish Balay } 385827bd09bSSatish Balay } 386827bd09bSSatish Balay } 387827bd09bSSatish Balay 388827bd09bSSatish Balay len -= (off - 1); 389827bd09bSSatish Balay 390db4deed7SKarl Rupp if (len > 0) { 391db4deed7SKarl Rupp if ((xxt_nnz + len) > xxt_max_nnz) { 3929566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "increasing space for X by 2x!\n")); 393827bd09bSSatish Balay xxt_max_nnz *= 2; 394a501084fSBarry Smith x_ptr = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar)); 395ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr, x, xxt_nnz); 396a501084fSBarry Smith free(x); 397827bd09bSSatish Balay x = x_ptr; 398827bd09bSSatish Balay x_ptr += xxt_nnz; 399827bd09bSSatish Balay } 400827bd09bSSatish Balay xxt_nnz += len; 401ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr, v + off, len); 402827bd09bSSatish Balay 403827bd09bSSatish Balay col_indices[2 * i] = off; 404827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = len; 405827bd09bSSatish Balay col_vals[i] = x_ptr; 4069371c9d4SSatish Balay } else { 407827bd09bSSatish Balay col_indices[2 * i] = 0; 408827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = 0; 409827bd09bSSatish Balay col_vals[i] = x_ptr; 410827bd09bSSatish Balay } 411827bd09bSSatish Balay } 412827bd09bSSatish Balay 413827bd09bSSatish Balay /* close off stages for execution phase */ 414db4deed7SKarl Rupp while (dim != level) { 415827bd09bSSatish Balay stages[dim++] = i; 41663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level)); 417827bd09bSSatish Balay } 418827bd09bSSatish Balay stages[dim] = i; 419827bd09bSSatish Balay 420827bd09bSSatish Balay xxt_handle->info->n = xxt_handle->mvi->n; 421827bd09bSSatish Balay xxt_handle->info->m = m; 422827bd09bSSatish Balay xxt_handle->info->nnz = xxt_nnz; 423827bd09bSSatish Balay xxt_handle->info->max_nnz = xxt_max_nnz; 424827bd09bSSatish Balay xxt_handle->info->msg_buf_sz = stages[level] - stages[0]; 425a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 426a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 427827bd09bSSatish Balay xxt_handle->info->x = x; 428827bd09bSSatish Balay xxt_handle->info->col_vals = col_vals; 429827bd09bSSatish Balay xxt_handle->info->col_sz = col_sz; 430827bd09bSSatish Balay xxt_handle->info->col_indices = col_indices; 431827bd09bSSatish Balay xxt_handle->info->stages = stages; 432827bd09bSSatish Balay xxt_handle->info->nsolves = 0; 433827bd09bSSatish Balay xxt_handle->info->tot_solve_time = 0.0; 434827bd09bSSatish Balay 435a501084fSBarry Smith free(segs); 436a501084fSBarry Smith free(u); 437a501084fSBarry Smith free(v); 438a501084fSBarry Smith free(uu); 439a501084fSBarry Smith free(z); 440a501084fSBarry Smith free(w); 441827bd09bSSatish Balay 442827bd09bSSatish Balay return (0); 443827bd09bSSatish Balay } 444827bd09bSSatish Balay 4459371c9d4SSatish Balay static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc) { 44652f87cdaSBarry Smith PetscInt off, len, *iptr; 44752f87cdaSBarry Smith PetscInt level = xxt_handle->level; 44852f87cdaSBarry Smith PetscInt n = xxt_handle->info->n; 44952f87cdaSBarry Smith PetscInt m = xxt_handle->info->m; 45052f87cdaSBarry Smith PetscInt *stages = xxt_handle->info->stages; 45152f87cdaSBarry Smith PetscInt *col_indices = xxt_handle->info->col_indices; 452a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr; 453a501084fSBarry Smith PetscScalar *solve_uu = xxt_handle->info->solve_uu; 454a501084fSBarry Smith PetscScalar *solve_w = xxt_handle->info->solve_w; 455a501084fSBarry Smith PetscScalar *x = xxt_handle->info->x; 45671a0148aSBarry Smith PetscBLASInt i1 = 1, dlen; 457827bd09bSSatish Balay 4580924e98cSBarry Smith PetscFunctionBegin; 459827bd09bSSatish Balay uu_ptr = solve_uu; 460ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr, m); 461827bd09bSSatish Balay 462827bd09bSSatish Balay /* x = X.Y^T.b */ 463827bd09bSSatish Balay /* uu = Y^T.b */ 464db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) { 465c5df96a5SBarry Smith off = *iptr++; 466c5df96a5SBarry Smith len = *iptr++; 4679566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 468792fecdfSBarry Smith PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, x_ptr, &i1)); 469827bd09bSSatish Balay } 470827bd09bSSatish Balay 471827bd09bSSatish Balay /* comunication of beta */ 472827bd09bSSatish Balay uu_ptr = solve_uu; 4739566063dSJacob Faibussowitsch if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages)); 474827bd09bSSatish Balay 475ca8e9878SJed Brown PCTFS_rvec_zero(uc, n); 476827bd09bSSatish Balay 477827bd09bSSatish Balay /* x = X.uu */ 478db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) { 479c5df96a5SBarry Smith off = *iptr++; 480c5df96a5SBarry Smith len = *iptr++; 4819566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 482792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1)); 483827bd09bSSatish Balay } 4840924e98cSBarry Smith PetscFunctionReturn(0); 485827bd09bSSatish Balay } 486827bd09bSSatish Balay 4879371c9d4SSatish Balay static PetscErrorCode check_handle(xxt_ADT xxt_handle) { 48852f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX}; 489827bd09bSSatish Balay 4900924e98cSBarry Smith PetscFunctionBegin; 49163a3b9bcSJacob Faibussowitsch PetscCheck(xxt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", xxt_handle); 492827bd09bSSatish Balay 493827bd09bSSatish Balay vals[0] = vals[1] = xxt_handle->id; 494dd39110bSPierre Jolivet PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op); 49563a3b9bcSJacob 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); 4960924e98cSBarry Smith PetscFunctionReturn(0); 497827bd09bSSatish Balay } 498827bd09bSSatish Balay 4999371c9d4SSatish Balay static PetscErrorCode det_separators(xxt_ADT xxt_handle) { 50052f87cdaSBarry Smith PetscInt i, ct, id; 50152f87cdaSBarry Smith PetscInt mask, edge, *iptr; 50252f87cdaSBarry Smith PetscInt *dir, *used; 50352f87cdaSBarry Smith PetscInt sum[4], w[4]; 504a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 50552f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 506a501084fSBarry Smith PetscScalar *lhs, *rhs; 50752f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo = 0; 508ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 50952f87cdaSBarry Smith PetscInt *local2global = xxt_handle->mvi->local2global; 51052f87cdaSBarry Smith PetscInt n = xxt_handle->mvi->n; 51152f87cdaSBarry Smith PetscInt m = xxt_handle->mvi->m; 51252f87cdaSBarry Smith PetscInt level = xxt_handle->level; 513ab824b78SBarry Smith PetscInt shared = 0; 514827bd09bSSatish Balay 5150924e98cSBarry Smith PetscFunctionBegin; 51652f87cdaSBarry Smith dir = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 51752f87cdaSBarry Smith nsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 51852f87cdaSBarry Smith lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 51952f87cdaSBarry Smith fo = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1)); 52052f87cdaSBarry Smith used = (PetscInt *)malloc(sizeof(PetscInt) * n); 521827bd09bSSatish Balay 522ca8e9878SJed Brown PCTFS_ivec_zero(dir, level + 1); 523ca8e9878SJed Brown PCTFS_ivec_zero(nsep, level + 1); 524ca8e9878SJed Brown PCTFS_ivec_zero(lnsep, level + 1); 525ca8e9878SJed Brown PCTFS_ivec_set(fo, -1, n + 1); 526ca8e9878SJed Brown PCTFS_ivec_zero(used, n); 527827bd09bSSatish Balay 5288cda6cd7SBarry Smith lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 5298cda6cd7SBarry Smith rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 530827bd09bSSatish Balay 531827bd09bSSatish Balay /* determine the # of unique dof */ 532ca8e9878SJed Brown PCTFS_rvec_zero(lhs, m); 533ca8e9878SJed Brown PCTFS_rvec_set(lhs, 1.0, n); 534ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level); 535ca8e9878SJed Brown PCTFS_rvec_zero(rsum, 2); 5363d3eaba7SBarry Smith for (i = 0; i < n; i++) { 5372fa5cd67SKarl Rupp if (lhs[i] != 0.0) { 5389371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 5399371c9d4SSatish Balay rsum[1] += lhs[i]; 5402fa5cd67SKarl Rupp } 541827bd09bSSatish Balay } 542b1c944f5SJed Brown PCTFS_grop_hc(rsum, rw, 2, op, level); 543827bd09bSSatish Balay rsum[0] += 0.1; 544827bd09bSSatish Balay rsum[1] += 0.1; 545827bd09bSSatish Balay 54677b4d14cSPeter Brune if (PetscAbsScalar(rsum[0] - rsum[1]) > EPS) shared = 1; 547827bd09bSSatish Balay 54852f87cdaSBarry Smith xxt_handle->info->n_global = xxt_handle->info->m_global = (PetscInt)rsum[0]; 54952f87cdaSBarry Smith xxt_handle->mvi->n_global = xxt_handle->mvi->m_global = (PetscInt)rsum[0]; 550827bd09bSSatish Balay 551827bd09bSSatish Balay /* determine separator sets top down */ 5522fa5cd67SKarl Rupp if (shared) { 553db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 554827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 555ca8e9878SJed Brown (id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m); 556ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge); 557827bd09bSSatish Balay 558827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 559ca8e9878SJed Brown (id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m); 560ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge); 561827bd09bSSatish Balay 562db4deed7SKarl Rupp for (i = 0; i < n; i++) { 563db4deed7SKarl Rupp if (id < mask) { 5642fa5cd67SKarl Rupp if (lhs[i] != 0.0) lhs[i] = 1.0; 565827bd09bSSatish Balay } 566db4deed7SKarl Rupp if (id >= mask) { 5672fa5cd67SKarl Rupp if (rhs[i] != 0.0) rhs[i] = 1.0; 568827bd09bSSatish Balay } 569827bd09bSSatish Balay } 570827bd09bSSatish Balay 5712fa5cd67SKarl Rupp if (id < mask) PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge - 1); 5722fa5cd67SKarl Rupp else PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge - 1); 573827bd09bSSatish Balay 574827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 575ca8e9878SJed Brown PCTFS_rvec_zero(rsum, 4); 576db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum, 4), ct = i = 0; i < n; i++) { 577db4deed7SKarl Rupp if (!used[i]) { 578827bd09bSSatish Balay /* number of unmarked dofs on node */ 579827bd09bSSatish Balay ct++; 580827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 581db4deed7SKarl Rupp if (id < mask) { 5829371c9d4SSatish Balay if (lhs[i] != 0.0) { 5839371c9d4SSatish Balay sum[0]++; 5849371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 5859371c9d4SSatish Balay } 586827bd09bSSatish Balay } 587827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 588db4deed7SKarl Rupp if (id >= mask) { 5899371c9d4SSatish Balay if (rhs[i] != 0.0) { 5909371c9d4SSatish Balay sum[1]++; 5919371c9d4SSatish Balay rsum[1] += 1.0 / rhs[i]; 5929371c9d4SSatish Balay } 593827bd09bSSatish Balay } 594827bd09bSSatish Balay } 595827bd09bSSatish Balay } 596827bd09bSSatish Balay 597827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 598827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 599827bd09bSSatish Balay (id < mask) ? (rsum[2] = ct) : (rsum[3] = ct); 600b1c944f5SJed Brown PCTFS_giop_hc(sum, w, 4, op, edge); 601b1c944f5SJed Brown PCTFS_grop_hc(rsum, rw, 4, op, edge); 6029371c9d4SSatish Balay rsum[0] += 0.1; 6039371c9d4SSatish Balay rsum[1] += 0.1; 6049371c9d4SSatish Balay rsum[2] += 0.1; 6059371c9d4SSatish Balay rsum[3] += 0.1; 606827bd09bSSatish Balay 607db4deed7SKarl Rupp if (id < mask) { 608827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 609db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 610db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 6119371c9d4SSatish Balay ct++; 6129371c9d4SSatish Balay nfo++; 613827bd09bSSatish Balay 61408401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n"); 615827bd09bSSatish Balay 616827bd09bSSatish Balay *--iptr = local2global[i]; 617827bd09bSSatish Balay used[i] = edge; 618827bd09bSSatish Balay } 619827bd09bSSatish Balay } 6202fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 621827bd09bSSatish Balay 622827bd09bSSatish Balay lnsep[edge] = ct; 62352f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[0]; 624827bd09bSSatish Balay dir[edge] = LEFT; 625827bd09bSSatish Balay } 626827bd09bSSatish Balay 627db4deed7SKarl Rupp if (id >= mask) { 628827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 629db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 630db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 6319371c9d4SSatish Balay ct++; 6329371c9d4SSatish Balay nfo++; 633827bd09bSSatish Balay 63408401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n"); 635827bd09bSSatish Balay 636827bd09bSSatish Balay *--iptr = local2global[i]; 637827bd09bSSatish Balay used[i] = edge; 638827bd09bSSatish Balay } 639827bd09bSSatish Balay } 6402fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 641827bd09bSSatish Balay 642827bd09bSSatish Balay lnsep[edge] = ct; 64352f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[1]; 644827bd09bSSatish Balay dir[edge] = RIGHT; 645827bd09bSSatish Balay } 646827bd09bSSatish Balay 647827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 648827bd09bSSatish Balay /* do we need full set of separators for this? */ 649827bd09bSSatish Balay 650827bd09bSSatish Balay /* fold rhs hc into lower */ 6512fa5cd67SKarl Rupp if (id >= mask) id -= mask; 652827bd09bSSatish Balay } 653db4deed7SKarl Rupp } else { 654db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 655827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 656ca8e9878SJed Brown (id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m); 657ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge); 658827bd09bSSatish Balay 659827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 660ca8e9878SJed Brown (id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m); 661ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge); 662827bd09bSSatish Balay 663827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 664db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum, 4), ct = i = 0; i < n; i++) { 665db4deed7SKarl Rupp if (!used[i]) { 666827bd09bSSatish Balay /* number of unmarked dofs on node */ 667827bd09bSSatish Balay ct++; 668827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6692fa5cd67SKarl Rupp if ((id < mask) && (lhs[i] != 0.0)) sum[0]++; 670827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6712fa5cd67SKarl Rupp if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++; 672827bd09bSSatish Balay } 673827bd09bSSatish Balay } 674827bd09bSSatish Balay 675827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 676827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 677b1c944f5SJed Brown PCTFS_giop_hc(sum, w, 4, op, edge); 678827bd09bSSatish Balay 679827bd09bSSatish Balay /* lhs hc wins */ 680db4deed7SKarl Rupp if (sum[2] >= sum[3]) { 681db4deed7SKarl Rupp if (id < mask) { 682827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 683db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 684db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 6859371c9d4SSatish Balay ct++; 6869371c9d4SSatish Balay nfo++; 687827bd09bSSatish Balay *--iptr = local2global[i]; 688827bd09bSSatish Balay used[i] = edge; 689827bd09bSSatish Balay } 690827bd09bSSatish Balay } 6912fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 692827bd09bSSatish Balay lnsep[edge] = ct; 693827bd09bSSatish Balay } 694827bd09bSSatish Balay nsep[edge] = sum[0]; 695827bd09bSSatish Balay dir[edge] = LEFT; 696db4deed7SKarl Rupp } else { /* rhs hc wins */ 697db4deed7SKarl Rupp if (id >= mask) { 698827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 699db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 700db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 7019371c9d4SSatish Balay ct++; 7029371c9d4SSatish Balay nfo++; 703827bd09bSSatish Balay *--iptr = local2global[i]; 704827bd09bSSatish Balay used[i] = edge; 705827bd09bSSatish Balay } 706827bd09bSSatish Balay } 7072fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 708827bd09bSSatish Balay lnsep[edge] = ct; 709827bd09bSSatish Balay } 710827bd09bSSatish Balay nsep[edge] = sum[1]; 711827bd09bSSatish Balay dir[edge] = RIGHT; 712827bd09bSSatish Balay } 713827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 714827bd09bSSatish Balay /* do we need full set of separators for this? */ 715827bd09bSSatish Balay 716827bd09bSSatish Balay /* fold rhs hc into lower */ 7172fa5cd67SKarl Rupp if (id >= mask) id -= mask; 718827bd09bSSatish Balay } 719827bd09bSSatish Balay } 720827bd09bSSatish Balay 721827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 722db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 723db4deed7SKarl Rupp if (!used[i]) { 7249371c9d4SSatish Balay ct++; 7259371c9d4SSatish Balay nfo++; 726827bd09bSSatish Balay *--iptr = local2global[i]; 727827bd09bSSatish Balay used[i] = edge; 728827bd09bSSatish Balay } 729827bd09bSSatish Balay } 7302fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 731827bd09bSSatish Balay lnsep[edge] = ct; 732827bd09bSSatish Balay nsep[edge] = ct; 733827bd09bSSatish Balay dir[edge] = LEFT; 734827bd09bSSatish Balay 735827bd09bSSatish Balay xxt_handle->info->nsep = nsep; 736827bd09bSSatish Balay xxt_handle->info->lnsep = lnsep; 737827bd09bSSatish Balay xxt_handle->info->fo = fo; 738827bd09bSSatish Balay xxt_handle->info->nfo = nfo; 739827bd09bSSatish Balay 740a501084fSBarry Smith free(dir); 741a501084fSBarry Smith free(lhs); 742a501084fSBarry Smith free(rhs); 743a501084fSBarry Smith free(used); 7440924e98cSBarry Smith PetscFunctionReturn(0); 745827bd09bSSatish Balay } 746827bd09bSSatish Balay 7479371c9d4SSatish Balay static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data) { 748827bd09bSSatish Balay mv_info *mvi; 749827bd09bSSatish Balay 750a501084fSBarry Smith mvi = (mv_info *)malloc(sizeof(mv_info)); 751827bd09bSSatish Balay mvi->n = n; 752827bd09bSSatish Balay mvi->m = m; 753827bd09bSSatish Balay mvi->n_global = -1; 754827bd09bSSatish Balay mvi->m_global = -1; 75552f87cdaSBarry Smith mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt)); 756ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global, local2global, m); 757827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7585c8f6a95SKarl Rupp mvi->matvec = matvec; 759827bd09bSSatish Balay mvi->grid_data = grid_data; 760827bd09bSSatish Balay 761827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */ 762ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 763827bd09bSSatish Balay 764827bd09bSSatish Balay return (mvi); 765827bd09bSSatish Balay } 766827bd09bSSatish Balay 7679371c9d4SSatish Balay static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) { 7680924e98cSBarry Smith PetscFunctionBegin; 769827bd09bSSatish Balay A->matvec((mv_info *)A->grid_data, v, u); 7700924e98cSBarry Smith PetscFunctionReturn(0); 771827bd09bSSatish Balay } 772