xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision 7ee088dc06ee5b2041d9ac07d4e2beb519944c90)
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 {
84b1c944f5SJed Brown   PCTFS_comm_init();
85827bd09bSSatish Balay   check_handle(xxt_handle);
86827bd09bSSatish Balay 
87827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
8863a3b9bcSJacob 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);
89827bd09bSSatish Balay 
90827bd09bSSatish Balay   /* space for X info */
91a501084fSBarry Smith   xxt_handle->info = (xxt_info *)malloc(sizeof(xxt_info));
92827bd09bSSatish Balay 
93827bd09bSSatish Balay   /* set up matvec handles */
945c8f6a95SKarl Rupp   xxt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode(*)(mv_info *, PetscScalar *, PetscScalar *))matvec, grid_data);
95827bd09bSSatish Balay 
96827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
97827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
98827bd09bSSatish Balay   xxt_handle->ns = 0;
99827bd09bSSatish Balay 
100827bd09bSSatish Balay   /* determine separators and generate firing order - NB xxt info set here */
101827bd09bSSatish Balay   det_separators(xxt_handle);
102827bd09bSSatish Balay 
103827bd09bSSatish Balay   return (do_xxt_factor(xxt_handle));
104827bd09bSSatish Balay }
105827bd09bSSatish Balay 
106d71ae5a4SJacob Faibussowitsch PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b)
107d71ae5a4SJacob Faibussowitsch {
108b1c944f5SJed Brown   PCTFS_comm_init();
109827bd09bSSatish Balay   check_handle(xxt_handle);
110827bd09bSSatish Balay 
111827bd09bSSatish Balay   /* need to copy b into x? */
1122fa5cd67SKarl Rupp   if (b) PCTFS_rvec_copy(x, b, xxt_handle->mvi->n);
11338dcbb09SBarry Smith   return do_xxt_solve(xxt_handle, x);
114827bd09bSSatish Balay }
115827bd09bSSatish Balay 
116d71ae5a4SJacob Faibussowitsch PetscInt XXT_free(xxt_ADT xxt_handle)
117d71ae5a4SJacob Faibussowitsch {
118b1c944f5SJed Brown   PCTFS_comm_init();
119827bd09bSSatish Balay   check_handle(xxt_handle);
120827bd09bSSatish Balay   n_xxt_handles--;
121827bd09bSSatish Balay 
122a501084fSBarry Smith   free(xxt_handle->info->nsep);
123a501084fSBarry Smith   free(xxt_handle->info->lnsep);
124a501084fSBarry Smith   free(xxt_handle->info->fo);
125a501084fSBarry Smith   free(xxt_handle->info->stages);
126a501084fSBarry Smith   free(xxt_handle->info->solve_uu);
127a501084fSBarry Smith   free(xxt_handle->info->solve_w);
128a501084fSBarry Smith   free(xxt_handle->info->x);
129a501084fSBarry Smith   free(xxt_handle->info->col_vals);
130a501084fSBarry Smith   free(xxt_handle->info->col_sz);
131a501084fSBarry Smith   free(xxt_handle->info->col_indices);
132a501084fSBarry Smith   free(xxt_handle->info);
133a501084fSBarry Smith   free(xxt_handle->mvi->local2global);
134ca8e9878SJed Brown   PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle);
135a501084fSBarry Smith   free(xxt_handle->mvi);
136a501084fSBarry Smith   free(xxt_handle);
137827bd09bSSatish Balay 
138827bd09bSSatish Balay   /* if the check fails we nuke */
139a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
140827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
141827bd09bSSatish Balay   return (0);
142827bd09bSSatish Balay }
143827bd09bSSatish Balay 
14411cc89d2SBarry Smith /* This function is currently unused */
145d71ae5a4SJacob Faibussowitsch PetscErrorCode XXT_stats(xxt_ADT xxt_handle)
146d71ae5a4SJacob Faibussowitsch {
14752f87cdaSBarry Smith   PetscInt    op[]  = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD};
14852f87cdaSBarry Smith   PetscInt    fop[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD};
14952f87cdaSBarry Smith   PetscInt    vals[9], work[9];
150a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
151827bd09bSSatish Balay 
15263a3b9bcSJacob Faibussowitsch   PetscFunctionBegin;
15363a3b9bcSJacob Faibussowitsch   PetscCall(PCTFS_comm_init());
15463a3b9bcSJacob Faibussowitsch   PetscCall(check_handle(xxt_handle));
155827bd09bSSatish Balay 
156827bd09bSSatish Balay   /* if factorization not done there are no stats */
157db4deed7SKarl Rupp   if (!xxt_handle->info || !xxt_handle->mvi) {
1589566063dSJacob Faibussowitsch     if (!PCTFS_my_id) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "XXT_stats() :: no stats available!\n"));
15911cc89d2SBarry Smith     PetscFunctionReturn(0);
160827bd09bSSatish Balay   }
161827bd09bSSatish Balay 
162827bd09bSSatish Balay   vals[0] = vals[1] = vals[2] = xxt_handle->info->nnz;
163827bd09bSSatish Balay   vals[3] = vals[4] = vals[5] = xxt_handle->mvi->n;
164827bd09bSSatish Balay   vals[6] = vals[7] = vals[8] = xxt_handle->info->msg_buf_sz;
165dd39110bSPierre Jolivet   PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op);
166827bd09bSSatish Balay 
1672fa5cd67SKarl Rupp   fvals[0] = fvals[1] = fvals[2] = xxt_handle->info->tot_solve_time / xxt_handle->info->nsolves++;
168dd39110bSPierre Jolivet   PCTFS_grop(fvals, fwork, PETSC_STATIC_ARRAY_LENGTH(fop) - 1, fop);
169827bd09bSSatish Balay 
170db4deed7SKarl Rupp   if (!PCTFS_my_id) {
17163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[0]));
17263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[1]));
17363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xxt_nnz=%g\n", PCTFS_my_id, (double)(1.0 * vals[2] / PCTFS_num_nodes)));
17463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot   xxt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[2]));
17563a3b9bcSJacob 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)))));
17663a3b9bcSJacob 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)))));
17763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xxt_n  =%" PetscInt_FMT "\n", PCTFS_my_id, vals[3]));
17863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xxt_n  =%" PetscInt_FMT "\n", PCTFS_my_id, vals[4]));
17963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xxt_n  =%g\n", PCTFS_my_id, (double)(1.0 * vals[5] / PCTFS_num_nodes)));
18063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot   xxt_n  =%" PetscInt_FMT "\n", PCTFS_my_id, vals[5]));
18163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xxt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[6]));
18263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xxt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[7]));
18363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xxt_buf=%g\n", PCTFS_my_id, (double)(1.0 * vals[8] / PCTFS_num_nodes)));
18463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[0])));
18563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[1])));
18663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xxt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[2] / PCTFS_num_nodes)));
187827bd09bSSatish Balay   }
18863a3b9bcSJacob Faibussowitsch   PetscFunctionReturn(0);
189827bd09bSSatish Balay }
190827bd09bSSatish Balay 
19111cc89d2SBarry Smith /*
192827bd09bSSatish Balay 
193827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
194827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
195827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
196827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
197827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
198827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
199ca8e9878SJed Brown    PCTFS_gs_init/gop).
200827bd09bSSatish Balay 
201827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
202827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
20311cc89d2SBarry Smith */
204d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle)
205d71ae5a4SJacob Faibussowitsch {
2067b1ae94cSBarry Smith   return xxt_generate(xxt_handle);
207827bd09bSSatish Balay }
208827bd09bSSatish Balay 
209d71ae5a4SJacob Faibussowitsch static PetscErrorCode xxt_generate(xxt_ADT xxt_handle)
210d71ae5a4SJacob Faibussowitsch {
21152f87cdaSBarry Smith   PetscInt      i, j, k, idex;
21252f87cdaSBarry Smith   PetscInt      dim, col;
213a501084fSBarry Smith   PetscScalar  *u, *uu, *v, *z, *w, alpha, alpha_w;
21452f87cdaSBarry Smith   PetscInt     *segs;
21552f87cdaSBarry Smith   PetscInt      op[] = {GL_ADD, 0};
21652f87cdaSBarry Smith   PetscInt      off, len;
217a501084fSBarry Smith   PetscScalar  *x_ptr;
21852f87cdaSBarry Smith   PetscInt     *iptr, flag;
21952f87cdaSBarry Smith   PetscInt      start = 0, end, work;
22052f87cdaSBarry Smith   PetscInt      op2[] = {GL_MIN, 0};
221ca8e9878SJed Brown   PCTFS_gs_ADT  PCTFS_gs_handle;
22252f87cdaSBarry Smith   PetscInt     *nsep, *lnsep, *fo;
22352f87cdaSBarry Smith   PetscInt      a_n            = xxt_handle->mvi->n;
22452f87cdaSBarry Smith   PetscInt      a_m            = xxt_handle->mvi->m;
22552f87cdaSBarry Smith   PetscInt     *a_local2global = xxt_handle->mvi->local2global;
22652f87cdaSBarry Smith   PetscInt      level;
22752f87cdaSBarry Smith   PetscInt      xxt_nnz = 0, xxt_max_nnz = 0;
22852f87cdaSBarry Smith   PetscInt      n, m;
22952f87cdaSBarry Smith   PetscInt     *col_sz, *col_indices, *stages;
230a501084fSBarry Smith   PetscScalar **col_vals, *x;
23152f87cdaSBarry Smith   PetscInt      n_global;
23271a0148aSBarry Smith   PetscBLASInt  i1  = 1, dlen;
233a501084fSBarry Smith   PetscScalar   dm1 = -1.0;
234827bd09bSSatish Balay 
235827bd09bSSatish Balay   n               = xxt_handle->mvi->n;
236827bd09bSSatish Balay   nsep            = xxt_handle->info->nsep;
237827bd09bSSatish Balay   lnsep           = xxt_handle->info->lnsep;
238827bd09bSSatish Balay   fo              = xxt_handle->info->fo;
239827bd09bSSatish Balay   end             = lnsep[0];
240827bd09bSSatish Balay   level           = xxt_handle->level;
241ca8e9878SJed Brown   PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
242827bd09bSSatish Balay 
243827bd09bSSatish Balay   /* is there a null space? */
244827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
2452fa5cd67SKarl Rupp   for (i = 0, j = 0; i <= level; i++) j += nsep[i];
246827bd09bSSatish Balay 
247827bd09bSSatish Balay   m = j - xxt_handle->ns;
24848a46eb9SPierre 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));
249827bd09bSSatish Balay 
250827bd09bSSatish Balay   /* get and initialize storage for x local         */
251827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
25252f87cdaSBarry Smith   col_sz      = (PetscInt *)malloc(m * sizeof(PetscInt));
25352f87cdaSBarry Smith   col_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt));
254a501084fSBarry Smith   col_vals    = (PetscScalar **)malloc(m * sizeof(PetscScalar *));
255db4deed7SKarl Rupp   for (i = j = 0; i < m; i++, j += 2) {
256827bd09bSSatish Balay     col_indices[j] = col_indices[j + 1] = col_sz[i] = -1;
257827bd09bSSatish Balay     col_vals[i]                                     = NULL;
258827bd09bSSatish Balay   }
259827bd09bSSatish Balay   col_indices[j] = -1;
260827bd09bSSatish Balay 
261827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
262827bd09bSSatish Balay   /* this looks like nsep[]=segments */
26352f87cdaSBarry Smith   stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt));
26452f87cdaSBarry Smith   segs   = (PetscInt *)malloc((level + 1) * sizeof(PetscInt));
265ca8e9878SJed Brown   PCTFS_ivec_zero(stages, level + 1);
266ca8e9878SJed Brown   PCTFS_ivec_copy(segs, nsep, level + 1);
2672fa5cd67SKarl Rupp   for (i = 0; i < level; i++) segs[i + 1] += segs[i];
268827bd09bSSatish Balay   stages[0] = segs[0];
269827bd09bSSatish Balay 
270827bd09bSSatish Balay   /* temporary vectors  */
271a501084fSBarry Smith   u  = (PetscScalar *)malloc(n * sizeof(PetscScalar));
272a501084fSBarry Smith   z  = (PetscScalar *)malloc(n * sizeof(PetscScalar));
273a501084fSBarry Smith   v  = (PetscScalar *)malloc(a_m * sizeof(PetscScalar));
274a501084fSBarry Smith   uu = (PetscScalar *)malloc(m * sizeof(PetscScalar));
275a501084fSBarry Smith   w  = (PetscScalar *)malloc(m * sizeof(PetscScalar));
276827bd09bSSatish Balay 
277827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
2782fa5cd67SKarl Rupp   for (i = 1, j = 0; i <= level; i++) j += nsep[i];
279827bd09bSSatish Balay 
280827bd09bSSatish Balay   /* storage for sparse x values */
281827bd09bSSatish Balay   n_global    = xxt_handle->info->n_global;
28285ec1a3cSBarry Smith   xxt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes;
283a501084fSBarry Smith   x           = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar));
284827bd09bSSatish Balay   xxt_nnz     = 0;
285827bd09bSSatish Balay 
286827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
287827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
2882fa5cd67SKarl Rupp   for (dim = i = j = 0; i < m; i++) {
289827bd09bSSatish Balay     /* time to move to the next level? */
290d4af0d30SBarry Smith     while (i == segs[dim]) {
29108401ef6SPierre Jolivet       PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level");
292827bd09bSSatish Balay       stages[dim++] = i;
293827bd09bSSatish Balay       end += lnsep[dim];
294827bd09bSSatish Balay     }
295827bd09bSSatish Balay     stages[dim] = i;
296827bd09bSSatish Balay 
297827bd09bSSatish Balay     /* which column are we firing? */
298827bd09bSSatish Balay     /* i.e. set v_l */
299827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
300827bd09bSSatish Balay     (start < end) ? (col = fo[start]) : (col = INT_MAX);
301b1c944f5SJed Brown     PCTFS_giop_hc(&col, &work, 1, op2, dim);
302827bd09bSSatish Balay 
303827bd09bSSatish Balay     /* shouldn't need this */
304db4deed7SKarl Rupp     if (col == INT_MAX) {
3059566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n"));
306827bd09bSSatish Balay       continue;
307827bd09bSSatish Balay     }
308827bd09bSSatish Balay 
309827bd09bSSatish Balay     /* do I own it? I should */
310ca8e9878SJed Brown     PCTFS_rvec_zero(v, a_m);
311db4deed7SKarl Rupp     if (col == fo[start]) {
312827bd09bSSatish Balay       start++;
313ca8e9878SJed Brown       idex = PCTFS_ivec_linear_search(col, a_local2global, a_n);
3140fdf79fbSJacob Faibussowitsch       PetscCheck(idex != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!");
3159371c9d4SSatish Balay       v[idex] = 1.0;
3169371c9d4SSatish Balay       j++;
317db4deed7SKarl Rupp     } else {
318ca8e9878SJed Brown       idex = PCTFS_ivec_linear_search(col, a_local2global, a_m);
3192fa5cd67SKarl Rupp       if (idex != -1) v[idex] = 1.0;
320827bd09bSSatish Balay     }
321827bd09bSSatish Balay 
322827bd09bSSatish Balay     /* perform u = A.v_l */
323ca8e9878SJed Brown     PCTFS_rvec_zero(u, n);
324827bd09bSSatish Balay     do_matvec(xxt_handle->mvi, v, u);
325827bd09bSSatish Balay 
326827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
327827bd09bSSatish Balay     /* technically only need to zero out first i entries */
328827bd09bSSatish Balay     /* later turn this into an XXT_solve call ? */
329ca8e9878SJed Brown     PCTFS_rvec_zero(uu, m);
330827bd09bSSatish Balay     x_ptr = x;
331827bd09bSSatish Balay     iptr  = col_indices;
332db4deed7SKarl Rupp     for (k = 0; k < i; k++) {
333827bd09bSSatish Balay       off = *iptr++;
334827bd09bSSatish Balay       len = *iptr++;
3359566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(len, &dlen));
336792fecdfSBarry Smith       PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, x_ptr, &i1));
337827bd09bSSatish Balay       x_ptr += len;
338827bd09bSSatish Balay     }
339827bd09bSSatish Balay 
340827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
3419566063dSJacob Faibussowitsch     PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages));
342827bd09bSSatish Balay 
343827bd09bSSatish Balay     /* z = X.uu */
344ca8e9878SJed Brown     PCTFS_rvec_zero(z, n);
345827bd09bSSatish Balay     x_ptr = x;
346827bd09bSSatish Balay     iptr  = col_indices;
347db4deed7SKarl Rupp     for (k = 0; k < i; k++) {
348827bd09bSSatish Balay       off = *iptr++;
349827bd09bSSatish Balay       len = *iptr++;
3509566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(len, &dlen));
351792fecdfSBarry Smith       PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1));
352827bd09bSSatish Balay       x_ptr += len;
353827bd09bSSatish Balay     }
354827bd09bSSatish Balay 
355827bd09bSSatish Balay     /* compute v_l = v_l - z */
356ca8e9878SJed Brown     PCTFS_rvec_zero(v + a_n, a_m - a_n);
3579566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &dlen));
358792fecdfSBarry Smith     PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1));
359827bd09bSSatish Balay 
360827bd09bSSatish Balay     /* compute u_l = A.v_l */
3612fa5cd67SKarl Rupp     if (a_n != a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim);
362ca8e9878SJed Brown     PCTFS_rvec_zero(u, n);
363827bd09bSSatish Balay     do_matvec(xxt_handle->mvi, v, u);
364827bd09bSSatish Balay 
365827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
3669566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &dlen));
367792fecdfSBarry Smith     PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, v, &i1));
368827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
369b1c944f5SJed Brown     PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
370827bd09bSSatish Balay 
3718f1a2a5eSBarry Smith     alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha);
372827bd09bSSatish Balay 
373827bd09bSSatish Balay     /* check for small alpha                             */
374827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
37563a3b9bcSJacob Faibussowitsch     PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha));
376827bd09bSSatish Balay 
377827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
378ca8e9878SJed Brown     PCTFS_rvec_scale(v, 1.0 / alpha, n);
379827bd09bSSatish Balay 
380827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
381827bd09bSSatish Balay     flag = 1;
382827bd09bSSatish Balay     off = len = 0;
383db4deed7SKarl Rupp     for (k = 0; k < n; k++) {
384db4deed7SKarl Rupp       if (v[k] != 0.0) {
385827bd09bSSatish Balay         len = k;
3869371c9d4SSatish Balay         if (flag) {
3879371c9d4SSatish Balay           off  = k;
3889371c9d4SSatish Balay           flag = 0;
3899371c9d4SSatish Balay         }
390827bd09bSSatish Balay       }
391827bd09bSSatish Balay     }
392827bd09bSSatish Balay 
393827bd09bSSatish Balay     len -= (off - 1);
394827bd09bSSatish Balay 
395db4deed7SKarl Rupp     if (len > 0) {
396db4deed7SKarl Rupp       if ((xxt_nnz + len) > xxt_max_nnz) {
3979566063dSJacob Faibussowitsch         PetscCall(PetscInfo(0, "increasing space for X by 2x!\n"));
398827bd09bSSatish Balay         xxt_max_nnz *= 2;
399a501084fSBarry Smith         x_ptr = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar));
400ca8e9878SJed Brown         PCTFS_rvec_copy(x_ptr, x, xxt_nnz);
401a501084fSBarry Smith         free(x);
402827bd09bSSatish Balay         x = x_ptr;
403827bd09bSSatish Balay         x_ptr += xxt_nnz;
404827bd09bSSatish Balay       }
405827bd09bSSatish Balay       xxt_nnz += len;
406ca8e9878SJed Brown       PCTFS_rvec_copy(x_ptr, v + off, len);
407827bd09bSSatish Balay 
408827bd09bSSatish Balay       col_indices[2 * i] = off;
409827bd09bSSatish Balay       col_sz[i] = col_indices[2 * i + 1] = len;
410827bd09bSSatish Balay       col_vals[i]                        = x_ptr;
4119371c9d4SSatish Balay     } else {
412827bd09bSSatish Balay       col_indices[2 * i] = 0;
413827bd09bSSatish Balay       col_sz[i] = col_indices[2 * i + 1] = 0;
414827bd09bSSatish Balay       col_vals[i]                        = x_ptr;
415827bd09bSSatish Balay     }
416827bd09bSSatish Balay   }
417827bd09bSSatish Balay 
418827bd09bSSatish Balay   /* close off stages for execution phase */
419db4deed7SKarl Rupp   while (dim != level) {
420827bd09bSSatish Balay     stages[dim++] = i;
42163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level));
422827bd09bSSatish Balay   }
423827bd09bSSatish Balay   stages[dim] = i;
424827bd09bSSatish Balay 
425827bd09bSSatish Balay   xxt_handle->info->n              = xxt_handle->mvi->n;
426827bd09bSSatish Balay   xxt_handle->info->m              = m;
427827bd09bSSatish Balay   xxt_handle->info->nnz            = xxt_nnz;
428827bd09bSSatish Balay   xxt_handle->info->max_nnz        = xxt_max_nnz;
429827bd09bSSatish Balay   xxt_handle->info->msg_buf_sz     = stages[level] - stages[0];
430a501084fSBarry Smith   xxt_handle->info->solve_uu       = (PetscScalar *)malloc(m * sizeof(PetscScalar));
431a501084fSBarry Smith   xxt_handle->info->solve_w        = (PetscScalar *)malloc(m * sizeof(PetscScalar));
432827bd09bSSatish Balay   xxt_handle->info->x              = x;
433827bd09bSSatish Balay   xxt_handle->info->col_vals       = col_vals;
434827bd09bSSatish Balay   xxt_handle->info->col_sz         = col_sz;
435827bd09bSSatish Balay   xxt_handle->info->col_indices    = col_indices;
436827bd09bSSatish Balay   xxt_handle->info->stages         = stages;
437827bd09bSSatish Balay   xxt_handle->info->nsolves        = 0;
438827bd09bSSatish Balay   xxt_handle->info->tot_solve_time = 0.0;
439827bd09bSSatish Balay 
440a501084fSBarry Smith   free(segs);
441a501084fSBarry Smith   free(u);
442a501084fSBarry Smith   free(v);
443a501084fSBarry Smith   free(uu);
444a501084fSBarry Smith   free(z);
445a501084fSBarry Smith   free(w);
446827bd09bSSatish Balay 
447827bd09bSSatish Balay   return (0);
448827bd09bSSatish Balay }
449827bd09bSSatish Balay 
450d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc)
451d71ae5a4SJacob Faibussowitsch {
45252f87cdaSBarry Smith   PetscInt     off, len, *iptr;
45352f87cdaSBarry Smith   PetscInt     level       = xxt_handle->level;
45452f87cdaSBarry Smith   PetscInt     n           = xxt_handle->info->n;
45552f87cdaSBarry Smith   PetscInt     m           = xxt_handle->info->m;
45652f87cdaSBarry Smith   PetscInt    *stages      = xxt_handle->info->stages;
45752f87cdaSBarry Smith   PetscInt    *col_indices = xxt_handle->info->col_indices;
458a501084fSBarry Smith   PetscScalar *x_ptr, *uu_ptr;
459a501084fSBarry Smith   PetscScalar *solve_uu = xxt_handle->info->solve_uu;
460a501084fSBarry Smith   PetscScalar *solve_w  = xxt_handle->info->solve_w;
461a501084fSBarry Smith   PetscScalar *x        = xxt_handle->info->x;
46271a0148aSBarry Smith   PetscBLASInt i1       = 1, dlen;
463827bd09bSSatish Balay 
4640924e98cSBarry Smith   PetscFunctionBegin;
465827bd09bSSatish Balay   uu_ptr = solve_uu;
466ca8e9878SJed Brown   PCTFS_rvec_zero(uu_ptr, m);
467827bd09bSSatish Balay 
468827bd09bSSatish Balay   /* x  = X.Y^T.b */
469827bd09bSSatish Balay   /* uu = Y^T.b */
470db4deed7SKarl Rupp   for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) {
471c5df96a5SBarry Smith     off = *iptr++;
472c5df96a5SBarry Smith     len = *iptr++;
4739566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(len, &dlen));
474792fecdfSBarry Smith     PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, x_ptr, &i1));
475827bd09bSSatish Balay   }
476827bd09bSSatish Balay 
477d5b43468SJose E. Roman   /* communication of beta */
478827bd09bSSatish Balay   uu_ptr = solve_uu;
4799566063dSJacob Faibussowitsch   if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages));
480827bd09bSSatish Balay 
481ca8e9878SJed Brown   PCTFS_rvec_zero(uc, n);
482827bd09bSSatish Balay 
483827bd09bSSatish Balay   /* x = X.uu */
484db4deed7SKarl Rupp   for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) {
485c5df96a5SBarry Smith     off = *iptr++;
486c5df96a5SBarry Smith     len = *iptr++;
4879566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(len, &dlen));
488792fecdfSBarry Smith     PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1));
489827bd09bSSatish Balay   }
4900924e98cSBarry Smith   PetscFunctionReturn(0);
491827bd09bSSatish Balay }
492827bd09bSSatish Balay 
493d71ae5a4SJacob Faibussowitsch static PetscErrorCode check_handle(xxt_ADT xxt_handle)
494d71ae5a4SJacob Faibussowitsch {
49552f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX};
496827bd09bSSatish Balay 
4970924e98cSBarry Smith   PetscFunctionBegin;
498*7ee088dcSPierre Jolivet   PetscCheck(xxt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", (void *)xxt_handle);
499827bd09bSSatish Balay 
500827bd09bSSatish Balay   vals[0] = vals[1] = xxt_handle->id;
501dd39110bSPierre Jolivet   PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op);
50263a3b9bcSJacob 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);
5030924e98cSBarry Smith   PetscFunctionReturn(0);
504827bd09bSSatish Balay }
505827bd09bSSatish Balay 
506d71ae5a4SJacob Faibussowitsch static PetscErrorCode det_separators(xxt_ADT xxt_handle)
507d71ae5a4SJacob Faibussowitsch {
50852f87cdaSBarry Smith   PetscInt     i, ct, id;
50952f87cdaSBarry Smith   PetscInt     mask, edge, *iptr;
51052f87cdaSBarry Smith   PetscInt    *dir, *used;
51152f87cdaSBarry Smith   PetscInt     sum[4], w[4];
512a501084fSBarry Smith   PetscScalar  rsum[4], rw[4];
51352f87cdaSBarry Smith   PetscInt     op[] = {GL_ADD, 0};
514a501084fSBarry Smith   PetscScalar *lhs, *rhs;
51552f87cdaSBarry Smith   PetscInt    *nsep, *lnsep, *fo, nfo = 0;
516ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
51752f87cdaSBarry Smith   PetscInt    *local2global    = xxt_handle->mvi->local2global;
51852f87cdaSBarry Smith   PetscInt     n               = xxt_handle->mvi->n;
51952f87cdaSBarry Smith   PetscInt     m               = xxt_handle->mvi->m;
52052f87cdaSBarry Smith   PetscInt     level           = xxt_handle->level;
521ab824b78SBarry Smith   PetscInt     shared          = 0;
522827bd09bSSatish Balay 
5230924e98cSBarry Smith   PetscFunctionBegin;
52452f87cdaSBarry Smith   dir   = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
52552f87cdaSBarry Smith   nsep  = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
52652f87cdaSBarry Smith   lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
52752f87cdaSBarry Smith   fo    = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1));
52852f87cdaSBarry Smith   used  = (PetscInt *)malloc(sizeof(PetscInt) * n);
529827bd09bSSatish Balay 
530ca8e9878SJed Brown   PCTFS_ivec_zero(dir, level + 1);
531ca8e9878SJed Brown   PCTFS_ivec_zero(nsep, level + 1);
532ca8e9878SJed Brown   PCTFS_ivec_zero(lnsep, level + 1);
533ca8e9878SJed Brown   PCTFS_ivec_set(fo, -1, n + 1);
534ca8e9878SJed Brown   PCTFS_ivec_zero(used, n);
535827bd09bSSatish Balay 
5368cda6cd7SBarry Smith   lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m);
5378cda6cd7SBarry Smith   rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m);
538827bd09bSSatish Balay 
539827bd09bSSatish Balay   /* determine the # of unique dof */
540ca8e9878SJed Brown   PCTFS_rvec_zero(lhs, m);
541ca8e9878SJed Brown   PCTFS_rvec_set(lhs, 1.0, n);
542ca8e9878SJed Brown   PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level);
543ca8e9878SJed Brown   PCTFS_rvec_zero(rsum, 2);
5443d3eaba7SBarry Smith   for (i = 0; i < n; i++) {
5452fa5cd67SKarl Rupp     if (lhs[i] != 0.0) {
5469371c9d4SSatish Balay       rsum[0] += 1.0 / lhs[i];
5479371c9d4SSatish Balay       rsum[1] += lhs[i];
5482fa5cd67SKarl Rupp     }
549827bd09bSSatish Balay   }
550b1c944f5SJed Brown   PCTFS_grop_hc(rsum, rw, 2, op, level);
551827bd09bSSatish Balay   rsum[0] += 0.1;
552827bd09bSSatish Balay   rsum[1] += 0.1;
553827bd09bSSatish Balay 
55477b4d14cSPeter Brune   if (PetscAbsScalar(rsum[0] - rsum[1]) > EPS) shared = 1;
555827bd09bSSatish Balay 
55652f87cdaSBarry Smith   xxt_handle->info->n_global = xxt_handle->info->m_global = (PetscInt)rsum[0];
55752f87cdaSBarry Smith   xxt_handle->mvi->n_global = xxt_handle->mvi->m_global = (PetscInt)rsum[0];
558827bd09bSSatish Balay 
559827bd09bSSatish Balay   /* determine separator sets top down */
5602fa5cd67SKarl Rupp   if (shared) {
561db4deed7SKarl Rupp     for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) {
562827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
563ca8e9878SJed Brown       (id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m);
564ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge);
565827bd09bSSatish Balay 
566827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
567ca8e9878SJed Brown       (id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m);
568ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge);
569827bd09bSSatish Balay 
570db4deed7SKarl Rupp       for (i = 0; i < n; i++) {
571db4deed7SKarl Rupp         if (id < mask) {
5722fa5cd67SKarl Rupp           if (lhs[i] != 0.0) lhs[i] = 1.0;
573827bd09bSSatish Balay         }
574db4deed7SKarl Rupp         if (id >= mask) {
5752fa5cd67SKarl Rupp           if (rhs[i] != 0.0) rhs[i] = 1.0;
576827bd09bSSatish Balay         }
577827bd09bSSatish Balay       }
578827bd09bSSatish Balay 
5792fa5cd67SKarl Rupp       if (id < mask) PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge - 1);
5802fa5cd67SKarl Rupp       else PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge - 1);
581827bd09bSSatish Balay 
582827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
583ca8e9878SJed Brown       PCTFS_rvec_zero(rsum, 4);
584db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum, 4), ct = i = 0; i < n; i++) {
585db4deed7SKarl Rupp         if (!used[i]) {
586827bd09bSSatish Balay           /* number of unmarked dofs on node */
587827bd09bSSatish Balay           ct++;
588827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
589db4deed7SKarl Rupp           if (id < mask) {
5909371c9d4SSatish Balay             if (lhs[i] != 0.0) {
5919371c9d4SSatish Balay               sum[0]++;
5929371c9d4SSatish Balay               rsum[0] += 1.0 / lhs[i];
5939371c9d4SSatish Balay             }
594827bd09bSSatish Balay           }
595827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
596db4deed7SKarl Rupp           if (id >= mask) {
5979371c9d4SSatish Balay             if (rhs[i] != 0.0) {
5989371c9d4SSatish Balay               sum[1]++;
5999371c9d4SSatish Balay               rsum[1] += 1.0 / rhs[i];
6009371c9d4SSatish Balay             }
601827bd09bSSatish Balay           }
602827bd09bSSatish Balay         }
603827bd09bSSatish Balay       }
604827bd09bSSatish Balay 
605827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
606827bd09bSSatish Balay       (id < mask) ? (sum[2] = ct) : (sum[3] = ct);
607827bd09bSSatish Balay       (id < mask) ? (rsum[2] = ct) : (rsum[3] = ct);
608b1c944f5SJed Brown       PCTFS_giop_hc(sum, w, 4, op, edge);
609b1c944f5SJed Brown       PCTFS_grop_hc(rsum, rw, 4, op, edge);
6109371c9d4SSatish Balay       rsum[0] += 0.1;
6119371c9d4SSatish Balay       rsum[1] += 0.1;
6129371c9d4SSatish Balay       rsum[2] += 0.1;
6139371c9d4SSatish Balay       rsum[3] += 0.1;
614827bd09bSSatish Balay 
615db4deed7SKarl Rupp       if (id < mask) {
616827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
617db4deed7SKarl Rupp         for (ct = i = 0; i < n; i++) {
618db4deed7SKarl Rupp           if ((!used[i]) && (lhs[i] != 0.0)) {
6199371c9d4SSatish Balay             ct++;
6209371c9d4SSatish Balay             nfo++;
621827bd09bSSatish Balay 
62208401ef6SPierre Jolivet             PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n");
623827bd09bSSatish Balay 
624827bd09bSSatish Balay             *--iptr = local2global[i];
625827bd09bSSatish Balay             used[i] = edge;
626827bd09bSSatish Balay           }
627827bd09bSSatish Balay         }
6282fa5cd67SKarl Rupp         if (ct > 1) PCTFS_ivec_sort(iptr, ct);
629827bd09bSSatish Balay 
630827bd09bSSatish Balay         lnsep[edge] = ct;
63152f87cdaSBarry Smith         nsep[edge]  = (PetscInt)rsum[0];
632827bd09bSSatish Balay         dir[edge]   = LEFT;
633827bd09bSSatish Balay       }
634827bd09bSSatish Balay 
635db4deed7SKarl Rupp       if (id >= mask) {
636827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
637db4deed7SKarl Rupp         for (ct = i = 0; i < n; i++) {
638db4deed7SKarl Rupp           if ((!used[i]) && (rhs[i] != 0.0)) {
6399371c9d4SSatish Balay             ct++;
6409371c9d4SSatish Balay             nfo++;
641827bd09bSSatish Balay 
64208401ef6SPierre Jolivet             PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n");
643827bd09bSSatish Balay 
644827bd09bSSatish Balay             *--iptr = local2global[i];
645827bd09bSSatish Balay             used[i] = edge;
646827bd09bSSatish Balay           }
647827bd09bSSatish Balay         }
6482fa5cd67SKarl Rupp         if (ct > 1) PCTFS_ivec_sort(iptr, ct);
649827bd09bSSatish Balay 
650827bd09bSSatish Balay         lnsep[edge] = ct;
65152f87cdaSBarry Smith         nsep[edge]  = (PetscInt)rsum[1];
652827bd09bSSatish Balay         dir[edge]   = RIGHT;
653827bd09bSSatish Balay       }
654827bd09bSSatish Balay 
655827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
656827bd09bSSatish Balay       /* do we need full set of separators for this?                */
657827bd09bSSatish Balay 
658827bd09bSSatish Balay       /* fold rhs hc into lower */
6592fa5cd67SKarl Rupp       if (id >= mask) id -= mask;
660827bd09bSSatish Balay     }
661db4deed7SKarl Rupp   } else {
662db4deed7SKarl Rupp     for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) {
663827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
664ca8e9878SJed Brown       (id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m);
665ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge);
666827bd09bSSatish Balay 
667827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
668ca8e9878SJed Brown       (id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m);
669ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge);
670827bd09bSSatish Balay 
671827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
672db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum, 4), ct = i = 0; i < n; i++) {
673db4deed7SKarl Rupp         if (!used[i]) {
674827bd09bSSatish Balay           /* number of unmarked dofs on node */
675827bd09bSSatish Balay           ct++;
676827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
6772fa5cd67SKarl Rupp           if ((id < mask) && (lhs[i] != 0.0)) sum[0]++;
678827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
6792fa5cd67SKarl Rupp           if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++;
680827bd09bSSatish Balay         }
681827bd09bSSatish Balay       }
682827bd09bSSatish Balay 
683827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
684827bd09bSSatish Balay       (id < mask) ? (sum[2] = ct) : (sum[3] = ct);
685b1c944f5SJed Brown       PCTFS_giop_hc(sum, w, 4, op, edge);
686827bd09bSSatish Balay 
687827bd09bSSatish Balay       /* lhs hc wins */
688db4deed7SKarl Rupp       if (sum[2] >= sum[3]) {
689db4deed7SKarl Rupp         if (id < mask) {
690827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
691db4deed7SKarl Rupp           for (ct = i = 0; i < n; i++) {
692db4deed7SKarl Rupp             if ((!used[i]) && (lhs[i] != 0.0)) {
6939371c9d4SSatish Balay               ct++;
6949371c9d4SSatish Balay               nfo++;
695827bd09bSSatish Balay               *--iptr = local2global[i];
696827bd09bSSatish Balay               used[i] = edge;
697827bd09bSSatish Balay             }
698827bd09bSSatish Balay           }
6992fa5cd67SKarl Rupp           if (ct > 1) PCTFS_ivec_sort(iptr, ct);
700827bd09bSSatish Balay           lnsep[edge] = ct;
701827bd09bSSatish Balay         }
702827bd09bSSatish Balay         nsep[edge] = sum[0];
703827bd09bSSatish Balay         dir[edge]  = LEFT;
704db4deed7SKarl Rupp       } else { /* rhs hc wins */
705db4deed7SKarl Rupp         if (id >= mask) {
706827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
707db4deed7SKarl Rupp           for (ct = i = 0; i < n; i++) {
708db4deed7SKarl Rupp             if ((!used[i]) && (rhs[i] != 0.0)) {
7099371c9d4SSatish Balay               ct++;
7109371c9d4SSatish Balay               nfo++;
711827bd09bSSatish Balay               *--iptr = local2global[i];
712827bd09bSSatish Balay               used[i] = edge;
713827bd09bSSatish Balay             }
714827bd09bSSatish Balay           }
7152fa5cd67SKarl Rupp           if (ct > 1) PCTFS_ivec_sort(iptr, ct);
716827bd09bSSatish Balay           lnsep[edge] = ct;
717827bd09bSSatish Balay         }
718827bd09bSSatish Balay         nsep[edge] = sum[1];
719827bd09bSSatish Balay         dir[edge]  = RIGHT;
720827bd09bSSatish Balay       }
721827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
722827bd09bSSatish Balay       /* do we need full set of separators for this?                */
723827bd09bSSatish Balay 
724827bd09bSSatish Balay       /* fold rhs hc into lower */
7252fa5cd67SKarl Rupp       if (id >= mask) id -= mask;
726827bd09bSSatish Balay     }
727827bd09bSSatish Balay   }
728827bd09bSSatish Balay 
729827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
730db4deed7SKarl Rupp   for (ct = i = 0; i < n; i++) {
731db4deed7SKarl Rupp     if (!used[i]) {
7329371c9d4SSatish Balay       ct++;
7339371c9d4SSatish Balay       nfo++;
734827bd09bSSatish Balay       *--iptr = local2global[i];
735827bd09bSSatish Balay       used[i] = edge;
736827bd09bSSatish Balay     }
737827bd09bSSatish Balay   }
7382fa5cd67SKarl Rupp   if (ct > 1) PCTFS_ivec_sort(iptr, ct);
739827bd09bSSatish Balay   lnsep[edge] = ct;
740827bd09bSSatish Balay   nsep[edge]  = ct;
741827bd09bSSatish Balay   dir[edge]   = LEFT;
742827bd09bSSatish Balay 
743827bd09bSSatish Balay   xxt_handle->info->nsep  = nsep;
744827bd09bSSatish Balay   xxt_handle->info->lnsep = lnsep;
745827bd09bSSatish Balay   xxt_handle->info->fo    = fo;
746827bd09bSSatish Balay   xxt_handle->info->nfo   = nfo;
747827bd09bSSatish Balay 
748a501084fSBarry Smith   free(dir);
749a501084fSBarry Smith   free(lhs);
750a501084fSBarry Smith   free(rhs);
751a501084fSBarry Smith   free(used);
7520924e98cSBarry Smith   PetscFunctionReturn(0);
753827bd09bSSatish Balay }
754827bd09bSSatish Balay 
755d71ae5a4SJacob Faibussowitsch static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data)
756d71ae5a4SJacob Faibussowitsch {
757827bd09bSSatish Balay   mv_info *mvi;
758827bd09bSSatish Balay 
759a501084fSBarry Smith   mvi               = (mv_info *)malloc(sizeof(mv_info));
760827bd09bSSatish Balay   mvi->n            = n;
761827bd09bSSatish Balay   mvi->m            = m;
762827bd09bSSatish Balay   mvi->n_global     = -1;
763827bd09bSSatish Balay   mvi->m_global     = -1;
76452f87cdaSBarry Smith   mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt));
765ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global, local2global, m);
766827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
7675c8f6a95SKarl Rupp   mvi->matvec          = matvec;
768827bd09bSSatish Balay   mvi->grid_data       = grid_data;
769827bd09bSSatish Balay 
770827bd09bSSatish Balay   /* set xxt communication handle to perform restricted matvec */
771ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
772827bd09bSSatish Balay 
773827bd09bSSatish Balay   return (mvi);
774827bd09bSSatish Balay }
775827bd09bSSatish Balay 
776d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
777d71ae5a4SJacob Faibussowitsch {
7780924e98cSBarry Smith   PetscFunctionBegin;
779827bd09bSSatish Balay   A->matvec((mv_info *)A->grid_data, v, u);
7800924e98cSBarry Smith   PetscFunctionReturn(0);
781827bd09bSSatish Balay }
782