1827bd09bSSatish Balay 211cc89d2SBarry Smith /* 3827bd09bSSatish Balay Module Name: xyt 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 xyt_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 *xcol_sz, *xcol_indices; 30a501084fSBarry Smith PetscScalar **xcol_vals, *x, *solve_uu, *solve_w; 3152f87cdaSBarry Smith PetscInt *ycol_sz, *ycol_indices; 32a501084fSBarry Smith PetscScalar **ycol_vals, *y; 3352f87cdaSBarry Smith PetscInt nsolves; 34a501084fSBarry Smith PetscScalar tot_solve_time; 35827bd09bSSatish Balay } xyt_info; 36827bd09bSSatish Balay 37827bd09bSSatish Balay typedef struct matvec_info { 3852f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 3952f87cdaSBarry Smith PetscInt *local2global; 40ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 41a501084fSBarry Smith PetscErrorCode (*matvec)(struct matvec_info *, PetscScalar *, PetscScalar *); 42827bd09bSSatish Balay void *grid_data; 43827bd09bSSatish Balay } mv_info; 44827bd09bSSatish Balay 45827bd09bSSatish Balay struct xyt_CDT { 4652f87cdaSBarry Smith PetscInt id; 4752f87cdaSBarry Smith PetscInt ns; 4852f87cdaSBarry Smith PetscInt level; 49827bd09bSSatish Balay xyt_info *info; 50827bd09bSSatish Balay mv_info *mvi; 51827bd09bSSatish Balay }; 52827bd09bSSatish Balay 5352f87cdaSBarry Smith static PetscInt n_xyt = 0; 5452f87cdaSBarry Smith static PetscInt n_xyt_handles = 0; 55827bd09bSSatish Balay 56827bd09bSSatish Balay /* prototypes */ 573fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs); 583fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle); 593fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle); 603fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 6138dcbb09SBarry Smith static PetscErrorCode xyt_generate(xyt_ADT xyt_handle); 6238dcbb09SBarry Smith static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle); 635c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data); 64827bd09bSSatish Balay 65d71ae5a4SJacob Faibussowitsch xyt_ADT XYT_new(void) 66d71ae5a4SJacob Faibussowitsch { 67827bd09bSSatish Balay xyt_ADT xyt_handle; 68827bd09bSSatish Balay 69827bd09bSSatish Balay /* rolling count on n_xyt ... pot. problem here */ 70827bd09bSSatish Balay n_xyt_handles++; 71a501084fSBarry Smith xyt_handle = (xyt_ADT)malloc(sizeof(struct xyt_CDT)); 72827bd09bSSatish Balay xyt_handle->id = ++n_xyt; 73827bd09bSSatish Balay xyt_handle->info = NULL; 74827bd09bSSatish Balay xyt_handle->mvi = NULL; 75827bd09bSSatish Balay 76827bd09bSSatish Balay return (xyt_handle); 77827bd09bSSatish Balay } 78827bd09bSSatish Balay 7938dcbb09SBarry Smith PetscErrorCode XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt handle */ 8052f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */ 8152f87cdaSBarry Smith PetscInt n, /* local num rows */ 8252f87cdaSBarry Smith PetscInt m, /* local num cols */ 835c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc */ 841147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */ 85827bd09bSSatish Balay { 86b1c944f5SJed Brown PCTFS_comm_init(); 87827bd09bSSatish Balay check_handle(xyt_handle); 88827bd09bSSatish Balay 89827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 9063a3b9bcSJacob Faibussowitsch PetscCheck((1 << (xyt_handle->level = PCTFS_i_log2_num_nodes)) == PCTFS_num_nodes, PETSC_COMM_SELF, PETSC_ERR_PLIB, "only 2^k for now and MPI_COMM_WORLD!!! %d != %d", 1 << PCTFS_i_log2_num_nodes, PCTFS_num_nodes); 91827bd09bSSatish Balay 92827bd09bSSatish Balay /* space for X info */ 93a501084fSBarry Smith xyt_handle->info = (xyt_info *)malloc(sizeof(xyt_info)); 94827bd09bSSatish Balay 95827bd09bSSatish Balay /* set up matvec handles */ 965c8f6a95SKarl Rupp xyt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode(*)(mv_info *, PetscScalar *, PetscScalar *))matvec, grid_data); 97827bd09bSSatish Balay 98827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 99827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 100827bd09bSSatish Balay xyt_handle->ns = 0; 101827bd09bSSatish Balay 102827bd09bSSatish Balay /* determine separators and generate firing order - NB xyt info set here */ 103827bd09bSSatish Balay det_separators(xyt_handle); 104827bd09bSSatish Balay 105827bd09bSSatish Balay return (do_xyt_factor(xyt_handle)); 106827bd09bSSatish Balay } 107827bd09bSSatish Balay 108d71ae5a4SJacob Faibussowitsch PetscErrorCode XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b) 109d71ae5a4SJacob Faibussowitsch { 110b1c944f5SJed Brown PCTFS_comm_init(); 111827bd09bSSatish Balay check_handle(xyt_handle); 112827bd09bSSatish Balay 113827bd09bSSatish Balay /* need to copy b into x? */ 11459bc5b24SSatish Balay if (b) PCTFS_rvec_copy(x, b, xyt_handle->mvi->n); 11538dcbb09SBarry Smith return do_xyt_solve(xyt_handle, x); 116827bd09bSSatish Balay } 117827bd09bSSatish Balay 118d71ae5a4SJacob Faibussowitsch PetscErrorCode XYT_free(xyt_ADT xyt_handle) 119d71ae5a4SJacob Faibussowitsch { 120b1c944f5SJed Brown PCTFS_comm_init(); 121827bd09bSSatish Balay check_handle(xyt_handle); 122827bd09bSSatish Balay n_xyt_handles--; 123827bd09bSSatish Balay 124a501084fSBarry Smith free(xyt_handle->info->nsep); 125a501084fSBarry Smith free(xyt_handle->info->lnsep); 126a501084fSBarry Smith free(xyt_handle->info->fo); 127a501084fSBarry Smith free(xyt_handle->info->stages); 128a501084fSBarry Smith free(xyt_handle->info->solve_uu); 129a501084fSBarry Smith free(xyt_handle->info->solve_w); 130a501084fSBarry Smith free(xyt_handle->info->x); 131a501084fSBarry Smith free(xyt_handle->info->xcol_vals); 132a501084fSBarry Smith free(xyt_handle->info->xcol_sz); 133a501084fSBarry Smith free(xyt_handle->info->xcol_indices); 134a501084fSBarry Smith free(xyt_handle->info->y); 135a501084fSBarry Smith free(xyt_handle->info->ycol_vals); 136a501084fSBarry Smith free(xyt_handle->info->ycol_sz); 137a501084fSBarry Smith free(xyt_handle->info->ycol_indices); 138a501084fSBarry Smith free(xyt_handle->info); 139a501084fSBarry Smith free(xyt_handle->mvi->local2global); 140ca8e9878SJed Brown PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle); 141a501084fSBarry Smith free(xyt_handle->mvi); 142a501084fSBarry Smith free(xyt_handle); 143827bd09bSSatish Balay 144827bd09bSSatish Balay /* if the check fails we nuke */ 145a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 146827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 147827bd09bSSatish Balay return (0); 148827bd09bSSatish Balay } 149827bd09bSSatish Balay 15011cc89d2SBarry Smith /* This function is currently not used */ 151d71ae5a4SJacob Faibussowitsch PetscErrorCode XYT_stats(xyt_ADT xyt_handle) 152d71ae5a4SJacob Faibussowitsch { 15352f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD}; 15452f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD}; 15552f87cdaSBarry Smith PetscInt vals[9], work[9]; 156a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 157827bd09bSSatish Balay 15863a3b9bcSJacob Faibussowitsch PetscFunctionBegin; 15963a3b9bcSJacob Faibussowitsch PetscCall(PCTFS_comm_init()); 16063a3b9bcSJacob Faibussowitsch PetscCall(check_handle(xyt_handle)); 161827bd09bSSatish Balay 162827bd09bSSatish Balay /* if factorization not done there are no stats */ 1637d0a6c19SBarry Smith if (!xyt_handle->info || !xyt_handle->mvi) { 16463a3b9bcSJacob Faibussowitsch if (!PCTFS_my_id) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "XYT_stats() :: no stats available!\n")); 16511cc89d2SBarry Smith PetscFunctionReturn(0); 166827bd09bSSatish Balay } 167827bd09bSSatish Balay 168827bd09bSSatish Balay vals[0] = vals[1] = vals[2] = xyt_handle->info->nnz; 169827bd09bSSatish Balay vals[3] = vals[4] = vals[5] = xyt_handle->mvi->n; 170827bd09bSSatish Balay vals[6] = vals[7] = vals[8] = xyt_handle->info->msg_buf_sz; 17163a3b9bcSJacob Faibussowitsch PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op)); 172827bd09bSSatish Balay 1732fa5cd67SKarl Rupp fvals[0] = fvals[1] = fvals[2] = xyt_handle->info->tot_solve_time / xyt_handle->info->nsolves++; 17463a3b9bcSJacob Faibussowitsch PetscCall(PCTFS_grop(fvals, fwork, PETSC_STATIC_ARRAY_LENGTH(fop) - 1, fop)); 175827bd09bSSatish Balay 176b1c944f5SJed Brown if (!PCTFS_my_id) { 17763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xyt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[0])); 17863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xyt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[1])); 17963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xyt_nnz=%g\n", PCTFS_my_id, (double)(1.0 * vals[2] / PCTFS_num_nodes))); 18063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot xyt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[2])); 18163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: xyt C(2d) =%g\n", PCTFS_my_id, (double)(vals[2] / (PetscPowReal(1.0 * vals[5], 1.5))))); 18263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: xyt C(3d) =%g\n", PCTFS_my_id, (double)(vals[2] / (PetscPowReal(1.0 * vals[5], 1.6667))))); 18363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xyt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[3])); 18463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xyt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[4])); 18563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xyt_n =%g\n", PCTFS_my_id, (double)(1.0 * vals[5] / PCTFS_num_nodes))); 18663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot xyt_n =%" PetscInt_FMT "\n", PCTFS_my_id, vals[5])); 18763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xyt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[6])); 18863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xyt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[7])); 18963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xyt_buf=%g\n", PCTFS_my_id, (double)(1.0 * vals[8] / PCTFS_num_nodes))); 19063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min xyt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[0]))); 19163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max xyt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[1]))); 19263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg xyt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[2] / PCTFS_num_nodes))); 193827bd09bSSatish Balay } 19463a3b9bcSJacob Faibussowitsch PetscFunctionReturn(0); 195827bd09bSSatish Balay } 196827bd09bSSatish Balay 19711cc89d2SBarry Smith /* 198827bd09bSSatish Balay 199827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 200827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 201827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 202827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 203827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 204827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 205ca8e9878SJed Brown PCTFS_gs_init/gop). 206827bd09bSSatish Balay 207827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 208827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 20911cc89d2SBarry Smith */ 210d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle) 211d71ae5a4SJacob Faibussowitsch { 2127b1ae94cSBarry Smith return xyt_generate(xyt_handle); 213827bd09bSSatish Balay } 214827bd09bSSatish Balay 215d71ae5a4SJacob Faibussowitsch static PetscErrorCode xyt_generate(xyt_ADT xyt_handle) 216d71ae5a4SJacob Faibussowitsch { 21752f87cdaSBarry Smith PetscInt i, j, k, idx; 21852f87cdaSBarry Smith PetscInt dim, col; 219a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 22052f87cdaSBarry Smith PetscInt *segs; 22152f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 22252f87cdaSBarry Smith PetscInt off, len; 223a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr; 22452f87cdaSBarry Smith PetscInt *iptr, flag; 22552f87cdaSBarry Smith PetscInt start = 0, end, work; 22652f87cdaSBarry Smith PetscInt op2[] = {GL_MIN, 0}; 227ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 22852f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 22952f87cdaSBarry Smith PetscInt a_n = xyt_handle->mvi->n; 23052f87cdaSBarry Smith PetscInt a_m = xyt_handle->mvi->m; 23152f87cdaSBarry Smith PetscInt *a_local2global = xyt_handle->mvi->local2global; 23252f87cdaSBarry Smith PetscInt level; 23352f87cdaSBarry Smith PetscInt n, m; 23452f87cdaSBarry Smith PetscInt *xcol_sz, *xcol_indices, *stages; 235a501084fSBarry Smith PetscScalar **xcol_vals, *x; 23652f87cdaSBarry Smith PetscInt *ycol_sz, *ycol_indices; 237a501084fSBarry Smith PetscScalar **ycol_vals, *y; 23852f87cdaSBarry Smith PetscInt n_global; 23952f87cdaSBarry Smith PetscInt xt_nnz = 0, xt_max_nnz = 0; 24052f87cdaSBarry Smith PetscInt yt_nnz = 0, yt_max_nnz = 0; 2414a0de3f6SBarry Smith PetscBLASInt i1 = 1, dlen; 242a501084fSBarry Smith PetscScalar dm1 = -1.0; 243827bd09bSSatish Balay 244827bd09bSSatish Balay n = xyt_handle->mvi->n; 245827bd09bSSatish Balay nsep = xyt_handle->info->nsep; 246827bd09bSSatish Balay lnsep = xyt_handle->info->lnsep; 247827bd09bSSatish Balay fo = xyt_handle->info->fo; 248827bd09bSSatish Balay end = lnsep[0]; 249827bd09bSSatish Balay level = xyt_handle->level; 250ca8e9878SJed Brown PCTFS_gs_handle = xyt_handle->mvi->PCTFS_gs_handle; 251827bd09bSSatish Balay 252827bd09bSSatish Balay /* is there a null space? */ 253827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2542fa5cd67SKarl Rupp for (i = 0, j = 0; i <= level; i++) j += nsep[i]; 255827bd09bSSatish Balay 256827bd09bSSatish Balay m = j - xyt_handle->ns; 25748a46eb9SPierre Jolivet if (m != j) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "xyt_generate() :: null space exists %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, j, xyt_handle->ns)); 258827bd09bSSatish Balay 25963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(0, "xyt_generate() :: X(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", n, m)); 260827bd09bSSatish Balay 261827bd09bSSatish Balay /* get and initialize storage for x local */ 262827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 26352f87cdaSBarry Smith xcol_sz = (PetscInt *)malloc(m * sizeof(PetscInt)); 26452f87cdaSBarry Smith xcol_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt)); 265a501084fSBarry Smith xcol_vals = (PetscScalar **)malloc(m * sizeof(PetscScalar *)); 266db4deed7SKarl Rupp for (i = j = 0; i < m; i++, j += 2) { 267827bd09bSSatish Balay xcol_indices[j] = xcol_indices[j + 1] = xcol_sz[i] = -1; 268827bd09bSSatish Balay xcol_vals[i] = NULL; 269827bd09bSSatish Balay } 270827bd09bSSatish Balay xcol_indices[j] = -1; 271827bd09bSSatish Balay 272827bd09bSSatish Balay /* get and initialize storage for y local */ 273827bd09bSSatish Balay /* note that y local is nxm and stored by columns */ 27452f87cdaSBarry Smith ycol_sz = (PetscInt *)malloc(m * sizeof(PetscInt)); 27552f87cdaSBarry Smith ycol_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt)); 276a501084fSBarry Smith ycol_vals = (PetscScalar **)malloc(m * sizeof(PetscScalar *)); 277db4deed7SKarl Rupp for (i = j = 0; i < m; i++, j += 2) { 278827bd09bSSatish Balay ycol_indices[j] = ycol_indices[j + 1] = ycol_sz[i] = -1; 279827bd09bSSatish Balay ycol_vals[i] = NULL; 280827bd09bSSatish Balay } 281827bd09bSSatish Balay ycol_indices[j] = -1; 282827bd09bSSatish Balay 283827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 284827bd09bSSatish Balay /* this looks like nsep[]=segments */ 28552f87cdaSBarry Smith stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 28652f87cdaSBarry Smith segs = (PetscInt *)malloc((level + 1) * sizeof(PetscInt)); 287ca8e9878SJed Brown PCTFS_ivec_zero(stages, level + 1); 288ca8e9878SJed Brown PCTFS_ivec_copy(segs, nsep, level + 1); 2892fa5cd67SKarl Rupp for (i = 0; i < level; i++) segs[i + 1] += segs[i]; 290827bd09bSSatish Balay stages[0] = segs[0]; 291827bd09bSSatish Balay 292827bd09bSSatish Balay /* temporary vectors */ 293a501084fSBarry Smith u = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 294a501084fSBarry Smith z = (PetscScalar *)malloc(n * sizeof(PetscScalar)); 295a501084fSBarry Smith v = (PetscScalar *)malloc(a_m * sizeof(PetscScalar)); 296a501084fSBarry Smith uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 297a501084fSBarry Smith w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 298827bd09bSSatish Balay 299827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 3002fa5cd67SKarl Rupp for (i = 1, j = 0; i <= level; i++) j += nsep[i]; 301827bd09bSSatish Balay 302827bd09bSSatish Balay /* storage for sparse x values */ 303827bd09bSSatish Balay n_global = xyt_handle->info->n_global; 30485ec1a3cSBarry Smith xt_max_nnz = yt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes; 305a501084fSBarry Smith x = (PetscScalar *)malloc(xt_max_nnz * sizeof(PetscScalar)); 306a501084fSBarry Smith y = (PetscScalar *)malloc(yt_max_nnz * sizeof(PetscScalar)); 307827bd09bSSatish Balay 308827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 309827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 310db4deed7SKarl Rupp for (dim = i = j = 0; i < m; i++) { 311827bd09bSSatish Balay /* time to move to the next level? */ 312d4af0d30SBarry Smith while (i == segs[dim]) { 31308401ef6SPierre Jolivet PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level"); 314827bd09bSSatish Balay stages[dim++] = i; 315827bd09bSSatish Balay end += lnsep[dim]; 316827bd09bSSatish Balay } 317827bd09bSSatish Balay stages[dim] = i; 318827bd09bSSatish Balay 319827bd09bSSatish Balay /* which column are we firing? */ 320827bd09bSSatish Balay /* i.e. set v_l */ 321827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 322827bd09bSSatish Balay (start < end) ? (col = fo[start]) : (col = INT_MAX); 323b1c944f5SJed Brown PCTFS_giop_hc(&col, &work, 1, op2, dim); 324827bd09bSSatish Balay 325827bd09bSSatish Balay /* shouldn't need this */ 326db4deed7SKarl Rupp if (col == INT_MAX) { 3279566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n")); 328827bd09bSSatish Balay continue; 329827bd09bSSatish Balay } 330827bd09bSSatish Balay 331827bd09bSSatish Balay /* do I own it? I should */ 332ca8e9878SJed Brown PCTFS_rvec_zero(v, a_m); 3332fa5cd67SKarl Rupp if (col == fo[start]) { 334827bd09bSSatish Balay start++; 335ca8e9878SJed Brown idx = PCTFS_ivec_linear_search(col, a_local2global, a_n); 3362fa5cd67SKarl Rupp if (idx != -1) { 3372fa5cd67SKarl Rupp v[idx] = 1.0; 3382fa5cd67SKarl Rupp j++; 339546078acSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!"); 340db4deed7SKarl Rupp } else { 341ca8e9878SJed Brown idx = PCTFS_ivec_linear_search(col, a_local2global, a_m); 3422fa5cd67SKarl Rupp if (idx != -1) v[idx] = 1.0; 343827bd09bSSatish Balay } 344827bd09bSSatish Balay 345827bd09bSSatish Balay /* perform u = A.v_l */ 346ca8e9878SJed Brown PCTFS_rvec_zero(u, n); 347827bd09bSSatish Balay do_matvec(xyt_handle->mvi, v, u); 348827bd09bSSatish Balay 349827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 350827bd09bSSatish Balay /* technically only need to zero out first i entries */ 351827bd09bSSatish Balay /* later turn this into an XYT_solve call ? */ 352ca8e9878SJed Brown PCTFS_rvec_zero(uu, m); 353827bd09bSSatish Balay y_ptr = y; 354827bd09bSSatish Balay iptr = ycol_indices; 355db4deed7SKarl Rupp for (k = 0; k < i; k++) { 356827bd09bSSatish Balay off = *iptr++; 357827bd09bSSatish Balay len = *iptr++; 3589566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 359792fecdfSBarry Smith PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, y_ptr, &i1)); 360827bd09bSSatish Balay y_ptr += len; 361827bd09bSSatish Balay } 362827bd09bSSatish Balay 363827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 3649566063dSJacob Faibussowitsch PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages)); 365827bd09bSSatish Balay 366827bd09bSSatish Balay /* z = X.uu */ 367ca8e9878SJed Brown PCTFS_rvec_zero(z, n); 368827bd09bSSatish Balay x_ptr = x; 369827bd09bSSatish Balay iptr = xcol_indices; 370db4deed7SKarl Rupp for (k = 0; k < i; k++) { 371827bd09bSSatish Balay off = *iptr++; 372827bd09bSSatish Balay len = *iptr++; 3739566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 374792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1)); 375827bd09bSSatish Balay x_ptr += len; 376827bd09bSSatish Balay } 377827bd09bSSatish Balay 378827bd09bSSatish Balay /* compute v_l = v_l - z */ 379ca8e9878SJed Brown PCTFS_rvec_zero(v + a_n, a_m - a_n); 3809566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 381792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1)); 382827bd09bSSatish Balay 383827bd09bSSatish Balay /* compute u_l = A.v_l */ 3842fa5cd67SKarl Rupp if (a_n != a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim); 385ca8e9878SJed Brown PCTFS_rvec_zero(u, n); 386827bd09bSSatish Balay do_matvec(xyt_handle->mvi, v, u); 387827bd09bSSatish Balay 388827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 3899566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen)); 390792fecdfSBarry Smith PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, u, &i1)); 391827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 392b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 393827bd09bSSatish Balay 3948f1a2a5eSBarry Smith alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha); 395827bd09bSSatish Balay 396827bd09bSSatish Balay /* check for small alpha */ 397827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 39863a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha)); 399827bd09bSSatish Balay 400827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 401ca8e9878SJed Brown PCTFS_rvec_scale(v, 1.0 / alpha, n); 402ca8e9878SJed Brown PCTFS_rvec_scale(u, 1.0 / alpha, n); 403827bd09bSSatish Balay 404827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 405827bd09bSSatish Balay flag = 1; 406827bd09bSSatish Balay off = len = 0; 407db4deed7SKarl Rupp for (k = 0; k < n; k++) { 408db4deed7SKarl Rupp if (v[k] != 0.0) { 409827bd09bSSatish Balay len = k; 4109371c9d4SSatish Balay if (flag) { 4119371c9d4SSatish Balay off = k; 4129371c9d4SSatish Balay flag = 0; 4139371c9d4SSatish Balay } 414827bd09bSSatish Balay } 415827bd09bSSatish Balay } 416827bd09bSSatish Balay 417827bd09bSSatish Balay len -= (off - 1); 418827bd09bSSatish Balay 419db4deed7SKarl Rupp if (len > 0) { 420db4deed7SKarl Rupp if ((xt_nnz + len) > xt_max_nnz) { 4219566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "increasing space for X by 2x!\n")); 422827bd09bSSatish Balay xt_max_nnz *= 2; 423a501084fSBarry Smith x_ptr = (PetscScalar *)malloc(xt_max_nnz * sizeof(PetscScalar)); 424ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr, x, xt_nnz); 425a501084fSBarry Smith free(x); 426827bd09bSSatish Balay x = x_ptr; 427827bd09bSSatish Balay x_ptr += xt_nnz; 428827bd09bSSatish Balay } 429827bd09bSSatish Balay xt_nnz += len; 430ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr, v + off, len); 431827bd09bSSatish Balay 432827bd09bSSatish Balay xcol_indices[2 * i] = off; 433827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2 * i + 1] = len; 434827bd09bSSatish Balay xcol_vals[i] = x_ptr; 435db4deed7SKarl Rupp } else { 436827bd09bSSatish Balay xcol_indices[2 * i] = 0; 437827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2 * i + 1] = 0; 438827bd09bSSatish Balay xcol_vals[i] = x_ptr; 439827bd09bSSatish Balay } 440827bd09bSSatish Balay 441827bd09bSSatish Balay /* add newly generated column, u_l, to Y */ 442827bd09bSSatish Balay flag = 1; 443827bd09bSSatish Balay off = len = 0; 444db4deed7SKarl Rupp for (k = 0; k < n; k++) { 445db4deed7SKarl Rupp if (u[k] != 0.0) { 446827bd09bSSatish Balay len = k; 4479371c9d4SSatish Balay if (flag) { 4489371c9d4SSatish Balay off = k; 4499371c9d4SSatish Balay flag = 0; 4509371c9d4SSatish Balay } 451827bd09bSSatish Balay } 452827bd09bSSatish Balay } 453827bd09bSSatish Balay 454827bd09bSSatish Balay len -= (off - 1); 455827bd09bSSatish Balay 456db4deed7SKarl Rupp if (len > 0) { 457db4deed7SKarl Rupp if ((yt_nnz + len) > yt_max_nnz) { 4589566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "increasing space for Y by 2x!\n")); 459827bd09bSSatish Balay yt_max_nnz *= 2; 460a501084fSBarry Smith y_ptr = (PetscScalar *)malloc(yt_max_nnz * sizeof(PetscScalar)); 461ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr, y, yt_nnz); 462a501084fSBarry Smith free(y); 463827bd09bSSatish Balay y = y_ptr; 464827bd09bSSatish Balay y_ptr += yt_nnz; 465827bd09bSSatish Balay } 466827bd09bSSatish Balay yt_nnz += len; 467ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr, u + off, len); 468827bd09bSSatish Balay 469827bd09bSSatish Balay ycol_indices[2 * i] = off; 470827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2 * i + 1] = len; 471827bd09bSSatish Balay ycol_vals[i] = y_ptr; 472db4deed7SKarl Rupp } else { 473827bd09bSSatish Balay ycol_indices[2 * i] = 0; 474827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2 * i + 1] = 0; 475827bd09bSSatish Balay ycol_vals[i] = y_ptr; 476827bd09bSSatish Balay } 477827bd09bSSatish Balay } 478827bd09bSSatish Balay 479827bd09bSSatish Balay /* close off stages for execution phase */ 480db4deed7SKarl Rupp while (dim != level) { 481827bd09bSSatish Balay stages[dim++] = i; 48263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level)); 483827bd09bSSatish Balay } 484827bd09bSSatish Balay stages[dim] = i; 485827bd09bSSatish Balay 486827bd09bSSatish Balay xyt_handle->info->n = xyt_handle->mvi->n; 487827bd09bSSatish Balay xyt_handle->info->m = m; 488827bd09bSSatish Balay xyt_handle->info->nnz = xt_nnz + yt_nnz; 489827bd09bSSatish Balay xyt_handle->info->max_nnz = xt_max_nnz + yt_max_nnz; 490827bd09bSSatish Balay xyt_handle->info->msg_buf_sz = stages[level] - stages[0]; 491a501084fSBarry Smith xyt_handle->info->solve_uu = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 492a501084fSBarry Smith xyt_handle->info->solve_w = (PetscScalar *)malloc(m * sizeof(PetscScalar)); 493827bd09bSSatish Balay xyt_handle->info->x = x; 494827bd09bSSatish Balay xyt_handle->info->xcol_vals = xcol_vals; 495827bd09bSSatish Balay xyt_handle->info->xcol_sz = xcol_sz; 496827bd09bSSatish Balay xyt_handle->info->xcol_indices = xcol_indices; 497827bd09bSSatish Balay xyt_handle->info->stages = stages; 498827bd09bSSatish Balay xyt_handle->info->y = y; 499827bd09bSSatish Balay xyt_handle->info->ycol_vals = ycol_vals; 500827bd09bSSatish Balay xyt_handle->info->ycol_sz = ycol_sz; 501827bd09bSSatish Balay xyt_handle->info->ycol_indices = ycol_indices; 502827bd09bSSatish Balay 503a501084fSBarry Smith free(segs); 504a501084fSBarry Smith free(u); 505a501084fSBarry Smith free(v); 506a501084fSBarry Smith free(uu); 507a501084fSBarry Smith free(z); 508a501084fSBarry Smith free(w); 509827bd09bSSatish Balay 510827bd09bSSatish Balay return (0); 511827bd09bSSatish Balay } 512827bd09bSSatish Balay 513d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 514d71ae5a4SJacob Faibussowitsch { 51552f87cdaSBarry Smith PetscInt off, len, *iptr; 51652f87cdaSBarry Smith PetscInt level = xyt_handle->level; 51752f87cdaSBarry Smith PetscInt n = xyt_handle->info->n; 51852f87cdaSBarry Smith PetscInt m = xyt_handle->info->m; 51952f87cdaSBarry Smith PetscInt *stages = xyt_handle->info->stages; 52052f87cdaSBarry Smith PetscInt *xcol_indices = xyt_handle->info->xcol_indices; 52152f87cdaSBarry Smith PetscInt *ycol_indices = xyt_handle->info->ycol_indices; 522a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr, *uu_ptr; 523a501084fSBarry Smith PetscScalar *solve_uu = xyt_handle->info->solve_uu; 524a501084fSBarry Smith PetscScalar *solve_w = xyt_handle->info->solve_w; 525a501084fSBarry Smith PetscScalar *x = xyt_handle->info->x; 526a501084fSBarry Smith PetscScalar *y = xyt_handle->info->y; 5274a0de3f6SBarry Smith PetscBLASInt i1 = 1, dlen; 528827bd09bSSatish Balay 5290924e98cSBarry Smith PetscFunctionBegin; 530827bd09bSSatish Balay uu_ptr = solve_uu; 531ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr, m); 532827bd09bSSatish Balay 533827bd09bSSatish Balay /* x = X.Y^T.b */ 534827bd09bSSatish Balay /* uu = Y^T.b */ 535947b95d8SBarry Smith for (y_ptr = y, iptr = ycol_indices; *iptr != -1; y_ptr += len) { 536c5df96a5SBarry Smith off = *iptr++; 537c5df96a5SBarry Smith len = *iptr++; 5389566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 539792fecdfSBarry Smith PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, y_ptr, &i1)); 540827bd09bSSatish Balay } 541827bd09bSSatish Balay 542d5b43468SJose E. Roman /* communication of beta */ 543827bd09bSSatish Balay uu_ptr = solve_uu; 5449566063dSJacob Faibussowitsch if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages)); 545ca8e9878SJed Brown PCTFS_rvec_zero(uc, n); 546827bd09bSSatish Balay 547827bd09bSSatish Balay /* x = X.uu */ 548db4deed7SKarl Rupp for (x_ptr = x, iptr = xcol_indices; *iptr != -1; x_ptr += len) { 549c5df96a5SBarry Smith off = *iptr++; 550c5df96a5SBarry Smith len = *iptr++; 5519566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen)); 552792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1)); 553827bd09bSSatish Balay } 5540924e98cSBarry Smith PetscFunctionReturn(0); 555827bd09bSSatish Balay } 556827bd09bSSatish Balay 557d71ae5a4SJacob Faibussowitsch static PetscErrorCode check_handle(xyt_ADT xyt_handle) 558d71ae5a4SJacob Faibussowitsch { 55952f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX}; 560827bd09bSSatish Balay 5610924e98cSBarry Smith PetscFunctionBegin; 56263a3b9bcSJacob Faibussowitsch PetscCheck(xyt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", xyt_handle); 563827bd09bSSatish Balay 564827bd09bSSatish Balay vals[0] = vals[1] = xyt_handle->id; 565dd39110bSPierre Jolivet PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op); 56663a3b9bcSJacob Faibussowitsch PetscCheck(!(vals[0] != vals[1]) && !(xyt_handle->id <= 0), PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: id mismatch min/max %" PetscInt_FMT "/%" PetscInt_FMT " %" PetscInt_FMT, vals[0], vals[1], xyt_handle->id); 5670924e98cSBarry Smith PetscFunctionReturn(0); 568827bd09bSSatish Balay } 569827bd09bSSatish Balay 570d71ae5a4SJacob Faibussowitsch static PetscErrorCode det_separators(xyt_ADT xyt_handle) 571d71ae5a4SJacob Faibussowitsch { 57252f87cdaSBarry Smith PetscInt i, ct, id; 57352f87cdaSBarry Smith PetscInt mask, edge, *iptr; 57452f87cdaSBarry Smith PetscInt *dir, *used; 57552f87cdaSBarry Smith PetscInt sum[4], w[4]; 576a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 57752f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0}; 578a501084fSBarry Smith PetscScalar *lhs, *rhs; 57952f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo = 0; 580ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle = xyt_handle->mvi->PCTFS_gs_handle; 58152f87cdaSBarry Smith PetscInt *local2global = xyt_handle->mvi->local2global; 58252f87cdaSBarry Smith PetscInt n = xyt_handle->mvi->n; 58352f87cdaSBarry Smith PetscInt m = xyt_handle->mvi->m; 58452f87cdaSBarry Smith PetscInt level = xyt_handle->level; 585ab824b78SBarry Smith PetscInt shared = 0; 586827bd09bSSatish Balay 5870924e98cSBarry Smith PetscFunctionBegin; 58852f87cdaSBarry Smith dir = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 58952f87cdaSBarry Smith nsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 59052f87cdaSBarry Smith lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1)); 59152f87cdaSBarry Smith fo = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1)); 59252f87cdaSBarry Smith used = (PetscInt *)malloc(sizeof(PetscInt) * n); 593827bd09bSSatish Balay 594ca8e9878SJed Brown PCTFS_ivec_zero(dir, level + 1); 595ca8e9878SJed Brown PCTFS_ivec_zero(nsep, level + 1); 596ca8e9878SJed Brown PCTFS_ivec_zero(lnsep, level + 1); 597ca8e9878SJed Brown PCTFS_ivec_set(fo, -1, n + 1); 598ca8e9878SJed Brown PCTFS_ivec_zero(used, n); 599827bd09bSSatish Balay 6008cda6cd7SBarry Smith lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 6018cda6cd7SBarry Smith rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m); 602827bd09bSSatish Balay 603827bd09bSSatish Balay /* determine the # of unique dof */ 604ca8e9878SJed Brown PCTFS_rvec_zero(lhs, m); 605ca8e9878SJed Brown PCTFS_rvec_set(lhs, 1.0, n); 606ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level); 6079566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "done first PCTFS_gs_gop_hc\n")); 608ca8e9878SJed Brown PCTFS_rvec_zero(rsum, 2); 609947b95d8SBarry Smith for (i = 0; i < n; i++) { 6109371c9d4SSatish Balay if (lhs[i] != 0.0) { 6119371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i]; 6129371c9d4SSatish Balay rsum[1] += lhs[i]; 6139371c9d4SSatish Balay } 6142fa5cd67SKarl Rupp if (lhs[i] != 1.0) shared = 1; 615827bd09bSSatish Balay } 616827bd09bSSatish Balay 617b1c944f5SJed Brown PCTFS_grop_hc(rsum, rw, 2, op, level); 618827bd09bSSatish Balay rsum[0] += 0.1; 619827bd09bSSatish Balay rsum[1] += 0.1; 620827bd09bSSatish Balay 62152f87cdaSBarry Smith xyt_handle->info->n_global = xyt_handle->info->m_global = (PetscInt)rsum[0]; 62252f87cdaSBarry Smith xyt_handle->mvi->n_global = xyt_handle->mvi->m_global = (PetscInt)rsum[0]; 623827bd09bSSatish Balay 624827bd09bSSatish Balay /* determine separator sets top down */ 625*0fdf79fbSJacob Faibussowitsch PetscCheck(!shared, PETSC_COMM_SELF, PETSC_ERR_PLIB, "shared dof separator determination not ready ... see hmt!!!"); 626827bd09bSSatish Balay /* solution is to do as in the symmetric shared case but then */ 627827bd09bSSatish Balay /* pick the sub-hc with the most free dofs and do a mat-vec */ 628827bd09bSSatish Balay /* and pick up the responses on the other sub-hc from the */ 629827bd09bSSatish Balay /* initial separator set obtained from the symm. shared case */ 6304e619c20SJed Brown /* [dead code deleted since it is unlikely to be completed] */ 631db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) { 632827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 633ca8e9878SJed Brown (id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m); 634ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge); 635827bd09bSSatish Balay 636827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 637ca8e9878SJed Brown (id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m); 638ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge); 639827bd09bSSatish Balay 640827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 641db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum, 4), ct = i = 0; i < n; i++) { 642db4deed7SKarl Rupp if (!used[i]) { 643827bd09bSSatish Balay /* number of unmarked dofs on node */ 644827bd09bSSatish Balay ct++; 645827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6462fa5cd67SKarl Rupp if ((id < mask) && (lhs[i] != 0.0)) sum[0]++; 647827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6482fa5cd67SKarl Rupp if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++; 649827bd09bSSatish Balay } 650827bd09bSSatish Balay } 651827bd09bSSatish Balay 652827bd09bSSatish Balay /* for the non-symmetric case we need separators of width 2 */ 653827bd09bSSatish Balay /* so take both sides */ 654827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct); 655b1c944f5SJed Brown PCTFS_giop_hc(sum, w, 4, op, edge); 656827bd09bSSatish Balay 657827bd09bSSatish Balay ct = 0; 658db4deed7SKarl Rupp if (id < mask) { 659827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 660db4deed7SKarl Rupp for (i = 0; i < n; i++) { 661db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) { 6629371c9d4SSatish Balay ct++; 6639371c9d4SSatish Balay nfo++; 664827bd09bSSatish Balay *--iptr = local2global[i]; 665827bd09bSSatish Balay used[i] = edge; 666827bd09bSSatish Balay } 667827bd09bSSatish Balay } 668827bd09bSSatish Balay /* LSH hc summation of ct should be sum[0] */ 669db4deed7SKarl Rupp } else { 670827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 671db4deed7SKarl Rupp for (i = 0; i < n; i++) { 672db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) { 6739371c9d4SSatish Balay ct++; 6749371c9d4SSatish Balay nfo++; 675827bd09bSSatish Balay *--iptr = local2global[i]; 676827bd09bSSatish Balay used[i] = edge; 677827bd09bSSatish Balay } 678827bd09bSSatish Balay } 679827bd09bSSatish Balay /* RSH hc summation of ct should be sum[1] */ 680827bd09bSSatish Balay } 681827bd09bSSatish Balay 6822fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 683827bd09bSSatish Balay lnsep[edge] = ct; 684827bd09bSSatish Balay nsep[edge] = sum[0] + sum[1]; 685827bd09bSSatish Balay dir[edge] = BOTH; 686827bd09bSSatish Balay 687827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 688827bd09bSSatish Balay /* do we need full set of separators for this? */ 689827bd09bSSatish Balay 690827bd09bSSatish Balay /* fold rhs hc into lower */ 6912fa5cd67SKarl Rupp if (id >= mask) id -= mask; 692827bd09bSSatish Balay } 693827bd09bSSatish Balay 694827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 695db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) { 696db4deed7SKarl Rupp if (!used[i]) { 6979371c9d4SSatish Balay ct++; 6989371c9d4SSatish Balay nfo++; 699827bd09bSSatish Balay *--iptr = local2global[i]; 700827bd09bSSatish Balay used[i] = edge; 701827bd09bSSatish Balay } 702827bd09bSSatish Balay } 7032fa5cd67SKarl Rupp if (ct > 1) PCTFS_ivec_sort(iptr, ct); 704827bd09bSSatish Balay lnsep[edge] = ct; 705827bd09bSSatish Balay nsep[edge] = ct; 706827bd09bSSatish Balay dir[edge] = BOTH; 707827bd09bSSatish Balay 708827bd09bSSatish Balay xyt_handle->info->nsep = nsep; 709827bd09bSSatish Balay xyt_handle->info->lnsep = lnsep; 710827bd09bSSatish Balay xyt_handle->info->fo = fo; 711827bd09bSSatish Balay xyt_handle->info->nfo = nfo; 712827bd09bSSatish Balay 713a501084fSBarry Smith free(dir); 714a501084fSBarry Smith free(lhs); 715a501084fSBarry Smith free(rhs); 716a501084fSBarry Smith free(used); 7170924e98cSBarry Smith PetscFunctionReturn(0); 718827bd09bSSatish Balay } 719827bd09bSSatish Balay 720d71ae5a4SJacob Faibussowitsch static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data) 721d71ae5a4SJacob Faibussowitsch { 722827bd09bSSatish Balay mv_info *mvi; 723827bd09bSSatish Balay 724a501084fSBarry Smith mvi = (mv_info *)malloc(sizeof(mv_info)); 725827bd09bSSatish Balay mvi->n = n; 726827bd09bSSatish Balay mvi->m = m; 727827bd09bSSatish Balay mvi->n_global = -1; 728827bd09bSSatish Balay mvi->m_global = -1; 72952f87cdaSBarry Smith mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt)); 7302fa5cd67SKarl Rupp 731ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global, local2global, m); 732827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7335c8f6a95SKarl Rupp mvi->matvec = matvec; 734827bd09bSSatish Balay mvi->grid_data = grid_data; 735827bd09bSSatish Balay 736827bd09bSSatish Balay /* set xyt communication handle to perform restricted matvec */ 737ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 738827bd09bSSatish Balay 739827bd09bSSatish Balay return (mvi); 740827bd09bSSatish Balay } 741827bd09bSSatish Balay 742d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 743d71ae5a4SJacob Faibussowitsch { 7440924e98cSBarry Smith PetscFunctionBegin; 745827bd09bSSatish Balay A->matvec((mv_info *)A->grid_data, v, u); 7460924e98cSBarry Smith PetscFunctionReturn(0); 747827bd09bSSatish Balay } 748