xref: /petsc/src/ksp/pc/impls/tfs/xyt.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
111cc89d2SBarry Smith /*
2827bd09bSSatish Balay Module Name: xyt
3827bd09bSSatish Balay Module Info:
4827bd09bSSatish Balay 
5827bd09bSSatish Balay author:  Henry M. Tufo III
6827bd09bSSatish Balay e-mail:  hmt@asci.uchicago.edu
7827bd09bSSatish Balay contact:
8827bd09bSSatish Balay +--------------------------------+--------------------------------+
9827bd09bSSatish Balay |MCS Division - Building 221     |Department of Computer Science  |
10827bd09bSSatish Balay |Argonne National Laboratory     |Ryerson 152                     |
11827bd09bSSatish Balay |9700 S. Cass Avenue             |The University of Chicago       |
12827bd09bSSatish Balay |Argonne, IL  60439              |Chicago, IL  60637              |
13827bd09bSSatish Balay |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
14827bd09bSSatish Balay +--------------------------------+--------------------------------+
15827bd09bSSatish Balay 
16827bd09bSSatish Balay Last Modification: 3.20.01
1711cc89d2SBarry Smith */
18c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
19827bd09bSSatish Balay 
20827bd09bSSatish Balay #define LEFT  -1
21827bd09bSSatish Balay #define RIGHT 1
22827bd09bSSatish Balay #define BOTH  0
23827bd09bSSatish Balay 
24827bd09bSSatish Balay typedef struct xyt_solver_info {
2552f87cdaSBarry Smith   PetscInt      n, m, n_global, m_global;
2652f87cdaSBarry Smith   PetscInt      nnz, max_nnz, msg_buf_sz;
2752f87cdaSBarry Smith   PetscInt     *nsep, *lnsep, *fo, nfo, *stages;
2852f87cdaSBarry Smith   PetscInt     *xcol_sz, *xcol_indices;
29a501084fSBarry Smith   PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
3052f87cdaSBarry Smith   PetscInt     *ycol_sz, *ycol_indices;
31a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
3252f87cdaSBarry Smith   PetscInt      nsolves;
33a501084fSBarry Smith   PetscScalar   tot_solve_time;
34827bd09bSSatish Balay } xyt_info;
35827bd09bSSatish Balay 
36827bd09bSSatish Balay typedef struct matvec_info {
3752f87cdaSBarry Smith   PetscInt     n, m, n_global, m_global;
3852f87cdaSBarry Smith   PetscInt    *local2global;
39ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle;
40a501084fSBarry Smith   PetscErrorCode (*matvec)(struct matvec_info *, PetscScalar *, PetscScalar *);
41827bd09bSSatish Balay   void *grid_data;
42827bd09bSSatish Balay } mv_info;
43827bd09bSSatish Balay 
44827bd09bSSatish Balay struct xyt_CDT {
4552f87cdaSBarry Smith   PetscInt  id;
4652f87cdaSBarry Smith   PetscInt  ns;
4752f87cdaSBarry Smith   PetscInt  level;
48827bd09bSSatish Balay   xyt_info *info;
49827bd09bSSatish Balay   mv_info  *mvi;
50827bd09bSSatish Balay };
51827bd09bSSatish Balay 
5252f87cdaSBarry Smith static PetscInt n_xyt         = 0;
5352f87cdaSBarry Smith static PetscInt n_xyt_handles = 0;
54827bd09bSSatish Balay 
55827bd09bSSatish Balay /* prototypes */
563fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
573fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle);
583fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle);
593fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
6038dcbb09SBarry Smith static PetscErrorCode xyt_generate(xyt_ADT xyt_handle);
6138dcbb09SBarry Smith static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle);
625c8f6a95SKarl Rupp static mv_info       *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data);
63827bd09bSSatish Balay 
XYT_new(void)64d71ae5a4SJacob Faibussowitsch xyt_ADT XYT_new(void)
65d71ae5a4SJacob Faibussowitsch {
66827bd09bSSatish Balay   xyt_ADT xyt_handle;
67827bd09bSSatish Balay 
68827bd09bSSatish Balay   /* rolling count on n_xyt ... pot. problem here */
69827bd09bSSatish Balay   n_xyt_handles++;
70a501084fSBarry Smith   xyt_handle       = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
71827bd09bSSatish Balay   xyt_handle->id   = ++n_xyt;
72827bd09bSSatish Balay   xyt_handle->info = NULL;
73827bd09bSSatish Balay   xyt_handle->mvi  = NULL;
74827bd09bSSatish Balay 
754ad8454bSPierre Jolivet   return xyt_handle;
76827bd09bSSatish Balay }
77827bd09bSSatish Balay 
XYT_factor(xyt_ADT xyt_handle,PetscInt * local2global,PetscInt n,PetscInt m,PetscErrorCode (* matvec)(void *,PetscScalar *,PetscScalar *),void * grid_data)7838dcbb09SBarry Smith PetscErrorCode XYT_factor(xyt_ADT   xyt_handle,                                           /* prev. allocated xyt  handle */
7952f87cdaSBarry Smith                           PetscInt *local2global,                                         /* global column mapping       */
8052f87cdaSBarry Smith                           PetscInt  n,                                                    /* local num rows              */
8152f87cdaSBarry Smith                           PetscInt  m,                                                    /* local num cols              */
825c8f6a95SKarl Rupp                           PetscErrorCode (*matvec)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc         */
831147fc2aSKarl Rupp                           void *grid_data)                                                /* grid data for matvec        */
84827bd09bSSatish Balay {
853ba16761SJacob Faibussowitsch   PetscFunctionBegin;
863ba16761SJacob Faibussowitsch   PetscCall(PCTFS_comm_init());
873ba16761SJacob Faibussowitsch   PetscCall(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 */
1033ba16761SJacob Faibussowitsch   PetscCall(det_separators(xyt_handle));
104827bd09bSSatish Balay 
1053ba16761SJacob Faibussowitsch   PetscCall(do_xyt_factor(xyt_handle));
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107827bd09bSSatish Balay }
108827bd09bSSatish Balay 
XYT_solve(xyt_ADT xyt_handle,PetscScalar * x,PetscScalar * b)109d71ae5a4SJacob Faibussowitsch PetscErrorCode XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b)
110d71ae5a4SJacob Faibussowitsch {
1113ba16761SJacob Faibussowitsch   PetscFunctionBegin;
1123ba16761SJacob Faibussowitsch   PetscCall(PCTFS_comm_init());
1133ba16761SJacob Faibussowitsch   PetscCall(check_handle(xyt_handle));
114827bd09bSSatish Balay 
115827bd09bSSatish Balay   /* need to copy b into x? */
1163ba16761SJacob Faibussowitsch   if (b) PetscCall(PCTFS_rvec_copy(x, b, xyt_handle->mvi->n));
1173ba16761SJacob Faibussowitsch   PetscCall(do_xyt_solve(xyt_handle, x));
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119827bd09bSSatish Balay }
120827bd09bSSatish Balay 
XYT_free(xyt_ADT xyt_handle)121d71ae5a4SJacob Faibussowitsch PetscErrorCode XYT_free(xyt_ADT xyt_handle)
122d71ae5a4SJacob Faibussowitsch {
1233ba16761SJacob Faibussowitsch   PetscFunctionBegin;
1243ba16761SJacob Faibussowitsch   PetscCall(PCTFS_comm_init());
1253ba16761SJacob Faibussowitsch   PetscCall(check_handle(xyt_handle));
126827bd09bSSatish Balay   n_xyt_handles--;
127827bd09bSSatish Balay 
128a501084fSBarry Smith   free(xyt_handle->info->nsep);
129a501084fSBarry Smith   free(xyt_handle->info->lnsep);
130a501084fSBarry Smith   free(xyt_handle->info->fo);
131a501084fSBarry Smith   free(xyt_handle->info->stages);
132a501084fSBarry Smith   free(xyt_handle->info->solve_uu);
133a501084fSBarry Smith   free(xyt_handle->info->solve_w);
134a501084fSBarry Smith   free(xyt_handle->info->x);
135a501084fSBarry Smith   free(xyt_handle->info->xcol_vals);
136a501084fSBarry Smith   free(xyt_handle->info->xcol_sz);
137a501084fSBarry Smith   free(xyt_handle->info->xcol_indices);
138a501084fSBarry Smith   free(xyt_handle->info->y);
139a501084fSBarry Smith   free(xyt_handle->info->ycol_vals);
140a501084fSBarry Smith   free(xyt_handle->info->ycol_sz);
141a501084fSBarry Smith   free(xyt_handle->info->ycol_indices);
142a501084fSBarry Smith   free(xyt_handle->info);
143a501084fSBarry Smith   free(xyt_handle->mvi->local2global);
1443ba16761SJacob Faibussowitsch   PetscCall(PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle));
145a501084fSBarry Smith   free(xyt_handle->mvi);
146a501084fSBarry Smith   free(xyt_handle);
147827bd09bSSatish Balay 
148827bd09bSSatish Balay   /* if the check fails we nuke */
149a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
150827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
1513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152827bd09bSSatish Balay }
153827bd09bSSatish Balay 
15411cc89d2SBarry Smith /* This function is currently not used */
XYT_stats(xyt_ADT xyt_handle)155d71ae5a4SJacob Faibussowitsch PetscErrorCode XYT_stats(xyt_ADT xyt_handle)
156d71ae5a4SJacob Faibussowitsch {
15752f87cdaSBarry Smith   PetscInt    op[]  = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD, GL_MIN, GL_MAX, GL_ADD};
15852f87cdaSBarry Smith   PetscInt    fop[] = {NON_UNIFORM, GL_MIN, GL_MAX, GL_ADD};
15952f87cdaSBarry Smith   PetscInt    vals[9], work[9];
160a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
161827bd09bSSatish Balay 
16263a3b9bcSJacob Faibussowitsch   PetscFunctionBegin;
16363a3b9bcSJacob Faibussowitsch   PetscCall(PCTFS_comm_init());
16463a3b9bcSJacob Faibussowitsch   PetscCall(check_handle(xyt_handle));
165827bd09bSSatish Balay 
166827bd09bSSatish Balay   /* if factorization not done there are no stats */
1677d0a6c19SBarry Smith   if (!xyt_handle->info || !xyt_handle->mvi) {
16863a3b9bcSJacob Faibussowitsch     if (!PCTFS_my_id) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "XYT_stats() :: no stats available!\n"));
1693ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
170827bd09bSSatish Balay   }
171827bd09bSSatish Balay 
172827bd09bSSatish Balay   vals[0] = vals[1] = vals[2] = xyt_handle->info->nnz;
173827bd09bSSatish Balay   vals[3] = vals[4] = vals[5] = xyt_handle->mvi->n;
174827bd09bSSatish Balay   vals[6] = vals[7] = vals[8] = xyt_handle->info->msg_buf_sz;
17563a3b9bcSJacob Faibussowitsch   PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op));
176827bd09bSSatish Balay 
1772fa5cd67SKarl Rupp   fvals[0] = fvals[1] = fvals[2] = xyt_handle->info->tot_solve_time / xyt_handle->info->nsolves++;
17863a3b9bcSJacob Faibussowitsch   PetscCall(PCTFS_grop(fvals, fwork, PETSC_STATIC_ARRAY_LENGTH(fop) - 1, fop));
179827bd09bSSatish Balay 
180b1c944f5SJed Brown   if (!PCTFS_my_id) {
18163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xyt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[0]));
18263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xyt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[1]));
18363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xyt_nnz=%g\n", PCTFS_my_id, (double)(1.0 * vals[2] / PCTFS_num_nodes)));
18463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot   xyt_nnz=%" PetscInt_FMT "\n", PCTFS_my_id, vals[2]));
18563a3b9bcSJacob 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)))));
18663a3b9bcSJacob 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)))));
18763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xyt_n  =%" PetscInt_FMT "\n", PCTFS_my_id, vals[3]));
18863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xyt_n  =%" PetscInt_FMT "\n", PCTFS_my_id, vals[4]));
18963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xyt_n  =%g\n", PCTFS_my_id, (double)(1.0 * vals[5] / PCTFS_num_nodes)));
19063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: tot   xyt_n  =%" PetscInt_FMT "\n", PCTFS_my_id, vals[5]));
19163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xyt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[6]));
19263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xyt_buf=%" PetscInt_FMT "\n", PCTFS_my_id, vals[7]));
19363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xyt_buf=%g\n", PCTFS_my_id, (double)(1.0 * vals[8] / PCTFS_num_nodes)));
19463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: min   xyt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[0])));
19563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: max   xyt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[1])));
19663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d :: avg   xyt_slv=%g\n", PCTFS_my_id, (double)PetscRealPart(fvals[2] / PCTFS_num_nodes)));
197827bd09bSSatish Balay   }
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
199827bd09bSSatish Balay }
200827bd09bSSatish Balay 
20111cc89d2SBarry Smith /*
202827bd09bSSatish Balay 
203827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
204827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
205827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
206827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
207827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
208827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
209ca8e9878SJed Brown    PCTFS_gs_init/gop).
210827bd09bSSatish Balay 
211827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
212827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
21311cc89d2SBarry Smith */
do_xyt_factor(xyt_ADT xyt_handle)214d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle)
215d71ae5a4SJacob Faibussowitsch {
2167b1ae94cSBarry Smith   return xyt_generate(xyt_handle);
217827bd09bSSatish Balay }
218827bd09bSSatish Balay 
xyt_generate(xyt_ADT xyt_handle)219d71ae5a4SJacob Faibussowitsch static PetscErrorCode xyt_generate(xyt_ADT xyt_handle)
220d71ae5a4SJacob Faibussowitsch {
22152f87cdaSBarry Smith   PetscInt      i, j, k, idx;
22252f87cdaSBarry Smith   PetscInt      dim, col;
223a501084fSBarry Smith   PetscScalar  *u, *uu, *v, *z, *w, alpha, alpha_w;
22452f87cdaSBarry Smith   PetscInt     *segs;
22552f87cdaSBarry Smith   PetscInt      op[] = {GL_ADD, 0};
22652f87cdaSBarry Smith   PetscInt      off, len;
227a501084fSBarry Smith   PetscScalar  *x_ptr, *y_ptr;
22852f87cdaSBarry Smith   PetscInt     *iptr, flag;
22952f87cdaSBarry Smith   PetscInt      start = 0, end, work;
23052f87cdaSBarry Smith   PetscInt      op2[] = {GL_MIN, 0};
231ca8e9878SJed Brown   PCTFS_gs_ADT  PCTFS_gs_handle;
23252f87cdaSBarry Smith   PetscInt     *nsep, *lnsep, *fo;
23352f87cdaSBarry Smith   PetscInt      a_n            = xyt_handle->mvi->n;
23452f87cdaSBarry Smith   PetscInt      a_m            = xyt_handle->mvi->m;
23552f87cdaSBarry Smith   PetscInt     *a_local2global = xyt_handle->mvi->local2global;
23652f87cdaSBarry Smith   PetscInt      level;
23752f87cdaSBarry Smith   PetscInt      n, m;
23852f87cdaSBarry Smith   PetscInt     *xcol_sz, *xcol_indices, *stages;
239a501084fSBarry Smith   PetscScalar **xcol_vals, *x;
24052f87cdaSBarry Smith   PetscInt     *ycol_sz, *ycol_indices;
241a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
24252f87cdaSBarry Smith   PetscInt      n_global;
24352f87cdaSBarry Smith   PetscInt      xt_nnz = 0, xt_max_nnz = 0;
24452f87cdaSBarry Smith   PetscInt      yt_nnz = 0, yt_max_nnz = 0;
2454a0de3f6SBarry Smith   PetscBLASInt  i1  = 1, dlen;
246a501084fSBarry Smith   PetscScalar   dm1 = -1.0;
247827bd09bSSatish Balay 
248827bd09bSSatish Balay   n               = xyt_handle->mvi->n;
249827bd09bSSatish Balay   nsep            = xyt_handle->info->nsep;
250827bd09bSSatish Balay   lnsep           = xyt_handle->info->lnsep;
251827bd09bSSatish Balay   fo              = xyt_handle->info->fo;
252827bd09bSSatish Balay   end             = lnsep[0];
253827bd09bSSatish Balay   level           = xyt_handle->level;
254ca8e9878SJed Brown   PCTFS_gs_handle = xyt_handle->mvi->PCTFS_gs_handle;
255827bd09bSSatish Balay 
256827bd09bSSatish Balay   /* is there a null space? */
257827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
2582fa5cd67SKarl Rupp   for (i = 0, j = 0; i <= level; i++) j += nsep[i];
259827bd09bSSatish Balay 
260827bd09bSSatish Balay   m = j - xyt_handle->ns;
26148a46eb9SPierre 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));
262827bd09bSSatish Balay 
26363a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(0, "xyt_generate() :: X(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", n, m));
264827bd09bSSatish Balay 
265827bd09bSSatish Balay   /* get and initialize storage for x local         */
266827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
26752f87cdaSBarry Smith   xcol_sz      = (PetscInt *)malloc(m * sizeof(PetscInt));
26852f87cdaSBarry Smith   xcol_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt));
269a501084fSBarry Smith   xcol_vals    = (PetscScalar **)malloc(m * sizeof(PetscScalar *));
270db4deed7SKarl Rupp   for (i = j = 0; i < m; i++, j += 2) {
271827bd09bSSatish Balay     xcol_indices[j] = xcol_indices[j + 1] = xcol_sz[i] = -1;
272827bd09bSSatish Balay     xcol_vals[i]                                       = NULL;
273827bd09bSSatish Balay   }
274827bd09bSSatish Balay   xcol_indices[j] = -1;
275827bd09bSSatish Balay 
276827bd09bSSatish Balay   /* get and initialize storage for y local         */
277827bd09bSSatish Balay   /* note that y local is nxm and stored by columns */
27852f87cdaSBarry Smith   ycol_sz      = (PetscInt *)malloc(m * sizeof(PetscInt));
27952f87cdaSBarry Smith   ycol_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt));
280a501084fSBarry Smith   ycol_vals    = (PetscScalar **)malloc(m * sizeof(PetscScalar *));
281db4deed7SKarl Rupp   for (i = j = 0; i < m; i++, j += 2) {
282827bd09bSSatish Balay     ycol_indices[j] = ycol_indices[j + 1] = ycol_sz[i] = -1;
283827bd09bSSatish Balay     ycol_vals[i]                                       = NULL;
284827bd09bSSatish Balay   }
285827bd09bSSatish Balay   ycol_indices[j] = -1;
286827bd09bSSatish Balay 
287827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
288827bd09bSSatish Balay   /* this looks like nsep[]=segments */
28952f87cdaSBarry Smith   stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt));
29052f87cdaSBarry Smith   segs   = (PetscInt *)malloc((level + 1) * sizeof(PetscInt));
2913ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_zero(stages, level + 1));
292ca8e9878SJed Brown   PCTFS_ivec_copy(segs, nsep, level + 1);
2932fa5cd67SKarl Rupp   for (i = 0; i < level; i++) segs[i + 1] += segs[i];
294827bd09bSSatish Balay   stages[0] = segs[0];
295827bd09bSSatish Balay 
296827bd09bSSatish Balay   /* temporary vectors  */
297a501084fSBarry Smith   u  = (PetscScalar *)malloc(n * sizeof(PetscScalar));
298a501084fSBarry Smith   z  = (PetscScalar *)malloc(n * sizeof(PetscScalar));
299a501084fSBarry Smith   v  = (PetscScalar *)malloc(a_m * sizeof(PetscScalar));
300a501084fSBarry Smith   uu = (PetscScalar *)malloc(m * sizeof(PetscScalar));
301a501084fSBarry Smith   w  = (PetscScalar *)malloc(m * sizeof(PetscScalar));
302827bd09bSSatish Balay 
303827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
3042fa5cd67SKarl Rupp   for (i = 1, j = 0; i <= level; i++) j += nsep[i];
305827bd09bSSatish Balay 
306827bd09bSSatish Balay   /* storage for sparse x values */
307827bd09bSSatish Balay   n_global   = xyt_handle->info->n_global;
30885ec1a3cSBarry Smith   xt_max_nnz = yt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes;
309a501084fSBarry Smith   x                       = (PetscScalar *)malloc(xt_max_nnz * sizeof(PetscScalar));
310a501084fSBarry Smith   y                       = (PetscScalar *)malloc(yt_max_nnz * sizeof(PetscScalar));
311827bd09bSSatish Balay 
312827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
313827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
314db4deed7SKarl Rupp   for (dim = i = j = 0; i < m; i++) {
315827bd09bSSatish Balay     /* time to move to the next level? */
316d4af0d30SBarry Smith     while (i == segs[dim]) {
31708401ef6SPierre Jolivet       PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level");
318827bd09bSSatish Balay       stages[dim++] = i;
319827bd09bSSatish Balay       end += lnsep[dim];
320827bd09bSSatish Balay     }
321827bd09bSSatish Balay     stages[dim] = i;
322827bd09bSSatish Balay 
323827bd09bSSatish Balay     /* which column are we firing? */
324827bd09bSSatish Balay     /* i.e. set v_l */
325827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
326827bd09bSSatish Balay     (start < end) ? (col = fo[start]) : (col = INT_MAX);
3273ba16761SJacob Faibussowitsch     PetscCall(PCTFS_giop_hc(&col, &work, 1, op2, dim));
328827bd09bSSatish Balay 
329827bd09bSSatish Balay     /* shouldn't need this */
330db4deed7SKarl Rupp     if (col == INT_MAX) {
3319566063dSJacob Faibussowitsch       PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n"));
332827bd09bSSatish Balay       continue;
333827bd09bSSatish Balay     }
334827bd09bSSatish Balay 
335827bd09bSSatish Balay     /* do I own it? I should */
3363ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_zero(v, a_m));
3372fa5cd67SKarl Rupp     if (col == fo[start]) {
338827bd09bSSatish Balay       start++;
339ca8e9878SJed Brown       idx = PCTFS_ivec_linear_search(col, a_local2global, a_n);
340*966bd95aSPierre Jolivet       PetscCheck(idx != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!");
3412fa5cd67SKarl Rupp       v[idx] = 1.0;
3422fa5cd67SKarl Rupp       j++;
343db4deed7SKarl Rupp     } else {
344ca8e9878SJed Brown       idx = PCTFS_ivec_linear_search(col, a_local2global, a_m);
3452fa5cd67SKarl Rupp       if (idx != -1) v[idx] = 1.0;
346827bd09bSSatish Balay     }
347827bd09bSSatish Balay 
348827bd09bSSatish Balay     /* perform u = A.v_l */
3493ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_zero(u, n));
3503ba16761SJacob Faibussowitsch     PetscCall(do_matvec(xyt_handle->mvi, v, u));
351827bd09bSSatish Balay 
352827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
353827bd09bSSatish Balay     /* technically only need to zero out first i entries */
354827bd09bSSatish Balay     /* later turn this into an XYT_solve call ? */
3553ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_zero(uu, m));
356827bd09bSSatish Balay     y_ptr = y;
357827bd09bSSatish Balay     iptr  = ycol_indices;
358db4deed7SKarl Rupp     for (k = 0; k < i; k++) {
359827bd09bSSatish Balay       off = *iptr++;
360827bd09bSSatish Balay       len = *iptr++;
3619566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(len, &dlen));
362792fecdfSBarry Smith       PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, y_ptr, &i1));
363827bd09bSSatish Balay       y_ptr += len;
364827bd09bSSatish Balay     }
365827bd09bSSatish Balay 
366827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
3679566063dSJacob Faibussowitsch     PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages));
368827bd09bSSatish Balay 
369827bd09bSSatish Balay     /* z = X.uu */
3703ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_zero(z, n));
371827bd09bSSatish Balay     x_ptr = x;
372827bd09bSSatish Balay     iptr  = xcol_indices;
373db4deed7SKarl Rupp     for (k = 0; k < i; k++) {
374827bd09bSSatish Balay       off = *iptr++;
375827bd09bSSatish Balay       len = *iptr++;
3769566063dSJacob Faibussowitsch       PetscCall(PetscBLASIntCast(len, &dlen));
377792fecdfSBarry Smith       PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1));
378827bd09bSSatish Balay       x_ptr += len;
379827bd09bSSatish Balay     }
380827bd09bSSatish Balay 
381827bd09bSSatish Balay     /* compute v_l = v_l - z */
3823ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_zero(v + a_n, a_m - a_n));
3839566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &dlen));
384792fecdfSBarry Smith     PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1));
385827bd09bSSatish Balay 
386827bd09bSSatish Balay     /* compute u_l = A.v_l */
3873ba16761SJacob Faibussowitsch     if (a_n != a_m) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim));
3883ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_zero(u, n));
3893ba16761SJacob Faibussowitsch     PetscCall(do_matvec(xyt_handle->mvi, v, u));
390827bd09bSSatish Balay 
391827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
3929566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &dlen));
393792fecdfSBarry Smith     PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, u, &i1));
394827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
3953ba16761SJacob Faibussowitsch     PetscCall(PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim));
396827bd09bSSatish Balay 
3978f1a2a5eSBarry Smith     alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha);
398827bd09bSSatish Balay 
399827bd09bSSatish Balay     /* check for small alpha                             */
400827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
40163a3b9bcSJacob Faibussowitsch     PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha));
402827bd09bSSatish Balay 
403827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
4043ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_scale(v, 1.0 / alpha, n));
4053ba16761SJacob Faibussowitsch     PetscCall(PCTFS_rvec_scale(u, 1.0 / alpha, n));
406827bd09bSSatish Balay 
407827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
408827bd09bSSatish Balay     flag = 1;
409827bd09bSSatish Balay     off = len = 0;
410db4deed7SKarl Rupp     for (k = 0; k < n; k++) {
411db4deed7SKarl Rupp       if (v[k] != 0.0) {
412827bd09bSSatish Balay         len = k;
4139371c9d4SSatish Balay         if (flag) {
4149371c9d4SSatish Balay           off  = k;
4159371c9d4SSatish Balay           flag = 0;
4169371c9d4SSatish Balay         }
417827bd09bSSatish Balay       }
418827bd09bSSatish Balay     }
419827bd09bSSatish Balay 
420827bd09bSSatish Balay     len -= (off - 1);
421827bd09bSSatish Balay 
422db4deed7SKarl Rupp     if (len > 0) {
423db4deed7SKarl Rupp       if ((xt_nnz + len) > xt_max_nnz) {
4249566063dSJacob Faibussowitsch         PetscCall(PetscInfo(0, "increasing space for X by 2x!\n"));
425827bd09bSSatish Balay         xt_max_nnz *= 2;
426a501084fSBarry Smith         x_ptr = (PetscScalar *)malloc(xt_max_nnz * sizeof(PetscScalar));
4273ba16761SJacob Faibussowitsch         PetscCall(PCTFS_rvec_copy(x_ptr, x, xt_nnz));
428a501084fSBarry Smith         free(x);
429827bd09bSSatish Balay         x = x_ptr;
430827bd09bSSatish Balay         x_ptr += xt_nnz;
431827bd09bSSatish Balay       }
432827bd09bSSatish Balay       xt_nnz += len;
4333ba16761SJacob Faibussowitsch       PetscCall(PCTFS_rvec_copy(x_ptr, v + off, len));
434827bd09bSSatish Balay 
435827bd09bSSatish Balay       xcol_indices[2 * i] = off;
436827bd09bSSatish Balay       xcol_sz[i] = xcol_indices[2 * i + 1] = len;
437827bd09bSSatish Balay       xcol_vals[i]                         = x_ptr;
438db4deed7SKarl Rupp     } else {
439827bd09bSSatish Balay       xcol_indices[2 * i] = 0;
440827bd09bSSatish Balay       xcol_sz[i] = xcol_indices[2 * i + 1] = 0;
441827bd09bSSatish Balay       xcol_vals[i]                         = x_ptr;
442827bd09bSSatish Balay     }
443827bd09bSSatish Balay 
444827bd09bSSatish Balay     /* add newly generated column, u_l, to Y */
445827bd09bSSatish Balay     flag = 1;
446827bd09bSSatish Balay     off = len = 0;
447db4deed7SKarl Rupp     for (k = 0; k < n; k++) {
448db4deed7SKarl Rupp       if (u[k] != 0.0) {
449827bd09bSSatish Balay         len = k;
4509371c9d4SSatish Balay         if (flag) {
4519371c9d4SSatish Balay           off  = k;
4529371c9d4SSatish Balay           flag = 0;
4539371c9d4SSatish Balay         }
454827bd09bSSatish Balay       }
455827bd09bSSatish Balay     }
456827bd09bSSatish Balay 
457827bd09bSSatish Balay     len -= (off - 1);
458827bd09bSSatish Balay 
459db4deed7SKarl Rupp     if (len > 0) {
460db4deed7SKarl Rupp       if ((yt_nnz + len) > yt_max_nnz) {
4619566063dSJacob Faibussowitsch         PetscCall(PetscInfo(0, "increasing space for Y by 2x!\n"));
462827bd09bSSatish Balay         yt_max_nnz *= 2;
463a501084fSBarry Smith         y_ptr = (PetscScalar *)malloc(yt_max_nnz * sizeof(PetscScalar));
4643ba16761SJacob Faibussowitsch         PetscCall(PCTFS_rvec_copy(y_ptr, y, yt_nnz));
465a501084fSBarry Smith         free(y);
466827bd09bSSatish Balay         y = y_ptr;
467827bd09bSSatish Balay         y_ptr += yt_nnz;
468827bd09bSSatish Balay       }
469827bd09bSSatish Balay       yt_nnz += len;
4703ba16761SJacob Faibussowitsch       PetscCall(PCTFS_rvec_copy(y_ptr, u + off, len));
471827bd09bSSatish Balay 
472827bd09bSSatish Balay       ycol_indices[2 * i] = off;
473827bd09bSSatish Balay       ycol_sz[i] = ycol_indices[2 * i + 1] = len;
474827bd09bSSatish Balay       ycol_vals[i]                         = y_ptr;
475db4deed7SKarl Rupp     } else {
476827bd09bSSatish Balay       ycol_indices[2 * i] = 0;
477827bd09bSSatish Balay       ycol_sz[i] = ycol_indices[2 * i + 1] = 0;
478827bd09bSSatish Balay       ycol_vals[i]                         = y_ptr;
479827bd09bSSatish Balay     }
480827bd09bSSatish Balay   }
481827bd09bSSatish Balay 
482827bd09bSSatish Balay   /* close off stages for execution phase */
483db4deed7SKarl Rupp   while (dim != level) {
484827bd09bSSatish Balay     stages[dim++] = i;
48563a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level));
486827bd09bSSatish Balay   }
487827bd09bSSatish Balay   stages[dim] = i;
488827bd09bSSatish Balay 
489827bd09bSSatish Balay   xyt_handle->info->n            = xyt_handle->mvi->n;
490827bd09bSSatish Balay   xyt_handle->info->m            = m;
491827bd09bSSatish Balay   xyt_handle->info->nnz          = xt_nnz + yt_nnz;
492827bd09bSSatish Balay   xyt_handle->info->max_nnz      = xt_max_nnz + yt_max_nnz;
493827bd09bSSatish Balay   xyt_handle->info->msg_buf_sz   = stages[level] - stages[0];
494a501084fSBarry Smith   xyt_handle->info->solve_uu     = (PetscScalar *)malloc(m * sizeof(PetscScalar));
495a501084fSBarry Smith   xyt_handle->info->solve_w      = (PetscScalar *)malloc(m * sizeof(PetscScalar));
496827bd09bSSatish Balay   xyt_handle->info->x            = x;
497827bd09bSSatish Balay   xyt_handle->info->xcol_vals    = xcol_vals;
498827bd09bSSatish Balay   xyt_handle->info->xcol_sz      = xcol_sz;
499827bd09bSSatish Balay   xyt_handle->info->xcol_indices = xcol_indices;
500827bd09bSSatish Balay   xyt_handle->info->stages       = stages;
501827bd09bSSatish Balay   xyt_handle->info->y            = y;
502827bd09bSSatish Balay   xyt_handle->info->ycol_vals    = ycol_vals;
503827bd09bSSatish Balay   xyt_handle->info->ycol_sz      = ycol_sz;
504827bd09bSSatish Balay   xyt_handle->info->ycol_indices = ycol_indices;
505827bd09bSSatish Balay 
506a501084fSBarry Smith   free(segs);
507a501084fSBarry Smith   free(u);
508a501084fSBarry Smith   free(v);
509a501084fSBarry Smith   free(uu);
510a501084fSBarry Smith   free(z);
511a501084fSBarry Smith   free(w);
512827bd09bSSatish Balay 
5133ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
514827bd09bSSatish Balay }
515827bd09bSSatish Balay 
do_xyt_solve(xyt_ADT xyt_handle,PetscScalar * uc)516d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc)
517d71ae5a4SJacob Faibussowitsch {
51852f87cdaSBarry Smith   PetscInt     off, len, *iptr;
51952f87cdaSBarry Smith   PetscInt     level        = xyt_handle->level;
52052f87cdaSBarry Smith   PetscInt     n            = xyt_handle->info->n;
52152f87cdaSBarry Smith   PetscInt     m            = xyt_handle->info->m;
52252f87cdaSBarry Smith   PetscInt    *stages       = xyt_handle->info->stages;
52352f87cdaSBarry Smith   PetscInt    *xcol_indices = xyt_handle->info->xcol_indices;
52452f87cdaSBarry Smith   PetscInt    *ycol_indices = xyt_handle->info->ycol_indices;
525a501084fSBarry Smith   PetscScalar *x_ptr, *y_ptr, *uu_ptr;
526a501084fSBarry Smith   PetscScalar *solve_uu = xyt_handle->info->solve_uu;
527a501084fSBarry Smith   PetscScalar *solve_w  = xyt_handle->info->solve_w;
528a501084fSBarry Smith   PetscScalar *x        = xyt_handle->info->x;
529a501084fSBarry Smith   PetscScalar *y        = xyt_handle->info->y;
5304a0de3f6SBarry Smith   PetscBLASInt i1       = 1, dlen;
531827bd09bSSatish Balay 
5320924e98cSBarry Smith   PetscFunctionBegin;
533827bd09bSSatish Balay   uu_ptr = solve_uu;
5343ba16761SJacob Faibussowitsch   PetscCall(PCTFS_rvec_zero(uu_ptr, m));
535827bd09bSSatish Balay 
536827bd09bSSatish Balay   /* x  = X.Y^T.b */
537827bd09bSSatish Balay   /* uu = Y^T.b */
538947b95d8SBarry Smith   for (y_ptr = y, iptr = ycol_indices; *iptr != -1; y_ptr += len) {
539c5df96a5SBarry Smith     off = *iptr++;
540c5df96a5SBarry Smith     len = *iptr++;
5419566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(len, &dlen));
542792fecdfSBarry Smith     PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, y_ptr, &i1));
543827bd09bSSatish Balay   }
544827bd09bSSatish Balay 
545d5b43468SJose E. Roman   /* communication of beta */
546827bd09bSSatish Balay   uu_ptr = solve_uu;
5479566063dSJacob Faibussowitsch   if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages));
5483ba16761SJacob Faibussowitsch   PetscCall(PCTFS_rvec_zero(uc, n));
549827bd09bSSatish Balay 
550827bd09bSSatish Balay   /* x = X.uu */
551db4deed7SKarl Rupp   for (x_ptr = x, iptr = xcol_indices; *iptr != -1; x_ptr += len) {
552c5df96a5SBarry Smith     off = *iptr++;
553c5df96a5SBarry Smith     len = *iptr++;
5549566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(len, &dlen));
555792fecdfSBarry Smith     PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1));
556827bd09bSSatish Balay   }
5573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
558827bd09bSSatish Balay }
559827bd09bSSatish Balay 
check_handle(xyt_ADT xyt_handle)560d71ae5a4SJacob Faibussowitsch static PetscErrorCode check_handle(xyt_ADT xyt_handle)
561d71ae5a4SJacob Faibussowitsch {
56252f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX};
563827bd09bSSatish Balay 
5640924e98cSBarry Smith   PetscFunctionBegin;
5653274405dSPierre Jolivet   PetscCheck(xyt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", (void *)xyt_handle);
566827bd09bSSatish Balay 
567827bd09bSSatish Balay   vals[0] = vals[1] = xyt_handle->id;
5683ba16761SJacob Faibussowitsch   PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op));
56963a3b9bcSJacob 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);
5703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
571827bd09bSSatish Balay }
572827bd09bSSatish Balay 
det_separators(xyt_ADT xyt_handle)573d71ae5a4SJacob Faibussowitsch static PetscErrorCode det_separators(xyt_ADT xyt_handle)
574d71ae5a4SJacob Faibussowitsch {
57552f87cdaSBarry Smith   PetscInt     i, ct, id;
57652f87cdaSBarry Smith   PetscInt     mask, edge, *iptr;
57752f87cdaSBarry Smith   PetscInt    *dir, *used;
57852f87cdaSBarry Smith   PetscInt     sum[4], w[4];
579a501084fSBarry Smith   PetscScalar  rsum[4], rw[4];
58052f87cdaSBarry Smith   PetscInt     op[] = {GL_ADD, 0};
581a501084fSBarry Smith   PetscScalar *lhs, *rhs;
58252f87cdaSBarry Smith   PetscInt    *nsep, *lnsep, *fo, nfo = 0;
583ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle = xyt_handle->mvi->PCTFS_gs_handle;
58452f87cdaSBarry Smith   PetscInt    *local2global    = xyt_handle->mvi->local2global;
58552f87cdaSBarry Smith   PetscInt     n               = xyt_handle->mvi->n;
58652f87cdaSBarry Smith   PetscInt     m               = xyt_handle->mvi->m;
58752f87cdaSBarry Smith   PetscInt     level           = xyt_handle->level;
588ab824b78SBarry Smith   PetscInt     shared          = 0;
589827bd09bSSatish Balay 
5900924e98cSBarry Smith   PetscFunctionBegin;
59152f87cdaSBarry Smith   dir   = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
59252f87cdaSBarry Smith   nsep  = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
59352f87cdaSBarry Smith   lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
59452f87cdaSBarry Smith   fo    = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1));
59552f87cdaSBarry Smith   used  = (PetscInt *)malloc(sizeof(PetscInt) * n);
596827bd09bSSatish Balay 
5973ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_zero(dir, level + 1));
5983ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_zero(nsep, level + 1));
5993ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_zero(lnsep, level + 1));
6003ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_set(fo, -1, n + 1));
6013ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_zero(used, n));
602827bd09bSSatish Balay 
6038cda6cd7SBarry Smith   lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m);
6048cda6cd7SBarry Smith   rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m);
605827bd09bSSatish Balay 
606827bd09bSSatish Balay   /* determine the # of unique dof */
6073ba16761SJacob Faibussowitsch   PetscCall(PCTFS_rvec_zero(lhs, m));
6083ba16761SJacob Faibussowitsch   PetscCall(PCTFS_rvec_set(lhs, 1.0, n));
6093ba16761SJacob Faibussowitsch   PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level));
6109566063dSJacob Faibussowitsch   PetscCall(PetscInfo(0, "done first PCTFS_gs_gop_hc\n"));
6113ba16761SJacob Faibussowitsch   PetscCall(PCTFS_rvec_zero(rsum, 2));
612947b95d8SBarry Smith   for (i = 0; i < n; i++) {
6139371c9d4SSatish Balay     if (lhs[i] != 0.0) {
6149371c9d4SSatish Balay       rsum[0] += 1.0 / lhs[i];
6159371c9d4SSatish Balay       rsum[1] += lhs[i];
6169371c9d4SSatish Balay     }
6172fa5cd67SKarl Rupp     if (lhs[i] != 1.0) shared = 1;
618827bd09bSSatish Balay   }
619827bd09bSSatish Balay 
6203ba16761SJacob Faibussowitsch   PetscCall(PCTFS_grop_hc(rsum, rw, 2, op, level));
621827bd09bSSatish Balay   rsum[0] += 0.1;
622827bd09bSSatish Balay   rsum[1] += 0.1;
623827bd09bSSatish Balay 
62452f87cdaSBarry Smith   xyt_handle->info->n_global = xyt_handle->info->m_global = (PetscInt)rsum[0];
62552f87cdaSBarry Smith   xyt_handle->mvi->n_global = xyt_handle->mvi->m_global = (PetscInt)rsum[0];
626827bd09bSSatish Balay 
627827bd09bSSatish Balay   /* determine separator sets top down */
6280fdf79fbSJacob Faibussowitsch   PetscCheck(!shared, PETSC_COMM_SELF, PETSC_ERR_PLIB, "shared dof separator determination not ready ... see hmt!!!");
629827bd09bSSatish Balay   /* solution is to do as in the symmetric shared case but then */
630827bd09bSSatish Balay   /* pick the sub-hc with the most free dofs and do a mat-vec   */
631827bd09bSSatish Balay   /* and pick up the responses on the other sub-hc from the     */
632827bd09bSSatish Balay   /* initial separator set obtained from the symm. shared case  */
6334e619c20SJed Brown   /* [dead code deleted since it is unlikely to be completed] */
634db4deed7SKarl Rupp   for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) {
635827bd09bSSatish Balay     /* set rsh of hc, fire, and collect lhs responses */
6363ba16761SJacob Faibussowitsch     PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m));
6373ba16761SJacob Faibussowitsch     PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge));
638827bd09bSSatish Balay 
639827bd09bSSatish Balay     /* set lsh of hc, fire, and collect rhs responses */
6403ba16761SJacob Faibussowitsch     PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m));
6413ba16761SJacob Faibussowitsch     PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge));
642827bd09bSSatish Balay 
643827bd09bSSatish Balay     /* count number of dofs I own that have signal and not in sep set */
6443ba16761SJacob Faibussowitsch     PetscCall(PCTFS_ivec_zero(sum, 4));
6453ba16761SJacob Faibussowitsch     for (ct = i = 0; i < n; i++) {
646db4deed7SKarl Rupp       if (!used[i]) {
647827bd09bSSatish Balay         /* number of unmarked dofs on node */
648827bd09bSSatish Balay         ct++;
649827bd09bSSatish Balay         /* number of dofs to be marked on lhs hc */
6502fa5cd67SKarl Rupp         if ((id < mask) && (lhs[i] != 0.0)) sum[0]++;
651827bd09bSSatish Balay         /* number of dofs to be marked on rhs hc */
6522fa5cd67SKarl Rupp         if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++;
653827bd09bSSatish Balay       }
654827bd09bSSatish Balay     }
655827bd09bSSatish Balay 
656827bd09bSSatish Balay     /* for the non-symmetric case we need separators of width 2 */
657827bd09bSSatish Balay     /* so take both sides */
658827bd09bSSatish Balay     (id < mask) ? (sum[2] = ct) : (sum[3] = ct);
6593ba16761SJacob Faibussowitsch     PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge));
660827bd09bSSatish Balay 
661827bd09bSSatish Balay     ct = 0;
662db4deed7SKarl Rupp     if (id < mask) {
663827bd09bSSatish Balay       /* mark dofs I own that have signal and not in sep set */
664db4deed7SKarl Rupp       for (i = 0; i < n; i++) {
665db4deed7SKarl Rupp         if ((!used[i]) && (lhs[i] != 0.0)) {
6669371c9d4SSatish Balay           ct++;
6679371c9d4SSatish Balay           nfo++;
668827bd09bSSatish Balay           *--iptr = local2global[i];
669827bd09bSSatish Balay           used[i] = edge;
670827bd09bSSatish Balay         }
671827bd09bSSatish Balay       }
672827bd09bSSatish Balay       /* LSH hc summation of ct should be sum[0] */
673db4deed7SKarl Rupp     } else {
674827bd09bSSatish Balay       /* mark dofs I own that have signal and not in sep set */
675db4deed7SKarl Rupp       for (i = 0; i < n; i++) {
676db4deed7SKarl Rupp         if ((!used[i]) && (rhs[i] != 0.0)) {
6779371c9d4SSatish Balay           ct++;
6789371c9d4SSatish Balay           nfo++;
679827bd09bSSatish Balay           *--iptr = local2global[i];
680827bd09bSSatish Balay           used[i] = edge;
681827bd09bSSatish Balay         }
682827bd09bSSatish Balay       }
683827bd09bSSatish Balay       /* RSH hc summation of ct should be sum[1] */
684827bd09bSSatish Balay     }
685827bd09bSSatish Balay 
6863ba16761SJacob Faibussowitsch     if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
687827bd09bSSatish Balay     lnsep[edge] = ct;
688827bd09bSSatish Balay     nsep[edge]  = sum[0] + sum[1];
689827bd09bSSatish Balay     dir[edge]   = BOTH;
690827bd09bSSatish Balay 
691827bd09bSSatish Balay     /* LATER or we can recur on these to order seps at this level */
692827bd09bSSatish Balay     /* do we need full set of separators for this?                */
693827bd09bSSatish Balay 
694827bd09bSSatish Balay     /* fold rhs hc into lower */
6952fa5cd67SKarl Rupp     if (id >= mask) id -= mask;
696827bd09bSSatish Balay   }
697827bd09bSSatish Balay 
698827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
699db4deed7SKarl Rupp   for (ct = i = 0; i < n; i++) {
700db4deed7SKarl Rupp     if (!used[i]) {
7019371c9d4SSatish Balay       ct++;
7029371c9d4SSatish Balay       nfo++;
703827bd09bSSatish Balay       *--iptr = local2global[i];
704827bd09bSSatish Balay       used[i] = edge;
705827bd09bSSatish Balay     }
706827bd09bSSatish Balay   }
7073ba16761SJacob Faibussowitsch   if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
708827bd09bSSatish Balay   lnsep[edge] = ct;
709827bd09bSSatish Balay   nsep[edge]  = ct;
710827bd09bSSatish Balay   dir[edge]   = BOTH;
711827bd09bSSatish Balay 
712827bd09bSSatish Balay   xyt_handle->info->nsep  = nsep;
713827bd09bSSatish Balay   xyt_handle->info->lnsep = lnsep;
714827bd09bSSatish Balay   xyt_handle->info->fo    = fo;
715827bd09bSSatish Balay   xyt_handle->info->nfo   = nfo;
716827bd09bSSatish Balay 
717a501084fSBarry Smith   free(dir);
718a501084fSBarry Smith   free(lhs);
719a501084fSBarry Smith   free(rhs);
720a501084fSBarry Smith   free(used);
7213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
722827bd09bSSatish Balay }
723827bd09bSSatish Balay 
set_mvi(PetscInt * local2global,PetscInt n,PetscInt m,PetscErrorCode (* matvec)(mv_info *,PetscScalar *,PetscScalar *),void * grid_data)724d71ae5a4SJacob Faibussowitsch static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data)
725d71ae5a4SJacob Faibussowitsch {
726827bd09bSSatish Balay   mv_info *mvi;
727827bd09bSSatish Balay 
728a501084fSBarry Smith   mvi               = (mv_info *)malloc(sizeof(mv_info));
729827bd09bSSatish Balay   mvi->n            = n;
730827bd09bSSatish Balay   mvi->m            = m;
731827bd09bSSatish Balay   mvi->n_global     = -1;
732827bd09bSSatish Balay   mvi->m_global     = -1;
73352f87cdaSBarry Smith   mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt));
7342fa5cd67SKarl Rupp 
735ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global, local2global, m);
736827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
7375c8f6a95SKarl Rupp   mvi->matvec          = matvec;
738827bd09bSSatish Balay   mvi->grid_data       = grid_data;
739827bd09bSSatish Balay 
740827bd09bSSatish Balay   /* set xyt communication handle to perform restricted matvec */
741ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
742827bd09bSSatish Balay 
7434ad8454bSPierre Jolivet   return mvi;
744827bd09bSSatish Balay }
745827bd09bSSatish Balay 
do_matvec(mv_info * A,PetscScalar * v,PetscScalar * u)746d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
747d71ae5a4SJacob Faibussowitsch {
7480924e98cSBarry Smith   PetscFunctionBegin;
7493ba16761SJacob Faibussowitsch   PetscCall(A->matvec((mv_info *)A->grid_data, v, u));
7503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
751827bd09bSSatish Balay }
752