111cc89d2SBarry Smith /*
2827bd09bSSatish Balay Module Name: xxt
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 xxt_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 *col_sz, *col_indices;
29a501084fSBarry Smith PetscScalar **col_vals, *x, *solve_uu, *solve_w;
3052f87cdaSBarry Smith PetscInt nsolves;
31a501084fSBarry Smith PetscScalar tot_solve_time;
32827bd09bSSatish Balay } xxt_info;
33827bd09bSSatish Balay
34827bd09bSSatish Balay typedef struct matvec_info {
3552f87cdaSBarry Smith PetscInt n, m, n_global, m_global;
3652f87cdaSBarry Smith PetscInt *local2global;
37ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle;
38a501084fSBarry Smith PetscErrorCode (*matvec)(struct matvec_info *, PetscScalar *, PetscScalar *);
39827bd09bSSatish Balay void *grid_data;
40827bd09bSSatish Balay } mv_info;
41827bd09bSSatish Balay
42827bd09bSSatish Balay struct xxt_CDT {
4352f87cdaSBarry Smith PetscInt id;
4452f87cdaSBarry Smith PetscInt ns;
4552f87cdaSBarry Smith PetscInt level;
46827bd09bSSatish Balay xxt_info *info;
47827bd09bSSatish Balay mv_info *mvi;
48827bd09bSSatish Balay };
49827bd09bSSatish Balay
5052f87cdaSBarry Smith static PetscInt n_xxt = 0;
5152f87cdaSBarry Smith static PetscInt n_xxt_handles = 0;
52827bd09bSSatish Balay
53827bd09bSSatish Balay /* prototypes */
543fdc5746SBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *rhs);
553fdc5746SBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle);
563fdc5746SBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle);
573fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
5838dcbb09SBarry Smith static PetscErrorCode xxt_generate(xxt_ADT xxt_handle);
5938dcbb09SBarry Smith static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle);
605c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data);
61a501084fSBarry Smith
XXT_new(void)62d71ae5a4SJacob Faibussowitsch xxt_ADT XXT_new(void)
63d71ae5a4SJacob Faibussowitsch {
64827bd09bSSatish Balay xxt_ADT xxt_handle;
65827bd09bSSatish Balay
66827bd09bSSatish Balay /* rolling count on n_xxt ... pot. problem here */
67827bd09bSSatish Balay n_xxt_handles++;
68a501084fSBarry Smith xxt_handle = (xxt_ADT)malloc(sizeof(struct xxt_CDT));
69827bd09bSSatish Balay xxt_handle->id = ++n_xxt;
709371c9d4SSatish Balay xxt_handle->info = NULL;
719371c9d4SSatish Balay xxt_handle->mvi = NULL;
72827bd09bSSatish Balay
73*4ad8454bSPierre Jolivet return xxt_handle;
74827bd09bSSatish Balay }
75827bd09bSSatish Balay
XXT_factor(xxt_ADT xxt_handle,PetscInt * local2global,PetscInt n,PetscInt m,PetscErrorCode (* matvec)(void *,PetscScalar *,PetscScalar *),void * grid_data)7638dcbb09SBarry Smith PetscErrorCode XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt handle */
7752f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */
7852f87cdaSBarry Smith PetscInt n, /* local num rows */
7952f87cdaSBarry Smith PetscInt m, /* local num cols */
805c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void *, PetscScalar *, PetscScalar *), /* b_loc=A_local.x_loc */
811147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */
82827bd09bSSatish Balay {
833ba16761SJacob Faibussowitsch PetscFunctionBegin;
843ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init());
853ba16761SJacob Faibussowitsch PetscCall(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 */
1013ba16761SJacob Faibussowitsch PetscCall(det_separators(xxt_handle));
102827bd09bSSatish Balay
1033ba16761SJacob Faibussowitsch PetscCall(do_xxt_factor(xxt_handle));
1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
105827bd09bSSatish Balay }
106827bd09bSSatish Balay
XXT_solve(xxt_ADT xxt_handle,PetscScalar * x,PetscScalar * b)107d71ae5a4SJacob Faibussowitsch PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b)
108d71ae5a4SJacob Faibussowitsch {
1093ba16761SJacob Faibussowitsch PetscFunctionBegin;
1103ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init());
1113ba16761SJacob Faibussowitsch PetscCall(check_handle(xxt_handle));
112827bd09bSSatish Balay
113827bd09bSSatish Balay /* need to copy b into x? */
1143ba16761SJacob Faibussowitsch if (b) PetscCall(PCTFS_rvec_copy(x, b, xxt_handle->mvi->n));
1153ba16761SJacob Faibussowitsch PetscCall(do_xxt_solve(xxt_handle, x));
1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
117827bd09bSSatish Balay }
118827bd09bSSatish Balay
XXT_free(xxt_ADT xxt_handle)1193ba16761SJacob Faibussowitsch PetscErrorCode XXT_free(xxt_ADT xxt_handle)
120d71ae5a4SJacob Faibussowitsch {
1213ba16761SJacob Faibussowitsch PetscFunctionBegin;
1223ba16761SJacob Faibussowitsch PetscCall(PCTFS_comm_init());
1233ba16761SJacob Faibussowitsch PetscCall(check_handle(xxt_handle));
124827bd09bSSatish Balay n_xxt_handles--;
125827bd09bSSatish Balay
126a501084fSBarry Smith free(xxt_handle->info->nsep);
127a501084fSBarry Smith free(xxt_handle->info->lnsep);
128a501084fSBarry Smith free(xxt_handle->info->fo);
129a501084fSBarry Smith free(xxt_handle->info->stages);
130a501084fSBarry Smith free(xxt_handle->info->solve_uu);
131a501084fSBarry Smith free(xxt_handle->info->solve_w);
132a501084fSBarry Smith free(xxt_handle->info->x);
133a501084fSBarry Smith free(xxt_handle->info->col_vals);
134a501084fSBarry Smith free(xxt_handle->info->col_sz);
135a501084fSBarry Smith free(xxt_handle->info->col_indices);
136a501084fSBarry Smith free(xxt_handle->info);
137a501084fSBarry Smith free(xxt_handle->mvi->local2global);
1383ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle));
139a501084fSBarry Smith free(xxt_handle->mvi);
140a501084fSBarry Smith free(xxt_handle);
141827bd09bSSatish Balay
142827bd09bSSatish Balay /* if the check fails we nuke */
143a501084fSBarry Smith /* if NULL pointer passed to free we nuke */
144827bd09bSSatish Balay /* if the calls to free fail that's not my problem */
1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
146827bd09bSSatish Balay }
147827bd09bSSatish Balay
14811cc89d2SBarry Smith /*
149827bd09bSSatish Balay
150827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
151827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
152827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid
153827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1)
154827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1)
155827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using
156ca8e9878SJed Brown PCTFS_gs_init/gop).
157827bd09bSSatish Balay
158827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
159827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
16011cc89d2SBarry Smith */
do_xxt_factor(xxt_ADT xxt_handle)161d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle)
162d71ae5a4SJacob Faibussowitsch {
1637b1ae94cSBarry Smith return xxt_generate(xxt_handle);
164827bd09bSSatish Balay }
165827bd09bSSatish Balay
xxt_generate(xxt_ADT xxt_handle)166d71ae5a4SJacob Faibussowitsch static PetscErrorCode xxt_generate(xxt_ADT xxt_handle)
167d71ae5a4SJacob Faibussowitsch {
16852f87cdaSBarry Smith PetscInt i, j, k, idex;
16952f87cdaSBarry Smith PetscInt dim, col;
170a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
17152f87cdaSBarry Smith PetscInt *segs;
17252f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0};
17352f87cdaSBarry Smith PetscInt off, len;
174a501084fSBarry Smith PetscScalar *x_ptr;
17552f87cdaSBarry Smith PetscInt *iptr, flag;
17652f87cdaSBarry Smith PetscInt start = 0, end, work;
17752f87cdaSBarry Smith PetscInt op2[] = {GL_MIN, 0};
178ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle;
17952f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo;
18052f87cdaSBarry Smith PetscInt a_n = xxt_handle->mvi->n;
18152f87cdaSBarry Smith PetscInt a_m = xxt_handle->mvi->m;
18252f87cdaSBarry Smith PetscInt *a_local2global = xxt_handle->mvi->local2global;
18352f87cdaSBarry Smith PetscInt level;
18452f87cdaSBarry Smith PetscInt xxt_nnz = 0, xxt_max_nnz = 0;
18552f87cdaSBarry Smith PetscInt n, m;
18652f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages;
187a501084fSBarry Smith PetscScalar **col_vals, *x;
18852f87cdaSBarry Smith PetscInt n_global;
18971a0148aSBarry Smith PetscBLASInt i1 = 1, dlen;
190a501084fSBarry Smith PetscScalar dm1 = -1.0;
191827bd09bSSatish Balay
192827bd09bSSatish Balay n = xxt_handle->mvi->n;
193827bd09bSSatish Balay nsep = xxt_handle->info->nsep;
194827bd09bSSatish Balay lnsep = xxt_handle->info->lnsep;
195827bd09bSSatish Balay fo = xxt_handle->info->fo;
196827bd09bSSatish Balay end = lnsep[0];
197827bd09bSSatish Balay level = xxt_handle->level;
198ca8e9878SJed Brown PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
199827bd09bSSatish Balay
200827bd09bSSatish Balay /* is there a null space? */
201827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */
2022fa5cd67SKarl Rupp for (i = 0, j = 0; i <= level; i++) j += nsep[i];
203827bd09bSSatish Balay
204827bd09bSSatish Balay m = j - xxt_handle->ns;
20548a46eb9SPierre 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));
206827bd09bSSatish Balay
207827bd09bSSatish Balay /* get and initialize storage for x local */
208827bd09bSSatish Balay /* note that x local is nxm and stored by columns */
20952f87cdaSBarry Smith col_sz = (PetscInt *)malloc(m * sizeof(PetscInt));
21052f87cdaSBarry Smith col_indices = (PetscInt *)malloc((2 * m + 1) * sizeof(PetscInt));
211a501084fSBarry Smith col_vals = (PetscScalar **)malloc(m * sizeof(PetscScalar *));
212db4deed7SKarl Rupp for (i = j = 0; i < m; i++, j += 2) {
213827bd09bSSatish Balay col_indices[j] = col_indices[j + 1] = col_sz[i] = -1;
214827bd09bSSatish Balay col_vals[i] = NULL;
215827bd09bSSatish Balay }
216827bd09bSSatish Balay col_indices[j] = -1;
217827bd09bSSatish Balay
218827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */
219827bd09bSSatish Balay /* this looks like nsep[]=segments */
22052f87cdaSBarry Smith stages = (PetscInt *)malloc((level + 1) * sizeof(PetscInt));
22152f87cdaSBarry Smith segs = (PetscInt *)malloc((level + 1) * sizeof(PetscInt));
2223ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(stages, level + 1));
223ca8e9878SJed Brown PCTFS_ivec_copy(segs, nsep, level + 1);
2242fa5cd67SKarl Rupp for (i = 0; i < level; i++) segs[i + 1] += segs[i];
225827bd09bSSatish Balay stages[0] = segs[0];
226827bd09bSSatish Balay
227827bd09bSSatish Balay /* temporary vectors */
228a501084fSBarry Smith u = (PetscScalar *)malloc(n * sizeof(PetscScalar));
229a501084fSBarry Smith z = (PetscScalar *)malloc(n * sizeof(PetscScalar));
230a501084fSBarry Smith v = (PetscScalar *)malloc(a_m * sizeof(PetscScalar));
231a501084fSBarry Smith uu = (PetscScalar *)malloc(m * sizeof(PetscScalar));
232a501084fSBarry Smith w = (PetscScalar *)malloc(m * sizeof(PetscScalar));
233827bd09bSSatish Balay
234827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */
2352fa5cd67SKarl Rupp for (i = 1, j = 0; i <= level; i++) j += nsep[i];
236827bd09bSSatish Balay
237827bd09bSSatish Balay /* storage for sparse x values */
238827bd09bSSatish Balay n_global = xxt_handle->info->n_global;
23985ec1a3cSBarry Smith xxt_max_nnz = (PetscInt)(2.5 * PetscPowReal(1.0 * n_global, 1.6667) + j * n / 2) / PCTFS_num_nodes;
240a501084fSBarry Smith x = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar));
241827bd09bSSatish Balay xxt_nnz = 0;
242827bd09bSSatish Balay
243827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */
244827bd09bSSatish Balay /* time to make the donuts - generate X factor */
2452fa5cd67SKarl Rupp for (dim = i = j = 0; i < m; i++) {
246827bd09bSSatish Balay /* time to move to the next level? */
247d4af0d30SBarry Smith while (i == segs[dim]) {
24808401ef6SPierre Jolivet PetscCheck(dim != level, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim about to exceed level");
249827bd09bSSatish Balay stages[dim++] = i;
250827bd09bSSatish Balay end += lnsep[dim];
251827bd09bSSatish Balay }
252827bd09bSSatish Balay stages[dim] = i;
253827bd09bSSatish Balay
254827bd09bSSatish Balay /* which column are we firing? */
255827bd09bSSatish Balay /* i.e. set v_l */
256827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */
257827bd09bSSatish Balay (start < end) ? (col = fo[start]) : (col = INT_MAX);
2583ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(&col, &work, 1, op2, dim));
259827bd09bSSatish Balay
260827bd09bSSatish Balay /* shouldn't need this */
261db4deed7SKarl Rupp if (col == INT_MAX) {
2629566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "hey ... col==INT_MAX??\n"));
263827bd09bSSatish Balay continue;
264827bd09bSSatish Balay }
265827bd09bSSatish Balay
266827bd09bSSatish Balay /* do I own it? I should */
2673ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(v, a_m));
268db4deed7SKarl Rupp if (col == fo[start]) {
269827bd09bSSatish Balay start++;
270ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_n);
2710fdf79fbSJacob Faibussowitsch PetscCheck(idex != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "NOT FOUND!");
2729371c9d4SSatish Balay v[idex] = 1.0;
2739371c9d4SSatish Balay j++;
274db4deed7SKarl Rupp } else {
275ca8e9878SJed Brown idex = PCTFS_ivec_linear_search(col, a_local2global, a_m);
2762fa5cd67SKarl Rupp if (idex != -1) v[idex] = 1.0;
277827bd09bSSatish Balay }
278827bd09bSSatish Balay
279827bd09bSSatish Balay /* perform u = A.v_l */
2803ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(u, n));
2813ba16761SJacob Faibussowitsch PetscCall(do_matvec(xxt_handle->mvi, v, u));
282827bd09bSSatish Balay
283827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */
284827bd09bSSatish Balay /* technically only need to zero out first i entries */
285827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */
2863ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uu, m));
287827bd09bSSatish Balay x_ptr = x;
288827bd09bSSatish Balay iptr = col_indices;
289db4deed7SKarl Rupp for (k = 0; k < i; k++) {
290827bd09bSSatish Balay off = *iptr++;
291827bd09bSSatish Balay len = *iptr++;
2929566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen));
293792fecdfSBarry Smith PetscCallBLAS("BLASdot", uu[k] = BLASdot_(&dlen, u + off, &i1, x_ptr, &i1));
294827bd09bSSatish Balay x_ptr += len;
295827bd09bSSatish Balay }
296827bd09bSSatish Balay
297827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */
2989566063dSJacob Faibussowitsch PetscCall(PCTFS_ssgl_radd(uu, w, dim, stages));
299827bd09bSSatish Balay
300827bd09bSSatish Balay /* z = X.uu */
3013ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(z, n));
302827bd09bSSatish Balay x_ptr = x;
303827bd09bSSatish Balay iptr = col_indices;
304db4deed7SKarl Rupp for (k = 0; k < i; k++) {
305827bd09bSSatish Balay off = *iptr++;
306827bd09bSSatish Balay len = *iptr++;
3079566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen));
308792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &uu[k], x_ptr, &i1, z + off, &i1));
309827bd09bSSatish Balay x_ptr += len;
310827bd09bSSatish Balay }
311827bd09bSSatish Balay
312827bd09bSSatish Balay /* compute v_l = v_l - z */
3133ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(v + a_n, a_m - a_n));
3149566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen));
315792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, &dm1, z, &i1, v, &i1));
316827bd09bSSatish Balay
317827bd09bSSatish Balay /* compute u_l = A.v_l */
3183ba16761SJacob Faibussowitsch if (a_n != a_m) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, v, "+\0", dim));
3193ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(u, n));
3203ba16761SJacob Faibussowitsch PetscCall(do_matvec(xxt_handle->mvi, v, u));
321827bd09bSSatish Balay
322827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
3239566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &dlen));
324792fecdfSBarry Smith PetscCallBLAS("BLASdot", alpha = BLASdot_(&dlen, u, &i1, v, &i1));
325827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
3263ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim));
327827bd09bSSatish Balay
3288f1a2a5eSBarry Smith alpha = (PetscScalar)PetscSqrtReal((PetscReal)alpha);
329827bd09bSSatish Balay
330827bd09bSSatish Balay /* check for small alpha */
331827bd09bSSatish Balay /* LATER use this to detect and determine null space */
33263a3b9bcSJacob Faibussowitsch PetscCheck(PetscAbsScalar(alpha) >= 1.0e-14, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bad alpha! %g", (double)PetscAbsScalar(alpha));
333827bd09bSSatish Balay
334827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */
3353ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_scale(v, 1.0 / alpha, n));
336827bd09bSSatish Balay
337827bd09bSSatish Balay /* add newly generated column, v_l, to X */
338827bd09bSSatish Balay flag = 1;
339827bd09bSSatish Balay off = len = 0;
340db4deed7SKarl Rupp for (k = 0; k < n; k++) {
341db4deed7SKarl Rupp if (v[k] != 0.0) {
342827bd09bSSatish Balay len = k;
3439371c9d4SSatish Balay if (flag) {
3449371c9d4SSatish Balay off = k;
3459371c9d4SSatish Balay flag = 0;
3469371c9d4SSatish Balay }
347827bd09bSSatish Balay }
348827bd09bSSatish Balay }
349827bd09bSSatish Balay
350827bd09bSSatish Balay len -= (off - 1);
351827bd09bSSatish Balay
352db4deed7SKarl Rupp if (len > 0) {
353db4deed7SKarl Rupp if ((xxt_nnz + len) > xxt_max_nnz) {
3549566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "increasing space for X by 2x!\n"));
355827bd09bSSatish Balay xxt_max_nnz *= 2;
356a501084fSBarry Smith x_ptr = (PetscScalar *)malloc(xxt_max_nnz * sizeof(PetscScalar));
3573ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_copy(x_ptr, x, xxt_nnz));
358a501084fSBarry Smith free(x);
359827bd09bSSatish Balay x = x_ptr;
360827bd09bSSatish Balay x_ptr += xxt_nnz;
361827bd09bSSatish Balay }
362827bd09bSSatish Balay xxt_nnz += len;
3633ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_copy(x_ptr, v + off, len));
364827bd09bSSatish Balay
365827bd09bSSatish Balay col_indices[2 * i] = off;
366827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = len;
367827bd09bSSatish Balay col_vals[i] = x_ptr;
3689371c9d4SSatish Balay } else {
369827bd09bSSatish Balay col_indices[2 * i] = 0;
370827bd09bSSatish Balay col_sz[i] = col_indices[2 * i + 1] = 0;
371827bd09bSSatish Balay col_vals[i] = x_ptr;
372827bd09bSSatish Balay }
373827bd09bSSatish Balay }
374827bd09bSSatish Balay
375827bd09bSSatish Balay /* close off stages for execution phase */
376db4deed7SKarl Rupp while (dim != level) {
377827bd09bSSatish Balay stages[dim++] = i;
37863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(0, "disconnected!!! dim(%" PetscInt_FMT ")!=level(%" PetscInt_FMT ")\n", dim, level));
379827bd09bSSatish Balay }
380827bd09bSSatish Balay stages[dim] = i;
381827bd09bSSatish Balay
382827bd09bSSatish Balay xxt_handle->info->n = xxt_handle->mvi->n;
383827bd09bSSatish Balay xxt_handle->info->m = m;
384827bd09bSSatish Balay xxt_handle->info->nnz = xxt_nnz;
385827bd09bSSatish Balay xxt_handle->info->max_nnz = xxt_max_nnz;
386827bd09bSSatish Balay xxt_handle->info->msg_buf_sz = stages[level] - stages[0];
387a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar *)malloc(m * sizeof(PetscScalar));
388a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar *)malloc(m * sizeof(PetscScalar));
389827bd09bSSatish Balay xxt_handle->info->x = x;
390827bd09bSSatish Balay xxt_handle->info->col_vals = col_vals;
391827bd09bSSatish Balay xxt_handle->info->col_sz = col_sz;
392827bd09bSSatish Balay xxt_handle->info->col_indices = col_indices;
393827bd09bSSatish Balay xxt_handle->info->stages = stages;
394827bd09bSSatish Balay xxt_handle->info->nsolves = 0;
395827bd09bSSatish Balay xxt_handle->info->tot_solve_time = 0.0;
396827bd09bSSatish Balay
397a501084fSBarry Smith free(segs);
398a501084fSBarry Smith free(u);
399a501084fSBarry Smith free(v);
400a501084fSBarry Smith free(uu);
401a501084fSBarry Smith free(z);
402a501084fSBarry Smith free(w);
403827bd09bSSatish Balay
4043ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
405827bd09bSSatish Balay }
406827bd09bSSatish Balay
do_xxt_solve(xxt_ADT xxt_handle,PetscScalar * uc)407d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc)
408d71ae5a4SJacob Faibussowitsch {
40952f87cdaSBarry Smith PetscInt off, len, *iptr;
41052f87cdaSBarry Smith PetscInt level = xxt_handle->level;
41152f87cdaSBarry Smith PetscInt n = xxt_handle->info->n;
41252f87cdaSBarry Smith PetscInt m = xxt_handle->info->m;
41352f87cdaSBarry Smith PetscInt *stages = xxt_handle->info->stages;
41452f87cdaSBarry Smith PetscInt *col_indices = xxt_handle->info->col_indices;
415a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr;
416a501084fSBarry Smith PetscScalar *solve_uu = xxt_handle->info->solve_uu;
417a501084fSBarry Smith PetscScalar *solve_w = xxt_handle->info->solve_w;
418a501084fSBarry Smith PetscScalar *x = xxt_handle->info->x;
41971a0148aSBarry Smith PetscBLASInt i1 = 1, dlen;
420827bd09bSSatish Balay
4210924e98cSBarry Smith PetscFunctionBegin;
422827bd09bSSatish Balay uu_ptr = solve_uu;
4233ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uu_ptr, m));
424827bd09bSSatish Balay
425827bd09bSSatish Balay /* x = X.Y^T.b */
426827bd09bSSatish Balay /* uu = Y^T.b */
427db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) {
428c5df96a5SBarry Smith off = *iptr++;
429c5df96a5SBarry Smith len = *iptr++;
4309566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen));
431792fecdfSBarry Smith PetscCallBLAS("BLASdot", *uu_ptr++ = BLASdot_(&dlen, uc + off, &i1, x_ptr, &i1));
432827bd09bSSatish Balay }
433827bd09bSSatish Balay
434d5b43468SJose E. Roman /* communication of beta */
435827bd09bSSatish Balay uu_ptr = solve_uu;
4369566063dSJacob Faibussowitsch if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages));
437827bd09bSSatish Balay
4383ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(uc, n));
439827bd09bSSatish Balay
440827bd09bSSatish Balay /* x = X.uu */
441db4deed7SKarl Rupp for (x_ptr = x, iptr = col_indices; *iptr != -1; x_ptr += len) {
442c5df96a5SBarry Smith off = *iptr++;
443c5df96a5SBarry Smith len = *iptr++;
4449566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len, &dlen));
445792fecdfSBarry Smith PetscCallBLAS("BLASaxpy", BLASaxpy_(&dlen, uu_ptr++, x_ptr, &i1, uc + off, &i1));
446827bd09bSSatish Balay }
4473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
448827bd09bSSatish Balay }
449827bd09bSSatish Balay
check_handle(xxt_ADT xxt_handle)450d71ae5a4SJacob Faibussowitsch static PetscErrorCode check_handle(xxt_ADT xxt_handle)
451d71ae5a4SJacob Faibussowitsch {
45252f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM, GL_MIN, GL_MAX};
453827bd09bSSatish Balay
4540924e98cSBarry Smith PetscFunctionBegin;
4557ee088dcSPierre Jolivet PetscCheck(xxt_handle, PETSC_COMM_SELF, PETSC_ERR_PLIB, "check_handle() :: bad handle :: NULL %p", (void *)xxt_handle);
456827bd09bSSatish Balay
457827bd09bSSatish Balay vals[0] = vals[1] = xxt_handle->id;
4583ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop(vals, work, PETSC_STATIC_ARRAY_LENGTH(op) - 1, op));
45963a3b9bcSJacob 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);
4603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
461827bd09bSSatish Balay }
462827bd09bSSatish Balay
det_separators(xxt_ADT xxt_handle)463d71ae5a4SJacob Faibussowitsch static PetscErrorCode det_separators(xxt_ADT xxt_handle)
464d71ae5a4SJacob Faibussowitsch {
46552f87cdaSBarry Smith PetscInt i, ct, id;
46652f87cdaSBarry Smith PetscInt mask, edge, *iptr;
46752f87cdaSBarry Smith PetscInt *dir, *used;
46852f87cdaSBarry Smith PetscInt sum[4], w[4];
469a501084fSBarry Smith PetscScalar rsum[4], rw[4];
47052f87cdaSBarry Smith PetscInt op[] = {GL_ADD, 0};
471a501084fSBarry Smith PetscScalar *lhs, *rhs;
47252f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo = 0;
473ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
47452f87cdaSBarry Smith PetscInt *local2global = xxt_handle->mvi->local2global;
47552f87cdaSBarry Smith PetscInt n = xxt_handle->mvi->n;
47652f87cdaSBarry Smith PetscInt m = xxt_handle->mvi->m;
47752f87cdaSBarry Smith PetscInt level = xxt_handle->level;
478ab824b78SBarry Smith PetscInt shared = 0;
479827bd09bSSatish Balay
4800924e98cSBarry Smith PetscFunctionBegin;
48152f87cdaSBarry Smith dir = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
48252f87cdaSBarry Smith nsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
48352f87cdaSBarry Smith lnsep = (PetscInt *)malloc(sizeof(PetscInt) * (level + 1));
48452f87cdaSBarry Smith fo = (PetscInt *)malloc(sizeof(PetscInt) * (n + 1));
48552f87cdaSBarry Smith used = (PetscInt *)malloc(sizeof(PetscInt) * n);
486827bd09bSSatish Balay
4873ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(dir, level + 1));
4883ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(nsep, level + 1));
4893ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(lnsep, level + 1));
4903ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_set(fo, -1, n + 1));
4913ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(used, n));
492827bd09bSSatish Balay
4938cda6cd7SBarry Smith lhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m);
4948cda6cd7SBarry Smith rhs = (PetscScalar *)malloc(sizeof(PetscScalar) * m);
495827bd09bSSatish Balay
496827bd09bSSatish Balay /* determine the # of unique dof */
4973ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(lhs, m));
4983ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_set(lhs, 1.0, n));
4993ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", level));
5003ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(rsum, 2));
5013d3eaba7SBarry Smith for (i = 0; i < n; i++) {
5022fa5cd67SKarl Rupp if (lhs[i] != 0.0) {
5039371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i];
5049371c9d4SSatish Balay rsum[1] += lhs[i];
5052fa5cd67SKarl Rupp }
506827bd09bSSatish Balay }
5073ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(rsum, rw, 2, op, level));
508827bd09bSSatish Balay rsum[0] += 0.1;
509827bd09bSSatish Balay rsum[1] += 0.1;
510827bd09bSSatish Balay
51177b4d14cSPeter Brune if (PetscAbsScalar(rsum[0] - rsum[1]) > EPS) shared = 1;
512827bd09bSSatish Balay
51352f87cdaSBarry Smith xxt_handle->info->n_global = xxt_handle->info->m_global = (PetscInt)rsum[0];
51452f87cdaSBarry Smith xxt_handle->mvi->n_global = xxt_handle->mvi->m_global = (PetscInt)rsum[0];
515827bd09bSSatish Balay
516827bd09bSSatish Balay /* determine separator sets top down */
5172fa5cd67SKarl Rupp if (shared) {
518db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) {
519827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */
5203ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m));
5213ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge));
522827bd09bSSatish Balay
523827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */
5243ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m));
5253ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge));
526827bd09bSSatish Balay
527db4deed7SKarl Rupp for (i = 0; i < n; i++) {
528db4deed7SKarl Rupp if (id < mask) {
5292fa5cd67SKarl Rupp if (lhs[i] != 0.0) lhs[i] = 1.0;
530827bd09bSSatish Balay }
531db4deed7SKarl Rupp if (id >= mask) {
5322fa5cd67SKarl Rupp if (rhs[i] != 0.0) rhs[i] = 1.0;
533827bd09bSSatish Balay }
534827bd09bSSatish Balay }
535827bd09bSSatish Balay
5363ba16761SJacob Faibussowitsch if (id < mask) PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge - 1));
5373ba16761SJacob Faibussowitsch else PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge - 1));
538827bd09bSSatish Balay
539827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */
5403ba16761SJacob Faibussowitsch PetscCall(PCTFS_rvec_zero(rsum, 4));
5413ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(sum, 4));
5423ba16761SJacob Faibussowitsch for (ct = i = 0; i < n; i++) {
543db4deed7SKarl Rupp if (!used[i]) {
544827bd09bSSatish Balay /* number of unmarked dofs on node */
545827bd09bSSatish Balay ct++;
546827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */
547db4deed7SKarl Rupp if (id < mask) {
5489371c9d4SSatish Balay if (lhs[i] != 0.0) {
5499371c9d4SSatish Balay sum[0]++;
5509371c9d4SSatish Balay rsum[0] += 1.0 / lhs[i];
5519371c9d4SSatish Balay }
552827bd09bSSatish Balay }
553827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */
554db4deed7SKarl Rupp if (id >= mask) {
5559371c9d4SSatish Balay if (rhs[i] != 0.0) {
5569371c9d4SSatish Balay sum[1]++;
5579371c9d4SSatish Balay rsum[1] += 1.0 / rhs[i];
5589371c9d4SSatish Balay }
559827bd09bSSatish Balay }
560827bd09bSSatish Balay }
561827bd09bSSatish Balay }
562827bd09bSSatish Balay
563827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */
564827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct);
565827bd09bSSatish Balay (id < mask) ? (rsum[2] = ct) : (rsum[3] = ct);
5663ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge));
5673ba16761SJacob Faibussowitsch PetscCall(PCTFS_grop_hc(rsum, rw, 4, op, edge));
5689371c9d4SSatish Balay rsum[0] += 0.1;
5699371c9d4SSatish Balay rsum[1] += 0.1;
5709371c9d4SSatish Balay rsum[2] += 0.1;
5719371c9d4SSatish Balay rsum[3] += 0.1;
572827bd09bSSatish Balay
573db4deed7SKarl Rupp if (id < mask) {
574827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */
575db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) {
576db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) {
5779371c9d4SSatish Balay ct++;
5789371c9d4SSatish Balay nfo++;
579827bd09bSSatish Balay
58008401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n");
581827bd09bSSatish Balay
582827bd09bSSatish Balay *--iptr = local2global[i];
583827bd09bSSatish Balay used[i] = edge;
584827bd09bSSatish Balay }
585827bd09bSSatish Balay }
5863ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
587827bd09bSSatish Balay
588827bd09bSSatish Balay lnsep[edge] = ct;
58952f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[0];
590827bd09bSSatish Balay dir[edge] = LEFT;
591827bd09bSSatish Balay }
592827bd09bSSatish Balay
593db4deed7SKarl Rupp if (id >= mask) {
594827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */
595db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) {
596db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) {
5979371c9d4SSatish Balay ct++;
5989371c9d4SSatish Balay nfo++;
599827bd09bSSatish Balay
60008401ef6SPierre Jolivet PetscCheck(nfo <= n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nfo about to exceed n");
601827bd09bSSatish Balay
602827bd09bSSatish Balay *--iptr = local2global[i];
603827bd09bSSatish Balay used[i] = edge;
604827bd09bSSatish Balay }
605827bd09bSSatish Balay }
6063ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
607827bd09bSSatish Balay
608827bd09bSSatish Balay lnsep[edge] = ct;
60952f87cdaSBarry Smith nsep[edge] = (PetscInt)rsum[1];
610827bd09bSSatish Balay dir[edge] = RIGHT;
611827bd09bSSatish Balay }
612827bd09bSSatish Balay
613827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */
614827bd09bSSatish Balay /* do we need full set of separators for this? */
615827bd09bSSatish Balay
616827bd09bSSatish Balay /* fold rhs hc into lower */
6172fa5cd67SKarl Rupp if (id >= mask) id -= mask;
618827bd09bSSatish Balay }
619db4deed7SKarl Rupp } else {
620db4deed7SKarl Rupp for (iptr = fo + n, id = PCTFS_my_id, mask = PCTFS_num_nodes >> 1, edge = level; edge > 0; edge--, mask >>= 1) {
621827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */
6223ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_zero(lhs, m) : PCTFS_rvec_set(lhs, 1.0, m));
6233ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, lhs, "+\0", edge));
624827bd09bSSatish Balay
625827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */
6263ba16761SJacob Faibussowitsch PetscCall((id < mask) ? PCTFS_rvec_set(rhs, 1.0, m) : PCTFS_rvec_zero(rhs, m));
6273ba16761SJacob Faibussowitsch PetscCall(PCTFS_gs_gop_hc(PCTFS_gs_handle, rhs, "+\0", edge));
628827bd09bSSatish Balay
629827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */
6303ba16761SJacob Faibussowitsch PetscCall(PCTFS_ivec_zero(sum, 4));
6313ba16761SJacob Faibussowitsch for (ct = i = 0; i < n; i++) {
632db4deed7SKarl Rupp if (!used[i]) {
633827bd09bSSatish Balay /* number of unmarked dofs on node */
634827bd09bSSatish Balay ct++;
635827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */
6362fa5cd67SKarl Rupp if ((id < mask) && (lhs[i] != 0.0)) sum[0]++;
637827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */
6382fa5cd67SKarl Rupp if ((id >= mask) && (rhs[i] != 0.0)) sum[1]++;
639827bd09bSSatish Balay }
640827bd09bSSatish Balay }
641827bd09bSSatish Balay
642827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */
643827bd09bSSatish Balay (id < mask) ? (sum[2] = ct) : (sum[3] = ct);
6443ba16761SJacob Faibussowitsch PetscCall(PCTFS_giop_hc(sum, w, 4, op, edge));
645827bd09bSSatish Balay
646827bd09bSSatish Balay /* lhs hc wins */
647db4deed7SKarl Rupp if (sum[2] >= sum[3]) {
648db4deed7SKarl Rupp if (id < mask) {
649827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */
650db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) {
651db4deed7SKarl Rupp if ((!used[i]) && (lhs[i] != 0.0)) {
6529371c9d4SSatish Balay ct++;
6539371c9d4SSatish Balay nfo++;
654827bd09bSSatish Balay *--iptr = local2global[i];
655827bd09bSSatish Balay used[i] = edge;
656827bd09bSSatish Balay }
657827bd09bSSatish Balay }
6583ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
659827bd09bSSatish Balay lnsep[edge] = ct;
660827bd09bSSatish Balay }
661827bd09bSSatish Balay nsep[edge] = sum[0];
662827bd09bSSatish Balay dir[edge] = LEFT;
663db4deed7SKarl Rupp } else { /* rhs hc wins */
664db4deed7SKarl Rupp if (id >= mask) {
665827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */
666db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) {
667db4deed7SKarl Rupp if ((!used[i]) && (rhs[i] != 0.0)) {
6689371c9d4SSatish Balay ct++;
6699371c9d4SSatish Balay nfo++;
670827bd09bSSatish Balay *--iptr = local2global[i];
671827bd09bSSatish Balay used[i] = edge;
672827bd09bSSatish Balay }
673827bd09bSSatish Balay }
6743ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
675827bd09bSSatish Balay lnsep[edge] = ct;
676827bd09bSSatish Balay }
677827bd09bSSatish Balay nsep[edge] = sum[1];
678827bd09bSSatish Balay dir[edge] = RIGHT;
679827bd09bSSatish Balay }
680827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */
681827bd09bSSatish Balay /* do we need full set of separators for this? */
682827bd09bSSatish Balay
683827bd09bSSatish Balay /* fold rhs hc into lower */
6842fa5cd67SKarl Rupp if (id >= mask) id -= mask;
685827bd09bSSatish Balay }
686827bd09bSSatish Balay }
687827bd09bSSatish Balay
688827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */
689db4deed7SKarl Rupp for (ct = i = 0; i < n; i++) {
690db4deed7SKarl Rupp if (!used[i]) {
6919371c9d4SSatish Balay ct++;
6929371c9d4SSatish Balay nfo++;
693827bd09bSSatish Balay *--iptr = local2global[i];
694827bd09bSSatish Balay used[i] = edge;
695827bd09bSSatish Balay }
696827bd09bSSatish Balay }
6973ba16761SJacob Faibussowitsch if (ct > 1) PetscCall(PCTFS_ivec_sort(iptr, ct));
698827bd09bSSatish Balay lnsep[edge] = ct;
699827bd09bSSatish Balay nsep[edge] = ct;
700827bd09bSSatish Balay dir[edge] = LEFT;
701827bd09bSSatish Balay
702827bd09bSSatish Balay xxt_handle->info->nsep = nsep;
703827bd09bSSatish Balay xxt_handle->info->lnsep = lnsep;
704827bd09bSSatish Balay xxt_handle->info->fo = fo;
705827bd09bSSatish Balay xxt_handle->info->nfo = nfo;
706827bd09bSSatish Balay
707a501084fSBarry Smith free(dir);
708a501084fSBarry Smith free(lhs);
709a501084fSBarry Smith free(rhs);
710a501084fSBarry Smith free(used);
7113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
712827bd09bSSatish Balay }
713827bd09bSSatish Balay
set_mvi(PetscInt * local2global,PetscInt n,PetscInt m,PetscErrorCode (* matvec)(mv_info *,PetscScalar *,PetscScalar *),void * grid_data)714d71ae5a4SJacob Faibussowitsch static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info *, PetscScalar *, PetscScalar *), void *grid_data)
715d71ae5a4SJacob Faibussowitsch {
716ff6a9541SJacob Faibussowitsch mv_info *mvi = (mv_info *)malloc(sizeof(mv_info));
717827bd09bSSatish Balay
718827bd09bSSatish Balay mvi->n = n;
719827bd09bSSatish Balay mvi->m = m;
720827bd09bSSatish Balay mvi->n_global = -1;
721827bd09bSSatish Balay mvi->m_global = -1;
72252f87cdaSBarry Smith mvi->local2global = (PetscInt *)malloc((m + 1) * sizeof(PetscInt));
723ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global, local2global, m);
724827bd09bSSatish Balay mvi->local2global[m] = INT_MAX;
7255c8f6a95SKarl Rupp mvi->matvec = matvec;
726827bd09bSSatish Balay mvi->grid_data = grid_data;
727827bd09bSSatish Balay
728827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */
729ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
730827bd09bSSatish Balay
731*4ad8454bSPierre Jolivet return mvi;
732827bd09bSSatish Balay }
733827bd09bSSatish Balay
do_matvec(mv_info * A,PetscScalar * v,PetscScalar * u)734d71ae5a4SJacob Faibussowitsch static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
735d71ae5a4SJacob Faibussowitsch {
7360924e98cSBarry Smith PetscFunctionBegin;
7373ba16761SJacob Faibussowitsch PetscCall(A->matvec((mv_info *)A->grid_data, v, u));
7383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
739827bd09bSSatish Balay }
740