xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /*************************************xxt.c************************************
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
18827bd09bSSatish Balay **************************************xxt.c***********************************/
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 
637b1ae94cSBarry Smith /**************************************xxt.c***********************************/
647b1ae94cSBarry Smith xxt_ADT XXT_new(void)
65827bd09bSSatish Balay {
66827bd09bSSatish Balay   xxt_ADT xxt_handle;
67827bd09bSSatish Balay 
68827bd09bSSatish Balay   /* rolling count on n_xxt ... pot. problem here */
69827bd09bSSatish Balay   n_xxt_handles++;
70a501084fSBarry Smith   xxt_handle       = (xxt_ADT)malloc(sizeof(struct xxt_CDT));
71827bd09bSSatish Balay   xxt_handle->id   = ++n_xxt;
72827bd09bSSatish Balay   xxt_handle->info = NULL; xxt_handle->mvi  = NULL;
73827bd09bSSatish Balay 
74827bd09bSSatish Balay   return(xxt_handle);
75827bd09bSSatish Balay }
76827bd09bSSatish Balay 
777b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7838dcbb09SBarry Smith PetscErrorCode XXT_factor(xxt_ADT xxt_handle,     /* prev. allocated xxt  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 {
85b1c944f5SJed Brown   PCTFS_comm_init();
86827bd09bSSatish Balay   check_handle(xxt_handle);
87827bd09bSSatish Balay 
88827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
892c71b3e2SJacob Faibussowitsch   PetscCheckFalse((1<<(xxt_handle->level=PCTFS_i_log2_num_nodes))!=PCTFS_num_nodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %D != %D",1<<PCTFS_i_log2_num_nodes,PCTFS_num_nodes);
90827bd09bSSatish Balay 
91827bd09bSSatish Balay   /* space for X info */
92a501084fSBarry Smith   xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info));
93827bd09bSSatish Balay 
94827bd09bSSatish Balay   /* set up matvec handles */
955c8f6a95SKarl Rupp   xxt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data);
96827bd09bSSatish Balay 
97827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
98827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
99827bd09bSSatish Balay   xxt_handle->ns=0;
100827bd09bSSatish Balay 
101827bd09bSSatish Balay   /* determine separators and generate firing order - NB xxt info set here */
102827bd09bSSatish Balay   det_separators(xxt_handle);
103827bd09bSSatish Balay 
104827bd09bSSatish Balay   return(do_xxt_factor(xxt_handle));
105827bd09bSSatish Balay }
106827bd09bSSatish Balay 
1077b1ae94cSBarry Smith /**************************************xxt.c***********************************/
10838dcbb09SBarry Smith PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b)
109827bd09bSSatish Balay {
110827bd09bSSatish Balay 
111b1c944f5SJed Brown   PCTFS_comm_init();
112827bd09bSSatish Balay   check_handle(xxt_handle);
113827bd09bSSatish Balay 
114827bd09bSSatish Balay   /* need to copy b into x? */
1152fa5cd67SKarl Rupp   if (b) PCTFS_rvec_copy(x,b,xxt_handle->mvi->n);
11638dcbb09SBarry Smith   return do_xxt_solve(xxt_handle,x);
117827bd09bSSatish Balay }
118827bd09bSSatish Balay 
1197b1ae94cSBarry Smith /**************************************xxt.c***********************************/
12052f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle)
121827bd09bSSatish Balay {
122827bd09bSSatish Balay 
123b1c944f5SJed Brown   PCTFS_comm_init();
124827bd09bSSatish Balay   check_handle(xxt_handle);
125827bd09bSSatish Balay   n_xxt_handles--;
126827bd09bSSatish Balay 
127a501084fSBarry Smith   free(xxt_handle->info->nsep);
128a501084fSBarry Smith   free(xxt_handle->info->lnsep);
129a501084fSBarry Smith   free(xxt_handle->info->fo);
130a501084fSBarry Smith   free(xxt_handle->info->stages);
131a501084fSBarry Smith   free(xxt_handle->info->solve_uu);
132a501084fSBarry Smith   free(xxt_handle->info->solve_w);
133a501084fSBarry Smith   free(xxt_handle->info->x);
134a501084fSBarry Smith   free(xxt_handle->info->col_vals);
135a501084fSBarry Smith   free(xxt_handle->info->col_sz);
136a501084fSBarry Smith   free(xxt_handle->info->col_indices);
137a501084fSBarry Smith   free(xxt_handle->info);
138a501084fSBarry Smith   free(xxt_handle->mvi->local2global);
139ca8e9878SJed Brown   PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle);
140a501084fSBarry Smith   free(xxt_handle->mvi);
141a501084fSBarry Smith   free(xxt_handle);
142827bd09bSSatish Balay 
143827bd09bSSatish Balay   /* if the check fails we nuke */
144a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
145827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
146827bd09bSSatish Balay   return(0);
147827bd09bSSatish Balay }
148827bd09bSSatish Balay 
1497b1ae94cSBarry Smith /**************************************xxt.c***********************************/
15038dcbb09SBarry Smith PetscErrorCode XXT_stats(xxt_ADT xxt_handle)
151827bd09bSSatish Balay {
15252f87cdaSBarry Smith   PetscInt       op[]  = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
15352f87cdaSBarry Smith   PetscInt       fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
15452f87cdaSBarry Smith   PetscInt       vals[9],  work[9];
155a501084fSBarry Smith   PetscScalar    fvals[3], fwork[3];
156827bd09bSSatish Balay 
157b1c944f5SJed Brown   PCTFS_comm_init();
158827bd09bSSatish Balay   check_handle(xxt_handle);
159827bd09bSSatish Balay 
160827bd09bSSatish Balay   /* if factorization not done there are no stats */
161db4deed7SKarl Rupp   if (!xxt_handle->info||!xxt_handle->mvi) {
1625f80ce2aSJacob Faibussowitsch     if (!PCTFS_my_id) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n"));
163827bd09bSSatish Balay     return 1;
164827bd09bSSatish Balay   }
165827bd09bSSatish Balay 
166827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
167827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
168827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
169b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
170827bd09bSSatish Balay 
1712fa5cd67SKarl Rupp   fvals[0]=fvals[1]=fvals[2] =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
172b1c944f5SJed Brown   PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
173827bd09bSSatish Balay 
174db4deed7SKarl Rupp   if (!PCTFS_my_id) {
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_nnz=%D\n",PCTFS_my_id,vals[0]));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_nnz=%D\n",PCTFS_my_id,vals[1]));
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes));
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_nnz=%D\n",PCTFS_my_id,vals[2]));
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(2d)  =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.5))));
1805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(3d)  =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.6667))));
1815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_n  =%D\n",PCTFS_my_id,vals[3]));
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_n  =%D\n",PCTFS_my_id,vals[4]));
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_n  =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes));
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_n  =%D\n",PCTFS_my_id,vals[5]));
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_buf=%D\n",PCTFS_my_id,vals[6]));
1865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_buf=%D\n",PCTFS_my_id,vals[7]));
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_slv=%g\n",PCTFS_my_id,fvals[0]));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_slv=%g\n",PCTFS_my_id,fvals[1]));
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes));
191827bd09bSSatish Balay   }
192827bd09bSSatish Balay 
193827bd09bSSatish Balay   return(0);
194827bd09bSSatish Balay }
195827bd09bSSatish Balay 
196827bd09bSSatish Balay /*************************************xxt.c************************************
197827bd09bSSatish Balay 
198827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
199827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
200827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
201827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
202827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
203827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
204ca8e9878SJed Brown    PCTFS_gs_init/gop).
205827bd09bSSatish Balay 
206827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
207827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
208827bd09bSSatish Balay **************************************xxt.c***********************************/
20938dcbb09SBarry Smith static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle)
210827bd09bSSatish Balay {
2117b1ae94cSBarry Smith   return xxt_generate(xxt_handle);
212827bd09bSSatish Balay }
213827bd09bSSatish Balay 
2147b1ae94cSBarry Smith /**************************************xxt.c***********************************/
21538dcbb09SBarry Smith static PetscErrorCode xxt_generate(xxt_ADT xxt_handle)
216827bd09bSSatish Balay {
21752f87cdaSBarry Smith   PetscInt       i,j,k,idex;
21852f87cdaSBarry Smith   PetscInt       dim, col;
219a501084fSBarry Smith   PetscScalar    *u, *uu, *v, *z, *w, alpha, alpha_w;
22052f87cdaSBarry Smith   PetscInt       *segs;
22152f87cdaSBarry Smith   PetscInt       op[] = {GL_ADD,0};
22252f87cdaSBarry Smith   PetscInt       off, len;
223a501084fSBarry Smith   PetscScalar    *x_ptr;
22452f87cdaSBarry Smith   PetscInt       *iptr, flag;
22552f87cdaSBarry Smith   PetscInt       start =0, end, work;
22652f87cdaSBarry Smith   PetscInt       op2[] = {GL_MIN,0};
227ca8e9878SJed Brown   PCTFS_gs_ADT   PCTFS_gs_handle;
22852f87cdaSBarry Smith   PetscInt       *nsep, *lnsep, *fo;
22952f87cdaSBarry Smith   PetscInt       a_n            =xxt_handle->mvi->n;
23052f87cdaSBarry Smith   PetscInt       a_m            =xxt_handle->mvi->m;
23152f87cdaSBarry Smith   PetscInt       *a_local2global=xxt_handle->mvi->local2global;
23252f87cdaSBarry Smith   PetscInt       level;
23352f87cdaSBarry Smith   PetscInt       xxt_nnz=0, xxt_max_nnz=0;
23452f87cdaSBarry Smith   PetscInt       n, m;
23552f87cdaSBarry Smith   PetscInt       *col_sz, *col_indices, *stages;
236a501084fSBarry Smith   PetscScalar    **col_vals, *x;
23752f87cdaSBarry Smith   PetscInt       n_global;
23852f87cdaSBarry Smith   PetscInt       xxt_zero_nnz  =0;
23952f87cdaSBarry Smith   PetscInt       xxt_zero_nnz_0=0;
24071a0148aSBarry Smith   PetscBLASInt   i1            = 1,dlen;
241a501084fSBarry Smith   PetscScalar    dm1           = -1.0;
242827bd09bSSatish Balay 
243827bd09bSSatish Balay   n               = xxt_handle->mvi->n;
244827bd09bSSatish Balay   nsep            = xxt_handle->info->nsep;
245827bd09bSSatish Balay   lnsep           = xxt_handle->info->lnsep;
246827bd09bSSatish Balay   fo              = xxt_handle->info->fo;
247827bd09bSSatish Balay   end             = lnsep[0];
248827bd09bSSatish Balay   level           = xxt_handle->level;
249ca8e9878SJed Brown   PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
250827bd09bSSatish Balay 
251827bd09bSSatish Balay   /* is there a null space? */
252827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
2532fa5cd67SKarl Rupp   for (i=0, j=0; i<=level; i++) j+=nsep[i];
254827bd09bSSatish Balay 
255827bd09bSSatish Balay   m = j-xxt_handle->ns;
25622d28d08SBarry Smith   if (m!=j) {
2575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns));
25822d28d08SBarry Smith   }
259827bd09bSSatish Balay 
260827bd09bSSatish Balay   /* get and initialize storage for x local         */
261827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
26252f87cdaSBarry Smith   col_sz      = (PetscInt*) malloc(m*sizeof(PetscInt));
26352f87cdaSBarry Smith   col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
264a501084fSBarry Smith   col_vals    = (PetscScalar**) malloc(m*sizeof(PetscScalar*));
265db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
266827bd09bSSatish Balay     col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
267827bd09bSSatish Balay     col_vals[i]   = NULL;
268827bd09bSSatish Balay   }
269827bd09bSSatish Balay   col_indices[j]=-1;
270827bd09bSSatish Balay 
271827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
272827bd09bSSatish Balay   /* this looks like nsep[]=segments */
27352f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
27452f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
275ca8e9878SJed Brown   PCTFS_ivec_zero(stages,level+1);
276ca8e9878SJed Brown   PCTFS_ivec_copy(segs,nsep,level+1);
2772fa5cd67SKarl Rupp   for (i=0; i<level; i++) segs[i+1] += segs[i];
278827bd09bSSatish Balay   stages[0] = segs[0];
279827bd09bSSatish Balay 
280827bd09bSSatish Balay   /* temporary vectors  */
281a501084fSBarry Smith   u  = (PetscScalar*) malloc(n*sizeof(PetscScalar));
282a501084fSBarry Smith   z  = (PetscScalar*) malloc(n*sizeof(PetscScalar));
283a501084fSBarry Smith   v  = (PetscScalar*) malloc(a_m*sizeof(PetscScalar));
284a501084fSBarry Smith   uu = (PetscScalar*) malloc(m*sizeof(PetscScalar));
285a501084fSBarry Smith   w  = (PetscScalar*) malloc(m*sizeof(PetscScalar));
286827bd09bSSatish Balay 
287827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
2882fa5cd67SKarl Rupp   for (i=1, j=0; i<=level; i++) j+=nsep[i];
289827bd09bSSatish Balay 
290827bd09bSSatish Balay   /* storage for sparse x values */
291827bd09bSSatish Balay   n_global    = xxt_handle->info->n_global;
29285ec1a3cSBarry Smith   xxt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes;
293a501084fSBarry Smith   x           = (PetscScalar*) malloc(xxt_max_nnz*sizeof(PetscScalar));
294827bd09bSSatish Balay   xxt_nnz     = 0;
295827bd09bSSatish Balay 
296827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
297827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
2982fa5cd67SKarl Rupp   for (dim=i=j=0; i<m; i++) {
299827bd09bSSatish Balay     /* time to move to the next level? */
300d4af0d30SBarry Smith     while (i==segs[dim]) {
3012c71b3e2SJacob Faibussowitsch       PetscCheckFalse(dim==level,PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level");
302827bd09bSSatish Balay       stages[dim++]=i;
303827bd09bSSatish Balay       end         +=lnsep[dim];
304827bd09bSSatish Balay     }
305827bd09bSSatish Balay     stages[dim]=i;
306827bd09bSSatish Balay 
307827bd09bSSatish Balay     /* which column are we firing? */
308827bd09bSSatish Balay     /* i.e. set v_l */
309827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
310827bd09bSSatish Balay     (start<end) ? (col=fo[start]) : (col=INT_MAX);
311b1c944f5SJed Brown     PCTFS_giop_hc(&col,&work,1,op2,dim);
312827bd09bSSatish Balay 
313827bd09bSSatish Balay     /* shouldn't need this */
314db4deed7SKarl Rupp     if (col==INT_MAX) {
3155f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(0,"hey ... col==INT_MAX??\n"));
316827bd09bSSatish Balay       continue;
317827bd09bSSatish Balay     }
318827bd09bSSatish Balay 
319827bd09bSSatish Balay     /* do I own it? I should */
320ca8e9878SJed Brown     PCTFS_rvec_zero(v,a_m);
321db4deed7SKarl Rupp     if (col==fo[start]) {
322827bd09bSSatish Balay       start++;
323ca8e9878SJed Brown       idex=PCTFS_ivec_linear_search(col, a_local2global, a_n);
324e7e72b3dSBarry Smith       if (idex!=-1) {
325e7e72b3dSBarry Smith         v[idex] = 1.0; j++;
326546078acSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!");
327db4deed7SKarl Rupp     } else {
328ca8e9878SJed Brown       idex=PCTFS_ivec_linear_search(col, a_local2global, a_m);
3292fa5cd67SKarl Rupp       if (idex!=-1) v[idex] = 1.0;
330827bd09bSSatish Balay     }
331827bd09bSSatish Balay 
332827bd09bSSatish Balay     /* perform u = A.v_l */
333ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
334827bd09bSSatish Balay     do_matvec(xxt_handle->mvi,v,u);
335827bd09bSSatish Balay 
336827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
337827bd09bSSatish Balay     /* technically only need to zero out first i entries */
338827bd09bSSatish Balay     /* later turn this into an XXT_solve call ? */
339ca8e9878SJed Brown     PCTFS_rvec_zero(uu,m);
340827bd09bSSatish Balay     x_ptr=x;
341827bd09bSSatish Balay     iptr = col_indices;
342db4deed7SKarl Rupp     for (k=0; k<i; k++) {
343827bd09bSSatish Balay       off   = *iptr++;
344827bd09bSSatish Balay       len   = *iptr++;
3455f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBLASIntCast(len,&dlen));
3468b83055fSJed Brown       PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1));
347827bd09bSSatish Balay       x_ptr+=len;
348827bd09bSSatish Balay     }
349827bd09bSSatish Balay 
350827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
3515f80ce2aSJacob Faibussowitsch     CHKERRQ(PCTFS_ssgl_radd  (uu, w, dim, stages));
352827bd09bSSatish Balay 
353827bd09bSSatish Balay     /* z = X.uu */
354ca8e9878SJed Brown     PCTFS_rvec_zero(z,n);
355827bd09bSSatish Balay     x_ptr=x;
356827bd09bSSatish Balay     iptr = col_indices;
357db4deed7SKarl Rupp     for (k=0; k<i; k++) {
358827bd09bSSatish Balay       off  = *iptr++;
359827bd09bSSatish Balay       len  = *iptr++;
3605f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBLASIntCast(len,&dlen));
3618b83055fSJed Brown       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1));
362827bd09bSSatish Balay       x_ptr+=len;
363827bd09bSSatish Balay     }
364827bd09bSSatish Balay 
365827bd09bSSatish Balay     /* compute v_l = v_l - z */
366ca8e9878SJed Brown     PCTFS_rvec_zero(v+a_n,a_m-a_n);
3675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(n,&dlen));
3688b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1));
369827bd09bSSatish Balay 
370827bd09bSSatish Balay     /* compute u_l = A.v_l */
3712fa5cd67SKarl Rupp     if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim);
372ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
373827bd09bSSatish Balay     do_matvec(xxt_handle->mvi,v,u);
374827bd09bSSatish Balay 
375827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
3765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(n,&dlen));
3778b83055fSJed Brown     PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,v,&i1));
378827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
379b1c944f5SJed Brown     PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
380827bd09bSSatish Balay 
3818f1a2a5eSBarry Smith     alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha);
382827bd09bSSatish Balay 
383827bd09bSSatish Balay     /* check for small alpha                             */
384827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
3852c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsScalar(alpha)<1.0e-14,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g",alpha);
386827bd09bSSatish Balay 
387827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
388ca8e9878SJed Brown     PCTFS_rvec_scale(v,1.0/alpha,n);
389827bd09bSSatish Balay 
390827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
391827bd09bSSatish Balay     flag = 1;
392827bd09bSSatish Balay     off=len=0;
393db4deed7SKarl Rupp     for (k=0; k<n; k++) {
394db4deed7SKarl Rupp       if (v[k]!=0.0) {
395827bd09bSSatish Balay         len=k;
396db4deed7SKarl Rupp         if (flag) { off=k; flag=0; }
397827bd09bSSatish Balay       }
398827bd09bSSatish Balay     }
399827bd09bSSatish Balay 
400827bd09bSSatish Balay     len -= (off-1);
401827bd09bSSatish Balay 
402db4deed7SKarl Rupp     if (len>0) {
403db4deed7SKarl Rupp       if ((xxt_nnz+len)>xxt_max_nnz) {
4045f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscInfo(0,"increasing space for X by 2x!\n"));
405827bd09bSSatish Balay         xxt_max_nnz *= 2;
406a501084fSBarry Smith         x_ptr        = (PetscScalar*) malloc(xxt_max_nnz*sizeof(PetscScalar));
407ca8e9878SJed Brown         PCTFS_rvec_copy(x_ptr,x,xxt_nnz);
408a501084fSBarry Smith         free(x);
409827bd09bSSatish Balay         x     = x_ptr;
410827bd09bSSatish Balay         x_ptr+=xxt_nnz;
411827bd09bSSatish Balay       }
412827bd09bSSatish Balay       xxt_nnz += len;
413ca8e9878SJed Brown       PCTFS_rvec_copy(x_ptr,v+off,len);
414827bd09bSSatish Balay 
415827bd09bSSatish Balay       /* keep track of number of zeros */
416db4deed7SKarl Rupp       if (dim) {
417db4deed7SKarl Rupp         for (k=0; k<len; k++) {
4182fa5cd67SKarl Rupp           if (x_ptr[k]==0.0) xxt_zero_nnz++;
419827bd09bSSatish Balay         }
420db4deed7SKarl Rupp       } else {
421db4deed7SKarl Rupp         for (k=0; k<len; k++) {
4222fa5cd67SKarl Rupp           if (x_ptr[k]==0.0) xxt_zero_nnz_0++;
423827bd09bSSatish Balay         }
424827bd09bSSatish Balay       }
425827bd09bSSatish Balay       col_indices[2*i] = off;
426827bd09bSSatish Balay       col_sz[i] = col_indices[2*i+1] = len;
427827bd09bSSatish Balay       col_vals[i] = x_ptr;
428827bd09bSSatish Balay     }
429db4deed7SKarl Rupp     else {
430827bd09bSSatish Balay       col_indices[2*i] = 0;
431827bd09bSSatish Balay       col_sz[i]        = col_indices[2*i+1] = 0;
432827bd09bSSatish Balay       col_vals[i]      = x_ptr;
433827bd09bSSatish Balay     }
434827bd09bSSatish Balay   }
435827bd09bSSatish Balay 
436827bd09bSSatish Balay   /* close off stages for execution phase */
437db4deed7SKarl Rupp   while (dim!=level) {
438827bd09bSSatish Balay     stages[dim++] = i;
4395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level));
440827bd09bSSatish Balay   }
441827bd09bSSatish Balay   stages[dim]=i;
442827bd09bSSatish Balay 
443827bd09bSSatish Balay   xxt_handle->info->n              = xxt_handle->mvi->n;
444827bd09bSSatish Balay   xxt_handle->info->m              = m;
445827bd09bSSatish Balay   xxt_handle->info->nnz            = xxt_nnz;
446827bd09bSSatish Balay   xxt_handle->info->max_nnz        = xxt_max_nnz;
447827bd09bSSatish Balay   xxt_handle->info->msg_buf_sz     = stages[level]-stages[0];
448a501084fSBarry Smith   xxt_handle->info->solve_uu       = (PetscScalar*) malloc(m*sizeof(PetscScalar));
449a501084fSBarry Smith   xxt_handle->info->solve_w        = (PetscScalar*) malloc(m*sizeof(PetscScalar));
450827bd09bSSatish Balay   xxt_handle->info->x              = x;
451827bd09bSSatish Balay   xxt_handle->info->col_vals       = col_vals;
452827bd09bSSatish Balay   xxt_handle->info->col_sz         = col_sz;
453827bd09bSSatish Balay   xxt_handle->info->col_indices    = col_indices;
454827bd09bSSatish Balay   xxt_handle->info->stages         = stages;
455827bd09bSSatish Balay   xxt_handle->info->nsolves        = 0;
456827bd09bSSatish Balay   xxt_handle->info->tot_solve_time = 0.0;
457827bd09bSSatish Balay 
458a501084fSBarry Smith   free(segs);
459a501084fSBarry Smith   free(u);
460a501084fSBarry Smith   free(v);
461a501084fSBarry Smith   free(uu);
462a501084fSBarry Smith   free(z);
463a501084fSBarry Smith   free(w);
464827bd09bSSatish Balay 
465827bd09bSSatish Balay   return(0);
466827bd09bSSatish Balay }
467827bd09bSSatish Balay 
4687b1ae94cSBarry Smith /**************************************xxt.c***********************************/
4690924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle,  PetscScalar *uc)
470827bd09bSSatish Balay {
47152f87cdaSBarry Smith   PetscInt       off, len, *iptr;
47252f87cdaSBarry Smith   PetscInt       level        = xxt_handle->level;
47352f87cdaSBarry Smith   PetscInt       n            = xxt_handle->info->n;
47452f87cdaSBarry Smith   PetscInt       m            = xxt_handle->info->m;
47552f87cdaSBarry Smith   PetscInt       *stages      = xxt_handle->info->stages;
47652f87cdaSBarry Smith   PetscInt       *col_indices = xxt_handle->info->col_indices;
477a501084fSBarry Smith   PetscScalar    *x_ptr, *uu_ptr;
478a501084fSBarry Smith   PetscScalar    *solve_uu = xxt_handle->info->solve_uu;
479a501084fSBarry Smith   PetscScalar    *solve_w  = xxt_handle->info->solve_w;
480a501084fSBarry Smith   PetscScalar    *x        = xxt_handle->info->x;
48171a0148aSBarry Smith   PetscBLASInt   i1        = 1,dlen;
482827bd09bSSatish Balay 
4830924e98cSBarry Smith   PetscFunctionBegin;
484827bd09bSSatish Balay   uu_ptr=solve_uu;
485ca8e9878SJed Brown   PCTFS_rvec_zero(uu_ptr,m);
486827bd09bSSatish Balay 
487827bd09bSSatish Balay   /* x  = X.Y^T.b */
488827bd09bSSatish Balay   /* uu = Y^T.b */
489db4deed7SKarl Rupp   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) {
490c5df96a5SBarry Smith     off       =*iptr++;
491c5df96a5SBarry Smith     len       =*iptr++;
4925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(len,&dlen));
4938b83055fSJed Brown     PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1));
494827bd09bSSatish Balay   }
495827bd09bSSatish Balay 
496827bd09bSSatish Balay   /* comunication of beta */
497827bd09bSSatish Balay   uu_ptr=solve_uu;
4985f80ce2aSJacob Faibussowitsch   if (level) CHKERRQ(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages));
499827bd09bSSatish Balay 
500ca8e9878SJed Brown   PCTFS_rvec_zero(uc,n);
501827bd09bSSatish Balay 
502827bd09bSSatish Balay   /* x = X.uu */
503db4deed7SKarl Rupp   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) {
504c5df96a5SBarry Smith     off  =*iptr++;
505c5df96a5SBarry Smith     len  =*iptr++;
5065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBLASIntCast(len,&dlen));
5078b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1));
508827bd09bSSatish Balay   }
5090924e98cSBarry Smith   PetscFunctionReturn(0);
510827bd09bSSatish Balay }
511827bd09bSSatish Balay 
5127b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5130924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle)
514827bd09bSSatish Balay {
51552f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
516827bd09bSSatish Balay 
5170924e98cSBarry Smith   PetscFunctionBegin;
518*28b400f6SJacob Faibussowitsch   PetscCheck(xxt_handle,PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D",xxt_handle);
519827bd09bSSatish Balay 
520827bd09bSSatish Balay   vals[0]=vals[1]=xxt_handle->id;
521b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
5222c71b3e2SJacob Faibussowitsch   PetscCheckFalse((vals[0]!=vals[1])||(xxt_handle->id<=0),PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %D/%D %D",vals[0],vals[1], xxt_handle->id);
5230924e98cSBarry Smith   PetscFunctionReturn(0);
524827bd09bSSatish Balay }
525827bd09bSSatish Balay 
5267b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5270924e98cSBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle)
528827bd09bSSatish Balay {
52952f87cdaSBarry Smith   PetscInt     i, ct, id;
53052f87cdaSBarry Smith   PetscInt     mask, edge, *iptr;
53152f87cdaSBarry Smith   PetscInt     *dir, *used;
53252f87cdaSBarry Smith   PetscInt     sum[4], w[4];
533a501084fSBarry Smith   PetscScalar  rsum[4], rw[4];
53452f87cdaSBarry Smith   PetscInt     op[] = {GL_ADD,0};
535a501084fSBarry Smith   PetscScalar  *lhs, *rhs;
53652f87cdaSBarry Smith   PetscInt     *nsep, *lnsep, *fo, nfo=0;
537ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
53852f87cdaSBarry Smith   PetscInt     *local2global   = xxt_handle->mvi->local2global;
53952f87cdaSBarry Smith   PetscInt     n               = xxt_handle->mvi->n;
54052f87cdaSBarry Smith   PetscInt     m               = xxt_handle->mvi->m;
54152f87cdaSBarry Smith   PetscInt     level           = xxt_handle->level;
542ab824b78SBarry Smith   PetscInt     shared          = 0;
543827bd09bSSatish Balay 
5440924e98cSBarry Smith   PetscFunctionBegin;
54552f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
54652f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
54752f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
54852f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
54952f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
550827bd09bSSatish Balay 
551ca8e9878SJed Brown   PCTFS_ivec_zero(dir,level+1);
552ca8e9878SJed Brown   PCTFS_ivec_zero(nsep,level+1);
553ca8e9878SJed Brown   PCTFS_ivec_zero(lnsep,level+1);
554ca8e9878SJed Brown   PCTFS_ivec_set (fo,-1,n+1);
555ca8e9878SJed Brown   PCTFS_ivec_zero(used,n);
556827bd09bSSatish Balay 
5578cda6cd7SBarry Smith   lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
5588cda6cd7SBarry Smith   rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
559827bd09bSSatish Balay 
560827bd09bSSatish Balay   /* determine the # of unique dof */
561ca8e9878SJed Brown   PCTFS_rvec_zero(lhs,m);
562ca8e9878SJed Brown   PCTFS_rvec_set(lhs,1.0,n);
563ca8e9878SJed Brown   PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level);
564ca8e9878SJed Brown   PCTFS_rvec_zero(rsum,2);
5653d3eaba7SBarry Smith   for (i=0; i<n; i++) {
5662fa5cd67SKarl Rupp     if (lhs[i]!=0.0) {
5672fa5cd67SKarl Rupp       rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];
5682fa5cd67SKarl Rupp     }
569827bd09bSSatish Balay   }
570b1c944f5SJed Brown   PCTFS_grop_hc(rsum,rw,2,op,level);
571827bd09bSSatish Balay   rsum[0]+=0.1;
572827bd09bSSatish Balay   rsum[1]+=0.1;
573827bd09bSSatish Balay 
57477b4d14cSPeter Brune   if (PetscAbsScalar(rsum[0]-rsum[1])>EPS) shared=1;
575827bd09bSSatish Balay 
57652f87cdaSBarry Smith   xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0];
57752f87cdaSBarry Smith   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0];
578827bd09bSSatish Balay 
579827bd09bSSatish Balay   /* determine separator sets top down */
5802fa5cd67SKarl Rupp   if (shared) {
581db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
582db4deed7SKarl Rupp 
583827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
584ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
585ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
586827bd09bSSatish Balay 
587827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
588ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
589ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
590827bd09bSSatish Balay 
591db4deed7SKarl Rupp       for (i=0;i<n;i++) {
592db4deed7SKarl Rupp         if (id< mask) {
5932fa5cd67SKarl Rupp           if (lhs[i]!=0.0) lhs[i]=1.0;
594827bd09bSSatish Balay         }
595db4deed7SKarl Rupp         if (id>=mask) {
5962fa5cd67SKarl Rupp           if (rhs[i]!=0.0) rhs[i]=1.0;
597827bd09bSSatish Balay         }
598827bd09bSSatish Balay       }
599827bd09bSSatish Balay 
6002fa5cd67SKarl Rupp       if (id< mask) PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1);
6012fa5cd67SKarl Rupp       else          PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1);
602827bd09bSSatish Balay 
603827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
604ca8e9878SJed Brown       PCTFS_rvec_zero(rsum,4);
605db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
606db4deed7SKarl Rupp         if (!used[i]) {
607827bd09bSSatish Balay           /* number of unmarked dofs on node */
608827bd09bSSatish Balay           ct++;
609827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
610db4deed7SKarl Rupp           if (id< mask) {
611db4deed7SKarl Rupp             if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; }
612827bd09bSSatish Balay           }
613827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
614db4deed7SKarl Rupp           if (id>=mask) {
615db4deed7SKarl Rupp             if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; }
616827bd09bSSatish Balay           }
617827bd09bSSatish Balay         }
618827bd09bSSatish Balay       }
619827bd09bSSatish Balay 
620827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
621827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
622827bd09bSSatish Balay       (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
623b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
624b1c944f5SJed Brown       PCTFS_grop_hc(rsum,rw,4,op,edge);
625827bd09bSSatish Balay       rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
626827bd09bSSatish Balay 
627db4deed7SKarl Rupp       if (id<mask) {
628827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
629db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
630db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
631827bd09bSSatish Balay             ct++; nfo++;
632827bd09bSSatish Balay 
6332c71b3e2SJacob Faibussowitsch             PetscCheckFalse(nfo>n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n");
634827bd09bSSatish Balay 
635827bd09bSSatish Balay             *--iptr = local2global[i];
636827bd09bSSatish Balay             used[i] = edge;
637827bd09bSSatish Balay           }
638827bd09bSSatish Balay         }
6392fa5cd67SKarl Rupp         if (ct>1) PCTFS_ivec_sort(iptr,ct);
640827bd09bSSatish Balay 
641827bd09bSSatish Balay         lnsep[edge]=ct;
64252f87cdaSBarry Smith         nsep[edge]=(PetscInt) rsum[0];
643827bd09bSSatish Balay         dir [edge]=LEFT;
644827bd09bSSatish Balay       }
645827bd09bSSatish Balay 
646db4deed7SKarl Rupp       if (id>=mask) {
647827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
648db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
649db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
650827bd09bSSatish Balay             ct++; nfo++;
651827bd09bSSatish Balay 
6522c71b3e2SJacob Faibussowitsch             PetscCheckFalse(nfo>n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n");
653827bd09bSSatish Balay 
654827bd09bSSatish Balay             *--iptr = local2global[i];
655827bd09bSSatish Balay             used[i] = edge;
656827bd09bSSatish Balay           }
657827bd09bSSatish Balay         }
6582fa5cd67SKarl Rupp         if (ct>1) PCTFS_ivec_sort(iptr,ct);
659827bd09bSSatish Balay 
660827bd09bSSatish Balay         lnsep[edge] = ct;
66152f87cdaSBarry Smith         nsep[edge]  = (PetscInt) rsum[1];
662827bd09bSSatish Balay         dir [edge]  = RIGHT;
663827bd09bSSatish Balay       }
664827bd09bSSatish Balay 
665827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
666827bd09bSSatish Balay       /* do we need full set of separators for this?                */
667827bd09bSSatish Balay 
668827bd09bSSatish Balay       /* fold rhs hc into lower */
6692fa5cd67SKarl Rupp       if (id>=mask) id-=mask;
670827bd09bSSatish Balay     }
671db4deed7SKarl Rupp   } else {
672db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
673827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
674ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
675ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
676827bd09bSSatish Balay 
677827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
678ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
679ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
680827bd09bSSatish Balay 
681827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
682db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
683db4deed7SKarl Rupp         if (!used[i]) {
684827bd09bSSatish Balay           /* number of unmarked dofs on node */
685827bd09bSSatish Balay           ct++;
686827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
6872fa5cd67SKarl Rupp           if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++;
688827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
6892fa5cd67SKarl Rupp           if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++;
690827bd09bSSatish Balay         }
691827bd09bSSatish Balay       }
692827bd09bSSatish Balay 
693827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
694827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
695b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
696827bd09bSSatish Balay 
697827bd09bSSatish Balay       /* lhs hc wins */
698db4deed7SKarl Rupp       if (sum[2]>=sum[3]) {
699db4deed7SKarl Rupp         if (id<mask) {
700827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
701db4deed7SKarl Rupp           for (ct=i=0;i<n;i++) {
702db4deed7SKarl Rupp             if ((!used[i])&&(lhs[i]!=0.0)) {
703827bd09bSSatish Balay               ct++; nfo++;
704827bd09bSSatish Balay               *--iptr = local2global[i];
705827bd09bSSatish Balay               used[i]=edge;
706827bd09bSSatish Balay             }
707827bd09bSSatish Balay           }
7082fa5cd67SKarl Rupp           if (ct>1) PCTFS_ivec_sort(iptr,ct);
709827bd09bSSatish Balay           lnsep[edge]=ct;
710827bd09bSSatish Balay         }
711827bd09bSSatish Balay         nsep[edge]=sum[0];
712827bd09bSSatish Balay         dir [edge]=LEFT;
713db4deed7SKarl Rupp       } else { /* rhs hc wins */
714db4deed7SKarl Rupp         if (id>=mask) {
715827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
716db4deed7SKarl Rupp           for (ct=i=0;i<n;i++) {
717db4deed7SKarl Rupp             if ((!used[i])&&(rhs[i]!=0.0)) {
718827bd09bSSatish Balay               ct++; nfo++;
719827bd09bSSatish Balay               *--iptr = local2global[i];
720827bd09bSSatish Balay               used[i]=edge;
721827bd09bSSatish Balay             }
722827bd09bSSatish Balay           }
7232fa5cd67SKarl Rupp           if (ct>1) PCTFS_ivec_sort(iptr,ct);
724827bd09bSSatish Balay           lnsep[edge]=ct;
725827bd09bSSatish Balay         }
726827bd09bSSatish Balay         nsep[edge]=sum[1];
727827bd09bSSatish Balay         dir [edge]=RIGHT;
728827bd09bSSatish Balay       }
729827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
730827bd09bSSatish Balay       /* do we need full set of separators for this?                */
731827bd09bSSatish Balay 
732827bd09bSSatish Balay       /* fold rhs hc into lower */
7332fa5cd67SKarl Rupp       if (id>=mask) id-=mask;
734827bd09bSSatish Balay     }
735827bd09bSSatish Balay   }
736827bd09bSSatish Balay 
737827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
738db4deed7SKarl Rupp   for (ct=i=0; i<n; i++) {
739db4deed7SKarl Rupp     if (!used[i]) {
740827bd09bSSatish Balay       ct++; nfo++;
741827bd09bSSatish Balay       *--iptr = local2global[i];
742827bd09bSSatish Balay       used[i] = edge;
743827bd09bSSatish Balay     }
744827bd09bSSatish Balay   }
7452fa5cd67SKarl Rupp   if (ct>1) PCTFS_ivec_sort(iptr,ct);
746827bd09bSSatish Balay   lnsep[edge]=ct;
747827bd09bSSatish Balay   nsep [edge]=ct;
748827bd09bSSatish Balay   dir  [edge]=LEFT;
749827bd09bSSatish Balay 
750827bd09bSSatish Balay   xxt_handle->info->nsep  = nsep;
751827bd09bSSatish Balay   xxt_handle->info->lnsep = lnsep;
752827bd09bSSatish Balay   xxt_handle->info->fo    = fo;
753827bd09bSSatish Balay   xxt_handle->info->nfo   = nfo;
754827bd09bSSatish Balay 
755a501084fSBarry Smith   free(dir);
756a501084fSBarry Smith   free(lhs);
757a501084fSBarry Smith   free(rhs);
758a501084fSBarry Smith   free(used);
7590924e98cSBarry Smith   PetscFunctionReturn(0);
760827bd09bSSatish Balay }
761827bd09bSSatish Balay 
7627b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7635c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global,PetscInt n,PetscInt m,PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*),void *grid_data)
764827bd09bSSatish Balay {
765827bd09bSSatish Balay   mv_info *mvi;
766827bd09bSSatish Balay 
767a501084fSBarry Smith   mvi               = (mv_info*)malloc(sizeof(mv_info));
768827bd09bSSatish Balay   mvi->n            = n;
769827bd09bSSatish Balay   mvi->m            = m;
770827bd09bSSatish Balay   mvi->n_global     = -1;
771827bd09bSSatish Balay   mvi->m_global     = -1;
77252f87cdaSBarry Smith   mvi->local2global = (PetscInt*)malloc((m+1)*sizeof(PetscInt));
773ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global,local2global,m);
774827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
7755c8f6a95SKarl Rupp   mvi->matvec          = matvec;
776827bd09bSSatish Balay   mvi->grid_data       = grid_data;
777827bd09bSSatish Balay 
778827bd09bSSatish Balay   /* set xxt communication handle to perform restricted matvec */
779ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
780827bd09bSSatish Balay 
781827bd09bSSatish Balay   return(mvi);
782827bd09bSSatish Balay }
783827bd09bSSatish Balay 
7847b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7850924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
786827bd09bSSatish Balay {
7870924e98cSBarry Smith   PetscFunctionBegin;
788827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
7890924e98cSBarry Smith   PetscFunctionReturn(0);
790827bd09bSSatish Balay }
791