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 { 84*3ba16761SJacob Faibussowitsch PetscFunctionBegin; 85*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 86*3ba16761SJacob 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 */ 102*3ba16761SJacob Faibussowitsch PetscCall(det_separators(xxt_handle)); 103827bd09bSSatish Balay 104*3ba16761SJacob Faibussowitsch PetscCall(do_xxt_factor(xxt_handle)); 105*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106827bd09bSSatish Balay } 107827bd09bSSatish Balay 108d71ae5a4SJacob Faibussowitsch PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b) 109d71ae5a4SJacob Faibussowitsch { 110*3ba16761SJacob Faibussowitsch PetscFunctionBegin; 111*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 112*3ba16761SJacob Faibussowitsch PetscCall(check_handle(xxt_handle)); 113827bd09bSSatish Balay 114827bd09bSSatish Balay /* need to copy b into x? */ 115*3ba16761SJacob Faibussowitsch if (b) PetscCall(PCTFS_rvec_copy(x, b, xxt_handle->mvi->n)); 116*3ba16761SJacob Faibussowitsch PetscCall(do_xxt_solve(xxt_handle, x)); 117*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118827bd09bSSatish Balay } 119827bd09bSSatish Balay 120*3ba16761SJacob Faibussowitsch PetscErrorCode XXT_free(xxt_ADT xxt_handle) 121d71ae5a4SJacob Faibussowitsch { 122*3ba16761SJacob Faibussowitsch PetscFunctionBegin; 123*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 124*3ba16761SJacob 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); 139*3ba16761SJacob 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 */ 146*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147827bd09bSSatish Balay } 148827bd09bSSatish Balay 14911cc89d2SBarry Smith /* This function is currently unused */ 150d71ae5a4SJacob Faibussowitsch PetscErrorCode XXT_stats(xxt_ADT xxt_handle) 151d71ae5a4SJacob Faibussowitsch { 15252f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD}; 15352f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD}; 15452f87cdaSBarry Smith PetscInt vals[9], work[9]; 155a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 156827bd09bSSatish Balay 15763a3b9bcSJacob Faibussowitsch PetscFunctionBegin; 15863a3b9bcSJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 15963a3b9bcSJacob Faibussowitsch PetscCall(check_handle(xxt_handle)); 160827bd09bSSatish Balay 161827bd09bSSatish Balay /* if factorization not done there are no stats */ 162db4deed7SKarl Rupp if (!xxt_handle->info || !xxt_handle->mvi) { 1639566063dSJacob Faibussowitsch if (!PCTFS_my_id) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "XXT_stats() :: no stats available!\n")); 164*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165827bd09bSSatish Balay } 166827bd09bSSatish Balay 167827bd09bSSatish Balay vals[0] = vals[1] = vals[2] = xxt_handle->info->nnz; 168827bd09bSSatish Balay vals[3] = vals[4] = vals[5] = xxt_handle->mvi->n; 169827bd09bSSatish Balay vals[6] = vals[7] = vals[8] = xxt_handle->info->msg_buf_sz; 170*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op)); 171827bd09bSSatish Balay 1722fa5cd67SKarl Rupp fvals[0] = fvals[1] = fvals[2] = xxt_handle->info->tot_solve_time / xxt_handle->info->nsolves++; 173*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop(fvals, fwork, PETSC_STATIC_ARRAY_LENGTH(fop) - 1, fop)); 174827bd09bSSatish Balay 175db4deed7SKarl Rupp if (!PCTFS_my_id) { 17663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[0])); 17763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[1])); 17863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_nnz=%g\n", PCTFS_my_id, (double)(1.0 * vals[2] / PCTFS_num_nodes))); 17963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[2])); 18063a3b9bcSJacob 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))))); 18163a3b9bcSJacob 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))))); 18263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[3])); 18363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[4])); 18463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_n =%g\n", PCTFS_my_id, (double)(1.0 * vals[5] / PCTFS_num_nodes))); 18563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot xxt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[5])); 18663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[6])); 18763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[7])); 18863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_buf=%g\n", PCTFS_my_id, (double)(1.0 * vals[8] / PCTFS_num_nodes))); 18963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[0]))); 19063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[1]))); 19163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[2] / PCTFS_num_nodes))); 192827bd09bSSatish Balay } 193*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 194827bd09bSSatish Balay } 195827bd09bSSatish Balay 19611cc89d2SBarry Smith /* 197827bd09bSSatish Balay 198827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 199827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 200827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 201827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 202827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 203827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 204ca8e9878SJed Brown PCTFS_gs_init/gop). 205827bd09bSSatish Balay 206827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 207827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 20811cc89d2SBarry Smith */ 209d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle) 210d71ae5a4SJacob Faibussowitsch { 2117b1ae94cSBarry Smith return xxt_generate(xxt_handle); 212827bd09bSSatish Balay } 213827bd09bSSatish Balay 214d71ae5a4SJacob Faibussowitsch static PetscErrorCode xxt_generate(xxt_ADT xxt_handle) 215d71ae5a4SJacob Faibussowitsch { 21652f87cdaSBarry Smith PetscInt i, j, k, idex; 21752f87cdaSBarry Smith PetscInt dim, col; 218a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 21952f87cdaSBarry Smith PetscInt *segs; 22052f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 22152f87cdaSBarry Smith PetscInt off, len; 222a501084fSBarry Smith PetscScalar *x_ptr; 22352f87cdaSBarry Smith PetscInt *iptr, flag; 22452f87cdaSBarry Smith PetscInt start = 0, end, work; 22552f87cdaSBarry Smith PetscInt op2[] = {GL_MIN, 0}; 226ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 22752f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 22852f87cdaSBarry Smith PetscInt a_n = xxt_handle->mvi->n; 22952f87cdaSBarry Smith PetscInt a_m = xxt_handle->mvi->m; 23052f87cdaSBarry Smith PetscInt *a_local2global = xxt_handle->mvi->local2global; 23152f87cdaSBarry Smith PetscInt level; 23252f87cdaSBarry Smith PetscInt xxt_nnz = 0, xxt_max_nnz = 0; 23352f87cdaSBarry Smith PetscInt n, m; 23452f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages; 235a501084fSBarry Smith PetscScalar **col_vals, *x; 23652f87cdaSBarry Smith PetscInt n_global; 23771a0148aSBarry Smith PetscBLASInt i1 = 1, dlen; 238a501084fSBarry Smith PetscScalar dm1 = -1.0; 239827bd09bSSatish Balay 240827bd09bSSatish Balay n = xxt_handle->mvi->n; 241827bd09bSSatish Balay nsep = xxt_handle->info->nsep; 242827bd09bSSatish Balay lnsep = xxt_handle->info->lnsep; 243827bd09bSSatish Balay fo = xxt_handle->info->fo; 244827bd09bSSatish Balay end = lnsep[0]; 245827bd09bSSatish Balay level = xxt_handle->level; 246ca8e9878SJed Brown PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 247827bd09bSSatish Balay 248827bd09bSSatish Balay /* is there a null space? */ 249827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2502fa5cd67SKarl Rupp for (i = 0, j = 0; i <= level; i++) j += nsep[i]; 251827bd09bSSatish Balay 252827bd09bSSatish Balay m = j - xxt_handle->ns; 25348a46eb9SPierre 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)); 254827bd09bSSatish Balay 255827bd09bSSatish Balay /* get and initialize storage for x local */ 256827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 25752f87cdaSBarry Smith col_sz = (PetscInt *)malloc(m * sizeof(PetscInt)); 25852f87cdaSBarry Smith col_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt)); 259a501084fSBarry Smith col_vals = (PetscScalar **)malloc(m * sizeof(PetscScalar *)); 260db4deed7SKarl Rupp for (i = j = 0; i < m; i++, j += 2) { 261827bd09bSSatish Balay col_indices[j] = col_indices[j + 1] = col_sz[i] = -1; 262827bd09bSSatish Balay col_vals[i] = NULL; 263827bd09bSSatish Balay } 264827bd09bSSatish Balay col_indices[j] = -1; 265827bd09bSSatish Balay 266827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 267827bd09bSSatish Balay /* this looks like nsep[]=segments */ 26852f87cdaSBarry Smith stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 26952f87cdaSBarry Smith segs = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 270*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(stages, level + 1)); 271ca8e9878SJed Brown PCTFS_ivec_copy(segs, nsep, level + 1); 2722fa5cd67SKarl Rupp for (i = 0; i < level; i++) segs[i + 1] += segs[i]; 273827bd09bSSatish Balay stages[0] = segs[0]; 274827bd09bSSatish Balay 275827bd09bSSatish Balay /* temporary vectors */ 276a501084fSBarry Smith u = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 277a501084fSBarry Smith z = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 278a501084fSBarry Smith v = (PetscScalar *)malloc(a_m * sizeof(PetscScalar)); 279a501084fSBarry Smith uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 280a501084fSBarry Smith w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 281827bd09bSSatish Balay 282827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 2832fa5cd67SKarl Rupp for (i = 1, j = 0; i <= level; i++) j += nsep[i]; 284827bd09bSSatish Balay 285827bd09bSSatish Balay /* storage for sparse x values */ 286827bd09bSSatish Balay n_global = xxt_handle->info->n_global; 28785ec1a3cSBarry Smith xxt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes; 288a501084fSBarry Smith x = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar)); 289827bd09bSSatish Balay xxt_nnz = 0; 290827bd09bSSatish Balay 291827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 292827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 2932fa5cd67SKarl Rupp for (dim = i = j = 0; i < m; i++) { 294827bd09bSSatish Balay /* time to move to the next level? */ 295d4af0d30SBarry Smith while (i == segs[dim]) { 29608401ef6SPierre Jolivet PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level"); 297827bd09bSSatish Balay stages[dim++] = i; 298827bd09bSSatish Balay end += lnsep[dim]; 299827bd09bSSatish Balay } 300827bd09bSSatish Balay stages[dim] = i; 301827bd09bSSatish Balay 302827bd09bSSatish Balay /* which column are we firing? */ 303827bd09bSSatish Balay /* i.e. set v_l */ 304827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 305827bd09bSSatish Balay (start < end) ? (col = fo[start]) : (col = INT_MAX); 306*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(&col, &work, 1, op2, dim)); 307827bd09bSSatish Balay 308827bd09bSSatish Balay /* shouldn't need this */ 309db4deed7SKarl Rupp if (col == INT_MAX) { 3109566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n")); 311827bd09bSSatish Balay continue; 312827bd09bSSatish Balay } 313827bd09bSSatish Balay 314827bd09bSSatish Balay /* do I own it? I should */ 315*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(v, a_m)); 316db4deed7SKarl Rupp if (col == fo[start]) { 317827bd09bSSatish Balay start++; 318ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_n); 3190fdf79fbSJacob Faibussowitsch PetscCheck(idex != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!"); 3209371c9d4SSatish Balay v[idex] = 1.0; 3219371c9d4SSatish Balay j++; 322db4deed7SKarl Rupp } else { 323ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_m); 3242fa5cd67SKarl Rupp if (idex != -1) v[idex] = 1.0; 325827bd09bSSatish Balay } 326827bd09bSSatish Balay 327827bd09bSSatish Balay /* perform u = A.v_l */ 328*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(u, n)); 329*3ba16761SJacob Faibussowitsch PetscCall(do_matvec(xxt_handle->mvi, v, u)); 330827bd09bSSatish Balay 331827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 332827bd09bSSatish Balay /* technically only need to zero out first i entries */ 333827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */ 334*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uu, m)); 335827bd09bSSatish Balay x_ptr = x; 336827bd09bSSatish Balay iptr = col_indices; 337db4deed7SKarl Rupp for (k = 0; k < i; k++) { 338827bd09bSSatish Balay off = *iptr++; 339827bd09bSSatish Balay len = *iptr++; 3409566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 341792fecdfSBarry Smith PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, x_ptr, &i1)); 342827bd09bSSatish Balay x_ptr += len; 343827bd09bSSatish Balay } 344827bd09bSSatish Balay 345827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 3469566063dSJacob Faibussowitsch PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages)); 347827bd09bSSatish Balay 348827bd09bSSatish Balay /* z = X.uu */ 349*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(z, n)); 350827bd09bSSatish Balay x_ptr = x; 351827bd09bSSatish Balay iptr = col_indices; 352db4deed7SKarl Rupp for (k = 0; k < i; k++) { 353827bd09bSSatish Balay off = *iptr++; 354827bd09bSSatish Balay len = *iptr++; 3559566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 356792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1)); 357827bd09bSSatish Balay x_ptr += len; 358827bd09bSSatish Balay } 359827bd09bSSatish Balay 360827bd09bSSatish Balay /* compute v_l = v_l - z */ 361*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(v + a_n, a_m - a_n)); 3629566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 363792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1)); 364827bd09bSSatish Balay 365827bd09bSSatish Balay /* compute u_l = A.v_l */ 366*3ba16761SJacob Faibussowitsch if (a_n != a_m) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim)); 367*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(u, n)); 368*3ba16761SJacob Faibussowitsch PetscCall(do_matvec(xxt_handle->mvi, v, u)); 369827bd09bSSatish Balay 370827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */ 3719566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 372792fecdfSBarry Smith PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, v, &i1)); 373827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */ 374*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim)); 375827bd09bSSatish Balay 3768f1a2a5eSBarry Smith alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha); 377827bd09bSSatish Balay 378827bd09bSSatish Balay /* check for small alpha */ 379827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 38063a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha)); 381827bd09bSSatish Balay 382827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 383*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_scale(v, 1.0 / alpha, n)); 384827bd09bSSatish Balay 385827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 386827bd09bSSatish Balay flag = 1; 387827bd09bSSatish Balay off = len = 0; 388db4deed7SKarl Rupp for (k = 0; k < n; k++) { 389db4deed7SKarl Rupp if (v[k] != 0.0) { 390827bd09bSSatish Balay len = k; 3919371c9d4SSatish Balay if (flag) { 3929371c9d4SSatish Balay off = k; 3939371c9d4SSatish Balay flag = 0; 3949371c9d4SSatish Balay } 395827bd09bSSatish Balay } 396827bd09bSSatish Balay } 397827bd09bSSatish Balay 398827bd09bSSatish Balay len -= (off - 1); 399827bd09bSSatish Balay 400db4deed7SKarl Rupp if (len > 0) { 401db4deed7SKarl Rupp if ((xxt_nnz + len) > xxt_max_nnz) { 4029566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "increasing space for X by 2x!\n")); 403827bd09bSSatish Balay xxt_max_nnz *= 2; 404a501084fSBarry Smith x_ptr = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar)); 405*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_copy(x_ptr, x, xxt_nnz)); 406a501084fSBarry Smith free(x); 407827bd09bSSatish Balay x = x_ptr; 408827bd09bSSatish Balay x_ptr += xxt_nnz; 409827bd09bSSatish Balay } 410827bd09bSSatish Balay xxt_nnz += len; 411*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_copy(x_ptr, v + off, len)); 412827bd09bSSatish Balay 413827bd09bSSatish Balay col_indices[2 * i] = off; 414827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = len; 415827bd09bSSatish Balay col_vals[i] = x_ptr; 4169371c9d4SSatish Balay } else { 417827bd09bSSatish Balay col_indices[2 * i] = 0; 418827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = 0; 419827bd09bSSatish Balay col_vals[i] = x_ptr; 420827bd09bSSatish Balay } 421827bd09bSSatish Balay } 422827bd09bSSatish Balay 423827bd09bSSatish Balay /* close off stages for execution phase */ 424db4deed7SKarl Rupp while (dim != level) { 425827bd09bSSatish Balay stages[dim++] = i; 42663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level)); 427827bd09bSSatish Balay } 428827bd09bSSatish Balay stages[dim] = i; 429827bd09bSSatish Balay 430827bd09bSSatish Balay xxt_handle->info->n = xxt_handle->mvi->n; 431827bd09bSSatish Balay xxt_handle->info->m = m; 432827bd09bSSatish Balay xxt_handle->info->nnz = xxt_nnz; 433827bd09bSSatish Balay xxt_handle->info->max_nnz = xxt_max_nnz; 434827bd09bSSatish Balay xxt_handle->info->msg_buf_sz = stages[level] - stages[0]; 435a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 436a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 437827bd09bSSatish Balay xxt_handle->info->x = x; 438827bd09bSSatish Balay xxt_handle->info->col_vals = col_vals; 439827bd09bSSatish Balay xxt_handle->info->col_sz = col_sz; 440827bd09bSSatish Balay xxt_handle->info->col_indices = col_indices; 441827bd09bSSatish Balay xxt_handle->info->stages = stages; 442827bd09bSSatish Balay xxt_handle->info->nsolves = 0; 443827bd09bSSatish Balay xxt_handle->info->tot_solve_time = 0.0; 444827bd09bSSatish Balay 445a501084fSBarry Smith free(segs); 446a501084fSBarry Smith free(u); 447a501084fSBarry Smith free(v); 448a501084fSBarry Smith free(uu); 449a501084fSBarry Smith free(z); 450a501084fSBarry Smith free(w); 451827bd09bSSatish Balay 452*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 453827bd09bSSatish Balay } 454827bd09bSSatish Balay 455d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc) 456d71ae5a4SJacob Faibussowitsch { 45752f87cdaSBarry Smith PetscInt off, len, *iptr; 45852f87cdaSBarry Smith PetscInt level = xxt_handle->level; 45952f87cdaSBarry Smith PetscInt n = xxt_handle->info->n; 46052f87cdaSBarry Smith PetscInt m = xxt_handle->info->m; 46152f87cdaSBarry Smith PetscInt *stages = xxt_handle->info->stages; 46252f87cdaSBarry Smith PetscInt *col_indices = xxt_handle->info->col_indices; 463a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr; 464a501084fSBarry Smith PetscScalar *solve_uu = xxt_handle->info->solve_uu; 465a501084fSBarry Smith PetscScalar *solve_w = xxt_handle->info->solve_w; 466a501084fSBarry Smith PetscScalar *x = xxt_handle->info->x; 46771a0148aSBarry Smith PetscBLASInt i1 = 1, dlen; 468827bd09bSSatish Balay 4690924e98cSBarry Smith PetscFunctionBegin; 470827bd09bSSatish Balay uu_ptr = solve_uu; 471*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uu_ptr, m)); 472827bd09bSSatish Balay 473827bd09bSSatish Balay /* x = X.Y^T.b */ 474827bd09bSSatish Balay /* uu = Y^T.b */ 475db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) { 476c5df96a5SBarry Smith off = *iptr++; 477c5df96a5SBarry Smith len = *iptr++; 4789566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 479792fecdfSBarry Smith PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, x_ptr, &i1)); 480827bd09bSSatish Balay } 481827bd09bSSatish Balay 482d5b43468SJose E. Roman /* communication of beta */ 483827bd09bSSatish Balay uu_ptr = solve_uu; 4849566063dSJacob Faibussowitsch if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages)); 485827bd09bSSatish Balay 486*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uc, n)); 487827bd09bSSatish Balay 488827bd09bSSatish Balay /* x = X.uu */ 489db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) { 490c5df96a5SBarry Smith off = *iptr++; 491c5df96a5SBarry Smith len = *iptr++; 4929566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 493792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1)); 494827bd09bSSatish Balay } 495*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 496827bd09bSSatish Balay } 497827bd09bSSatish Balay 498d71ae5a4SJacob Faibussowitsch static PetscErrorCode check_handle(xxt_ADT xxt_handle) 499d71ae5a4SJacob Faibussowitsch { 50052f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX}; 501827bd09bSSatish Balay 5020924e98cSBarry Smith PetscFunctionBegin; 5037ee088dcSPierre Jolivet PetscCheck(xxt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", (void *)xxt_handle); 504827bd09bSSatish Balay 505827bd09bSSatish Balay vals[0] = vals[1] = xxt_handle->id; 506*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op)); 50763a3b9bcSJacob 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); 508*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509827bd09bSSatish Balay } 510827bd09bSSatish Balay 511d71ae5a4SJacob Faibussowitsch static PetscErrorCode det_separators(xxt_ADT xxt_handle) 512d71ae5a4SJacob Faibussowitsch { 51352f87cdaSBarry Smith PetscInt i, ct, id; 51452f87cdaSBarry Smith PetscInt mask, edge, *iptr; 51552f87cdaSBarry Smith PetscInt *dir, *used; 51652f87cdaSBarry Smith PetscInt sum[4], w[4]; 517a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 51852f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 519a501084fSBarry Smith PetscScalar *lhs, *rhs; 52052f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo = 0; 521ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 52252f87cdaSBarry Smith PetscInt *local2global = xxt_handle->mvi->local2global; 52352f87cdaSBarry Smith PetscInt n = xxt_handle->mvi->n; 52452f87cdaSBarry Smith PetscInt m = xxt_handle->mvi->m; 52552f87cdaSBarry Smith PetscInt level = xxt_handle->level; 526ab824b78SBarry Smith PetscInt shared = 0; 527827bd09bSSatish Balay 5280924e98cSBarry Smith PetscFunctionBegin; 52952f87cdaSBarry Smith dir = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 53052f87cdaSBarry Smith nsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 53152f87cdaSBarry Smith lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 53252f87cdaSBarry Smith fo = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1)); 53352f87cdaSBarry Smith used = (PetscInt *)malloc(sizeof(PetscInt) * n); 534827bd09bSSatish Balay 535*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(dir, level + 1)); 536*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(nsep, level + 1)); 537*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(lnsep, level + 1)); 538*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_set(fo, -1, n + 1)); 539*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(used, n)); 540827bd09bSSatish Balay 5418cda6cd7SBarry Smith lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 5428cda6cd7SBarry Smith rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 543827bd09bSSatish Balay 544827bd09bSSatish Balay /* determine the # of unique dof */ 545*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(lhs, m)); 546*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_set(lhs, 1.0, n)); 547*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level)); 548*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(rsum, 2)); 5493d3eaba7SBarry Smith for (i = 0; i < n; i++) { 5502fa5cd67SKarl Rupp if (lhs[i] != 0.0) { 5519371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 5529371c9d4SSatish Balay rsum[1] += lhs[i]; 5532fa5cd67SKarl Rupp } 554827bd09bSSatish Balay } 555*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(rsum, rw, 2, op, level)); 556827bd09bSSatish Balay rsum[0] += 0.1; 557827bd09bSSatish Balay rsum[1] += 0.1; 558827bd09bSSatish Balay 55977b4d14cSPeter Brune if (PetscAbsScalar(rsum[0] - rsum[1]) > EPS) shared = 1; 560827bd09bSSatish Balay 56152f87cdaSBarry Smith xxt_handle->info->n_global = xxt_handle->info->m_global = (PetscInt)rsum[0]; 56252f87cdaSBarry Smith xxt_handle->mvi->n_global = xxt_handle->mvi->m_global = (PetscInt)rsum[0]; 563827bd09bSSatish Balay 564827bd09bSSatish Balay /* determine separator sets top down */ 5652fa5cd67SKarl Rupp if (shared) { 566db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 567827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 568*3ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m)); 569*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge)); 570827bd09bSSatish Balay 571827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 572*3ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m)); 573*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge)); 574827bd09bSSatish Balay 575db4deed7SKarl Rupp for (i = 0; i < n; i++) { 576db4deed7SKarl Rupp if (id < mask) { 5772fa5cd67SKarl Rupp if (lhs[i] != 0.0) lhs[i] = 1.0; 578827bd09bSSatish Balay } 579db4deed7SKarl Rupp if (id >= mask) { 5802fa5cd67SKarl Rupp if (rhs[i] != 0.0) rhs[i] = 1.0; 581827bd09bSSatish Balay } 582827bd09bSSatish Balay } 583827bd09bSSatish Balay 584*3ba16761SJacob Faibussowitsch if (id < mask) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge - 1)); 585*3ba16761SJacob Faibussowitsch else PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge - 1)); 586827bd09bSSatish Balay 587827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 588*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(rsum, 4)); 589*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(sum, 4)); 590*3ba16761SJacob Faibussowitsch for (ct = i = 0; i < n; i++) { 591db4deed7SKarl Rupp if (!used[i]) { 592827bd09bSSatish Balay /* number of unmarked dofs on node */ 593827bd09bSSatish Balay ct++; 594827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 595db4deed7SKarl Rupp if (id < mask) { 5969371c9d4SSatish Balay if (lhs[i] != 0.0) { 5979371c9d4SSatish Balay sum[0]++; 5989371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 5999371c9d4SSatish Balay } 600827bd09bSSatish Balay } 601827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 602db4deed7SKarl Rupp if (id >= mask) { 6039371c9d4SSatish Balay if (rhs[i] != 0.0) { 6049371c9d4SSatish Balay sum[1]++; 6059371c9d4SSatish Balay rsum[1] += 1.0 / rhs[i]; 6069371c9d4SSatish Balay } 607827bd09bSSatish Balay } 608827bd09bSSatish Balay } 609827bd09bSSatish Balay } 610827bd09bSSatish Balay 611827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 612827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 613827bd09bSSatish Balay (id < mask) ? (rsum[2] = ct) : (rsum[3] = ct); 614*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge)); 615*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(rsum, rw, 4, op, edge)); 6169371c9d4SSatish Balay rsum[0] += 0.1; 6179371c9d4SSatish Balay rsum[1] += 0.1; 6189371c9d4SSatish Balay rsum[2] += 0.1; 6199371c9d4SSatish Balay rsum[3] += 0.1; 620827bd09bSSatish Balay 621db4deed7SKarl Rupp if (id < mask) { 622827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 623db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 624db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 6259371c9d4SSatish Balay ct++; 6269371c9d4SSatish Balay nfo++; 627827bd09bSSatish Balay 62808401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n"); 629827bd09bSSatish Balay 630827bd09bSSatish Balay *--iptr = local2global[i]; 631827bd09bSSatish Balay used[i] = edge; 632827bd09bSSatish Balay } 633827bd09bSSatish Balay } 634*3ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 635827bd09bSSatish Balay 636827bd09bSSatish Balay lnsep[edge] = ct; 63752f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[0]; 638827bd09bSSatish Balay dir[edge] = LEFT; 639827bd09bSSatish Balay } 640827bd09bSSatish Balay 641db4deed7SKarl Rupp if (id >= mask) { 642827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 643db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 644db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 6459371c9d4SSatish Balay ct++; 6469371c9d4SSatish Balay nfo++; 647827bd09bSSatish Balay 64808401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n"); 649827bd09bSSatish Balay 650827bd09bSSatish Balay *--iptr = local2global[i]; 651827bd09bSSatish Balay used[i] = edge; 652827bd09bSSatish Balay } 653827bd09bSSatish Balay } 654*3ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 655827bd09bSSatish Balay 656827bd09bSSatish Balay lnsep[edge] = ct; 65752f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[1]; 658827bd09bSSatish Balay dir[edge] = RIGHT; 659827bd09bSSatish Balay } 660827bd09bSSatish Balay 661827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 662827bd09bSSatish Balay /* do we need full set of separators for this? */ 663827bd09bSSatish Balay 664827bd09bSSatish Balay /* fold rhs hc into lower */ 6652fa5cd67SKarl Rupp if (id >= mask) id -= mask; 666827bd09bSSatish Balay } 667db4deed7SKarl Rupp } else { 668db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 669827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 670*3ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m)); 671*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge)); 672827bd09bSSatish Balay 673827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 674*3ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m)); 675*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge)); 676827bd09bSSatish Balay 677827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 678*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(sum, 4)); 679*3ba16761SJacob Faibussowitsch for (ct = i = 0; i < n; i++) { 680db4deed7SKarl Rupp if (!used[i]) { 681827bd09bSSatish Balay /* number of unmarked dofs on node */ 682827bd09bSSatish Balay ct++; 683827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6842fa5cd67SKarl Rupp if ((id < mask) && (lhs[i] != 0.0)) sum[0]++; 685827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6862fa5cd67SKarl Rupp if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++; 687827bd09bSSatish Balay } 688827bd09bSSatish Balay } 689827bd09bSSatish Balay 690827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 691827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 692*3ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge)); 693827bd09bSSatish Balay 694827bd09bSSatish Balay /* lhs hc wins */ 695db4deed7SKarl Rupp if (sum[2] >= sum[3]) { 696db4deed7SKarl Rupp if (id < mask) { 697827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 698db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 699db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 7009371c9d4SSatish Balay ct++; 7019371c9d4SSatish Balay nfo++; 702827bd09bSSatish Balay *--iptr = local2global[i]; 703827bd09bSSatish Balay used[i] = edge; 704827bd09bSSatish Balay } 705827bd09bSSatish Balay } 706*3ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 707827bd09bSSatish Balay lnsep[edge] = ct; 708827bd09bSSatish Balay } 709827bd09bSSatish Balay nsep[edge] = sum[0]; 710827bd09bSSatish Balay dir[edge] = LEFT; 711db4deed7SKarl Rupp } else { /* rhs hc wins */ 712db4deed7SKarl Rupp if (id >= mask) { 713827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 714db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 715db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 7169371c9d4SSatish Balay ct++; 7179371c9d4SSatish Balay nfo++; 718827bd09bSSatish Balay *--iptr = local2global[i]; 719827bd09bSSatish Balay used[i] = edge; 720827bd09bSSatish Balay } 721827bd09bSSatish Balay } 722*3ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 723827bd09bSSatish Balay lnsep[edge] = ct; 724827bd09bSSatish Balay } 725827bd09bSSatish Balay nsep[edge] = sum[1]; 726827bd09bSSatish Balay dir[edge] = RIGHT; 727827bd09bSSatish Balay } 728827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 729827bd09bSSatish Balay /* do we need full set of separators for this? */ 730827bd09bSSatish Balay 731827bd09bSSatish Balay /* fold rhs hc into lower */ 7322fa5cd67SKarl Rupp if (id >= mask) id -= mask; 733827bd09bSSatish Balay } 734827bd09bSSatish Balay } 735827bd09bSSatish Balay 736827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 737db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 738db4deed7SKarl Rupp if (!used[i]) { 7399371c9d4SSatish Balay ct++; 7409371c9d4SSatish Balay nfo++; 741827bd09bSSatish Balay *--iptr = local2global[i]; 742827bd09bSSatish Balay used[i] = edge; 743827bd09bSSatish Balay } 744827bd09bSSatish Balay } 745*3ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct)); 746827bd09bSSatish Balay lnsep[edge] = ct; 747827bd09bSSatish Balay nsep[edge] = ct; 748827bd09bSSatish Balay dir[edge] = LEFT; 749827bd09bSSatish Balay 750827bd09bSSatish Balay xxt_handle->info->nsep = nsep; 751827bd09bSSatish Balay xxt_handle->info->lnsep = lnsep; 752827bd09bSSatish Balay xxt_handle->info->fo = fo; 753827bd09bSSatish Balay xxt_handle->info->nfo = nfo; 754827bd09bSSatish Balay 755a501084fSBarry Smith free(dir); 756a501084fSBarry Smith free(lhs); 757a501084fSBarry Smith free(rhs); 758a501084fSBarry Smith free(used); 759*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 760827bd09bSSatish Balay } 761827bd09bSSatish Balay 762d71ae5a4SJacob Faibussowitsch static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data) 763d71ae5a4SJacob Faibussowitsch { 764827bd09bSSatish Balay mv_info *mvi; 765827bd09bSSatish Balay 766a501084fSBarry Smith mvi = (mv_info *)malloc(sizeof(mv_info)); 767827bd09bSSatish Balay mvi->n = n; 768827bd09bSSatish Balay mvi->m = m; 769827bd09bSSatish Balay mvi->n_global = -1; 770827bd09bSSatish Balay mvi->m_global = -1; 77152f87cdaSBarry Smith mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt)); 772ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global, local2global, m); 773827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7745c8f6a95SKarl Rupp mvi->matvec = matvec; 775827bd09bSSatish Balay mvi->grid_data = grid_data; 776827bd09bSSatish Balay 777827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */ 778ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 779827bd09bSSatish Balay 780827bd09bSSatish Balay return (mvi); 781827bd09bSSatish Balay } 782827bd09bSSatish Balay 783d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 784d71ae5a4SJacob Faibussowitsch { 7850924e98cSBarry Smith PetscFunctionBegin; 786*3ba16761SJacob Faibussowitsch PetscCall(A->matvec((mv_info *)A->grid_data, v, u)); 787*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 788827bd09bSSatish Balay } 789