xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision d1528f56e7efe90e41719e8af1e4ce61b838432e)
1 #define PETSCKSP_DLL
2 
3 /*************************************xxt.c************************************
4 Module Name: xxt
5 Module Info:
6 
7 author:  Henry M. Tufo III
8 e-mail:  hmt@asci.uchicago.edu
9 contact:
10 +--------------------------------+--------------------------------+
11 |MCS Division - Building 221     |Department of Computer Science  |
12 |Argonne National Laboratory     |Ryerson 152                     |
13 |9700 S. Cass Avenue             |The University of Chicago       |
14 |Argonne, IL  60439              |Chicago, IL  60637              |
15 |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
16 +--------------------------------+--------------------------------+
17 
18 Last Modification: 3.20.01
19 **************************************xxt.c***********************************/
20 #include "src/ksp/pc/impls/tfs/tfs.h"
21 
22 #define LEFT  -1
23 #define RIGHT  1
24 #define BOTH   0
25 
26 typedef struct xxt_solver_info {
27   PetscInt n, m, n_global, m_global;
28   PetscInt nnz, max_nnz, msg_buf_sz;
29   PetscInt *nsep, *lnsep, *fo, nfo, *stages;
30   PetscInt *col_sz, *col_indices;
31   PetscScalar **col_vals, *x, *solve_uu, *solve_w;
32   PetscInt nsolves;
33   PetscScalar tot_solve_time;
34 } xxt_info;
35 
36 typedef struct matvec_info {
37   PetscInt n, m, n_global, m_global;
38   PetscInt *local2global;
39   gs_ADT gs_handle;
40   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
41   void *grid_data;
42 } mv_info;
43 
44 struct xxt_CDT{
45   PetscInt id;
46   PetscInt ns;
47   PetscInt level;
48   xxt_info *info;
49   mv_info  *mvi;
50 };
51 
52 static PetscInt n_xxt=0;
53 static PetscInt n_xxt_handles=0;
54 
55 /* prototypes */
56 static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *rhs);
57 static PetscErrorCode check_handle(xxt_ADT xxt_handle);
58 static PetscErrorCode det_separators(xxt_ADT xxt_handle);
59 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
60 static PetscInt xxt_generate(xxt_ADT xxt_handle);
61 static PetscInt do_xxt_factor(xxt_ADT xxt_handle);
62 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data);
63 
64 /**************************************xxt.c***********************************/
65 xxt_ADT XXT_new(void)
66 {
67   xxt_ADT xxt_handle;
68 
69   /* rolling count on n_xxt ... pot. problem here */
70   n_xxt_handles++;
71   xxt_handle       = (xxt_ADT)malloc(sizeof(struct xxt_CDT));
72   xxt_handle->id   = ++n_xxt;
73   xxt_handle->info = NULL; xxt_handle->mvi  = NULL;
74 
75   return(xxt_handle);
76 }
77 
78 /**************************************xxt.c***********************************/
79 PetscInt XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt  handle */
80 	   PetscInt *local2global,  /* global column mapping       */
81 	   PetscInt n,              /* local num rows              */
82 	   PetscInt m,              /* local num cols              */
83 	   void *matvec,       /* b_loc=A_local.x_loc         */
84 	   void *grid_data     /* grid data for matvec        */
85 	   )
86 {
87   comm_init();
88   check_handle(xxt_handle);
89 
90   /* only 2^k for now and all nodes participating */
91   if ((1<<(xxt_handle->level=i_log2_num_nodes))!=num_nodes)
92     {SETERRQ2(PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %d != %d\n",1<<i_log2_num_nodes,num_nodes);}
93 
94   /* space for X info */
95   xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info));
96 
97   /* set up matvec handles */
98   xxt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);
99 
100   /* matrix is assumed to be of full rank */
101   /* LATER we can reset to indicate rank def. */
102   xxt_handle->ns=0;
103 
104   /* determine separators and generate firing order - NB xxt info set here */
105   det_separators(xxt_handle);
106 
107   return(do_xxt_factor(xxt_handle));
108 }
109 
110 /**************************************xxt.c***********************************/
111 PetscInt XXT_solve(xxt_ADT xxt_handle, double *x, double *b)
112 {
113 
114   comm_init();
115   check_handle(xxt_handle);
116 
117   /* need to copy b into x? */
118   if (b)
119     {rvec_copy(x,b,xxt_handle->mvi->n);}
120   do_xxt_solve(xxt_handle,x);
121 
122   return(0);
123 }
124 
125 /**************************************xxt.c***********************************/
126 PetscInt XXT_free(xxt_ADT xxt_handle)
127 {
128 
129   comm_init();
130   check_handle(xxt_handle);
131   n_xxt_handles--;
132 
133   free(xxt_handle->info->nsep);
134   free(xxt_handle->info->lnsep);
135   free(xxt_handle->info->fo);
136   free(xxt_handle->info->stages);
137   free(xxt_handle->info->solve_uu);
138   free(xxt_handle->info->solve_w);
139   free(xxt_handle->info->x);
140   free(xxt_handle->info->col_vals);
141   free(xxt_handle->info->col_sz);
142   free(xxt_handle->info->col_indices);
143   free(xxt_handle->info);
144   free(xxt_handle->mvi->local2global);
145    gs_free(xxt_handle->mvi->gs_handle);
146   free(xxt_handle->mvi);
147   free(xxt_handle);
148 
149   /* if the check fails we nuke */
150   /* if NULL pointer passed to free we nuke */
151   /* if the calls to free fail that's not my problem */
152   return(0);
153 }
154 
155 /**************************************xxt.c***********************************/
156 PetscInt XXT_stats(xxt_ADT xxt_handle)
157 {
158   PetscInt  op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
159   PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
160   PetscInt   vals[9],  work[9];
161   PetscScalar fvals[3], fwork[3];
162 
163   comm_init();
164   check_handle(xxt_handle);
165 
166   /* if factorization not done there are no stats */
167   if (!xxt_handle->info||!xxt_handle->mvi)
168     {
169       if (!my_id)
170 	{printf("XXT_stats() :: no stats available!\n");}
171       return 1;
172     }
173 
174   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
175   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
176   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
177   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
178 
179   fvals[0]=fvals[1]=fvals[2]
180     =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
181   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
182 
183   if (!my_id)
184     {
185       printf("%d :: min   xxt_nnz=%d\n",my_id,vals[0]);
186       printf("%d :: max   xxt_nnz=%d\n",my_id,vals[1]);
187       printf("%d :: avg   xxt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
188       printf("%d :: tot   xxt_nnz=%d\n",my_id,vals[2]);
189       printf("%d :: xxt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
190       printf("%d :: xxt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
191       printf("%d :: min   xxt_n  =%d\n",my_id,vals[3]);
192       printf("%d :: max   xxt_n  =%d\n",my_id,vals[4]);
193       printf("%d :: avg   xxt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
194       printf("%d :: tot   xxt_n  =%d\n",my_id,vals[5]);
195       printf("%d :: min   xxt_buf=%d\n",my_id,vals[6]);
196       printf("%d :: max   xxt_buf=%d\n",my_id,vals[7]);
197       printf("%d :: avg   xxt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
198       printf("%d :: min   xxt_slv=%g\n",my_id,fvals[0]);
199       printf("%d :: max   xxt_slv=%g\n",my_id,fvals[1]);
200       printf("%d :: avg   xxt_slv=%g\n",my_id,fvals[2]/num_nodes);
201     }
202 
203   return(0);
204 }
205 
206 /*************************************xxt.c************************************
207 
208 Description: get A_local, local portion of global coarse matrix which
209 is a row dist. nxm matrix w/ n<m.
210    o my_ml holds address of ML struct associated w/A_local and coarse grid
211    o local2global holds global number of column i (i=0,...,m-1)
212    o local2global holds global number of row    i (i=0,...,n-1)
213    o mylocmatvec performs A_local . vec_local (note that gs is performed using
214    gs_init/gop).
215 
216 mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
217 mylocmatvec (void :: void *data, double *in, double *out)
218 **************************************xxt.c***********************************/
219 static PetscInt do_xxt_factor(xxt_ADT xxt_handle)
220 {
221   return xxt_generate(xxt_handle);
222 }
223 
224 /**************************************xxt.c***********************************/
225 static PetscInt xxt_generate(xxt_ADT xxt_handle)
226 {
227   PetscInt i,j,k,idex;
228   PetscInt dim, col;
229   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
230   PetscInt *segs;
231   PetscInt op[] = {GL_ADD,0};
232   PetscInt off, len;
233   PetscScalar *x_ptr;
234   PetscInt *iptr, flag;
235   PetscInt start=0, end, work;
236   PetscInt op2[] = {GL_MIN,0};
237   gs_ADT gs_handle;
238   PetscInt *nsep, *lnsep, *fo;
239   PetscInt a_n=xxt_handle->mvi->n;
240   PetscInt a_m=xxt_handle->mvi->m;
241   PetscInt *a_local2global=xxt_handle->mvi->local2global;
242   PetscInt level;
243   PetscInt xxt_nnz=0, xxt_max_nnz=0;
244   PetscInt n, m;
245   PetscInt *col_sz, *col_indices, *stages;
246   PetscScalar **col_vals, *x;
247   PetscInt n_global;
248   PetscInt xxt_zero_nnz=0;
249   PetscInt xxt_zero_nnz_0=0;
250   PetscBLASInt i1 = 1;
251   PetscScalar dm1 = -1.0;
252   PetscErrorCode ierr;
253 
254   n=xxt_handle->mvi->n;
255   nsep=xxt_handle->info->nsep;
256   lnsep=xxt_handle->info->lnsep;
257   fo=xxt_handle->info->fo;
258   end=lnsep[0];
259   level=xxt_handle->level;
260   gs_handle=xxt_handle->mvi->gs_handle;
261 
262   /* is there a null space? */
263   /* LATER add in ability to detect null space by checking alpha */
264   for (i=0, j=0; i<=level; i++)
265     {j+=nsep[i];}
266 
267   m = j-xxt_handle->ns;
268   if (m!=j)
269     {printf("xxt_generate() :: null space exists %d %d %d\n",m,j,xxt_handle->ns);}
270 
271   /* get and initialize storage for x local         */
272   /* note that x local is nxm and stored by columns */
273   col_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
274   col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
275   col_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
276   for (i=j=0; i<m; i++, j+=2)
277     {
278       col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
279       col_vals[i] = NULL;
280     }
281   col_indices[j]=-1;
282 
283   /* size of separators for each sub-hc working from bottom of tree to top */
284   /* this looks like nsep[]=segments */
285   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
286   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
287   ivec_zero(stages,level+1);
288   ivec_copy(segs,nsep,level+1);
289   for (i=0; i<level; i++)
290     {segs[i+1] += segs[i];}
291   stages[0] = segs[0];
292 
293   /* temporary vectors  */
294   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
295   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
296   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
297   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
298   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
299 
300   /* extra nnz due to replication of vertices across separators */
301   for (i=1, j=0; i<=level; i++)
302     {j+=nsep[i];}
303 
304   /* storage for sparse x values */
305   n_global = xxt_handle->info->n_global;
306   xxt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
307   x = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
308   xxt_nnz = 0;
309 
310   /* LATER - can embed next sep to fire in gs */
311   /* time to make the donuts - generate X factor */
312   for (dim=i=j=0;i<m;i++)
313     {
314       /* time to move to the next level? */
315       while (i==segs[dim])
316 	{
317 	  if (dim==level)
318 	    {SETERRQ(PETSC_ERR_PLIB,"dim about to exceed level\n"); break;}
319 
320 	  stages[dim++]=i;
321 	  end+=lnsep[dim];
322 	}
323       stages[dim]=i;
324 
325       /* which column are we firing? */
326       /* i.e. set v_l */
327       /* use new seps and do global min across hc to determine which one to fire */
328       (start<end) ? (col=fo[start]) : (col=INT_MAX);
329       giop_hc(&col,&work,1,op2,dim);
330 
331       /* shouldn't need this */
332       if (col==INT_MAX)
333 	{
334 	  ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");
335 	  continue;
336 	}
337 
338       /* do I own it? I should */
339       rvec_zero(v ,a_m);
340       if (col==fo[start])
341 	{
342 	  start++;
343 	  idex=ivec_linear_search(col, a_local2global, a_n);
344 	  if (idex!=-1)
345 	    {v[idex] = 1.0; j++;}
346 	  else
347 	    {SETERRQ(PETSC_ERR_PLIB,"NOT FOUND!\n");}
348 	}
349       else
350 	{
351 	  idex=ivec_linear_search(col, a_local2global, a_m);
352 	  if (idex!=-1)
353 	    {v[idex] = 1.0;}
354 	}
355 
356       /* perform u = A.v_l */
357       rvec_zero(u,n);
358       do_matvec(xxt_handle->mvi,v,u);
359 
360       /* uu =  X^T.u_l (local portion) */
361       /* technically only need to zero out first i entries */
362       /* later turn this into an XXT_solve call ? */
363       rvec_zero(uu,m);
364       x_ptr=x;
365       iptr = col_indices;
366       for (k=0; k<i; k++)
367 	{
368 	  off = *iptr++;
369 	  len = *iptr++;
370 
371 	  uu[k] = BLASdot_(&len,u+off,&i1,x_ptr,&i1);
372 	  x_ptr+=len;
373 	}
374 
375 
376       /* uu = X^T.u_l (comm portion) */
377       ssgl_radd  (uu, w, dim, stages);
378 
379       /* z = X.uu */
380       rvec_zero(z,n);
381       x_ptr=x;
382       iptr = col_indices;
383       for (k=0; k<i; k++)
384 	{
385 	  off = *iptr++;
386 	  len = *iptr++;
387 
388 	  BLASaxpy_(&len,&uu[k],x_ptr,&i1,z+off,&i1);
389 	  x_ptr+=len;
390 	}
391 
392       /* compute v_l = v_l - z */
393       rvec_zero(v+a_n,a_m-a_n);
394       BLASaxpy_(&n,&dm1,z,&i1,v,&i1);
395 
396       /* compute u_l = A.v_l */
397       if (a_n!=a_m)
398 	{gs_gop_hc(gs_handle,v,"+\0",dim);}
399       rvec_zero(u,n);
400       do_matvec(xxt_handle->mvi,v,u);
401 
402       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
403       alpha = BLASdot_(&n,u,&i1,v,&i1);
404       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
405       grop_hc(&alpha, &alpha_w, 1, op, dim);
406 
407       alpha = (PetscScalar) sqrt((double)alpha);
408 
409       /* check for small alpha                             */
410       /* LATER use this to detect and determine null space */
411       if (fabs(alpha)<1.0e-14)
412 	{SETERRQ1(PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);}
413 
414       /* compute v_l = v_l/sqrt(alpha) */
415       rvec_scale(v,1.0/alpha,n);
416 
417       /* add newly generated column, v_l, to X */
418       flag = 1;
419       off=len=0;
420       for (k=0; k<n; k++)
421 	{
422 	  if (v[k]!=0.0)
423 	    {
424 	      len=k;
425 	      if (flag)
426 		{off=k; flag=0;}
427 	    }
428 	}
429 
430       len -= (off-1);
431 
432       if (len>0)
433 	{
434 	  if ((xxt_nnz+len)>xxt_max_nnz)
435 	    {
436 	      ierr = PetscInfo(0,"increasing space for X by 2x!\n");
437 	      xxt_max_nnz *= 2;
438 	      x_ptr = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
439 	      rvec_copy(x_ptr,x,xxt_nnz);
440 	      free(x);
441 	      x = x_ptr;
442 	      x_ptr+=xxt_nnz;
443 	    }
444 	  xxt_nnz += len;
445 	  rvec_copy(x_ptr,v+off,len);
446 
447           /* keep track of number of zeros */
448 	  if (dim)
449 	    {
450 	      for (k=0; k<len; k++)
451 		{
452 		  if (x_ptr[k]==0.0)
453 		    {xxt_zero_nnz++;}
454 		}
455 	    }
456 	  else
457 	    {
458 	      for (k=0; k<len; k++)
459 		{
460 		  if (x_ptr[k]==0.0)
461 		    {xxt_zero_nnz_0++;}
462 		}
463 	    }
464 	  col_indices[2*i] = off;
465 	  col_sz[i] = col_indices[2*i+1] = len;
466 	  col_vals[i] = x_ptr;
467 	}
468       else
469 	{
470 	  col_indices[2*i] = 0;
471 	  col_sz[i] = col_indices[2*i+1] = 0;
472 	  col_vals[i] = x_ptr;
473 	}
474     }
475 
476   /* close off stages for execution phase */
477   while (dim!=level)
478     {
479       stages[dim++]=i;
480       ierr = PetscInfo2(0,"disconnected!!! dim(%d)!=level(%d)\n",dim,level);CHKERRQ(ierr);
481     }
482   stages[dim]=i;
483 
484   xxt_handle->info->n=xxt_handle->mvi->n;
485   xxt_handle->info->m=m;
486   xxt_handle->info->nnz=xxt_nnz;
487   xxt_handle->info->max_nnz=xxt_max_nnz;
488   xxt_handle->info->msg_buf_sz=stages[level]-stages[0];
489   xxt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
490   xxt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
491   xxt_handle->info->x=x;
492   xxt_handle->info->col_vals=col_vals;
493   xxt_handle->info->col_sz=col_sz;
494   xxt_handle->info->col_indices=col_indices;
495   xxt_handle->info->stages=stages;
496   xxt_handle->info->nsolves=0;
497   xxt_handle->info->tot_solve_time=0.0;
498 
499   free(segs);
500   free(u);
501   free(v);
502   free(uu);
503   free(z);
504   free(w);
505 
506   return(0);
507 }
508 
509 /**************************************xxt.c***********************************/
510 static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle,  PetscScalar *uc)
511 {
512    PetscInt off, len, *iptr;
513   PetscInt level       =xxt_handle->level;
514   PetscInt n           =xxt_handle->info->n;
515   PetscInt m           =xxt_handle->info->m;
516   PetscInt *stages     =xxt_handle->info->stages;
517   PetscInt *col_indices=xxt_handle->info->col_indices;
518   PetscScalar *x_ptr, *uu_ptr;
519   PetscScalar *solve_uu=xxt_handle->info->solve_uu;
520   PetscScalar *solve_w =xxt_handle->info->solve_w;
521   PetscScalar *x       =xxt_handle->info->x;
522   PetscBLASInt i1 = 1;
523 
524   PetscFunctionBegin;
525   uu_ptr=solve_uu;
526   rvec_zero(uu_ptr,m);
527 
528   /* x  = X.Y^T.b */
529   /* uu = Y^T.b */
530   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
531     {
532       off=*iptr++; len=*iptr++;
533       *uu_ptr++ = BLASdot_(&len,uc+off,&i1,x_ptr,&i1);
534     }
535 
536   /* comunication of beta */
537   uu_ptr=solve_uu;
538   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}
539 
540   rvec_zero(uc,n);
541 
542   /* x = X.uu */
543   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
544     {
545       off=*iptr++; len=*iptr++;
546       BLASaxpy_(&len,uu_ptr++,x_ptr,&i1,uc+off,&i1);
547     }
548   PetscFunctionReturn(0);
549 }
550 
551 /**************************************xxt.c***********************************/
552 static PetscErrorCode check_handle(xxt_ADT xxt_handle)
553 {
554   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
555 
556   PetscFunctionBegin;
557   if (xxt_handle==NULL)
558     {SETERRQ1(PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %d\n",xxt_handle);}
559 
560   vals[0]=vals[1]=xxt_handle->id;
561   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
562   if ((vals[0]!=vals[1])||(xxt_handle->id<=0))
563     {SETERRQ3(PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %d/%d %d\n",vals[0],vals[1], xxt_handle->id);}
564   PetscFunctionReturn(0);
565 }
566 
567 /**************************************xxt.c***********************************/
568 static  PetscErrorCode det_separators(xxt_ADT xxt_handle)
569 {
570   PetscInt i, ct, id;
571   PetscInt mask, edge, *iptr;
572   PetscInt *dir, *used;
573   PetscInt sum[4], w[4];
574   PetscScalar rsum[4], rw[4];
575   PetscInt op[] = {GL_ADD,0};
576   PetscScalar *lhs, *rhs;
577   PetscInt *nsep, *lnsep, *fo, nfo=0;
578   gs_ADT gs_handle=xxt_handle->mvi->gs_handle;
579   PetscInt *local2global=xxt_handle->mvi->local2global;
580   PetscInt  n=xxt_handle->mvi->n;
581   PetscInt  m=xxt_handle->mvi->m;
582   PetscInt level=xxt_handle->level;
583   PetscInt shared=FALSE;
584 
585   PetscFunctionBegin;
586   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
587   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
588   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
589   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
590   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
591 
592   ivec_zero(dir  ,level+1);
593   ivec_zero(nsep ,level+1);
594   ivec_zero(lnsep,level+1);
595   ivec_set (fo   ,-1,n+1);
596   ivec_zero(used,n);
597 
598   lhs  = (double*)malloc(sizeof(PetscScalar)*m);
599   rhs  = (double*)malloc(sizeof(PetscScalar)*m);
600 
601   /* determine the # of unique dof */
602   rvec_zero(lhs,m);
603   rvec_set(lhs,1.0,n);
604   gs_gop_hc(gs_handle,lhs,"+\0",level);
605   rvec_zero(rsum,2);
606   for (ct=i=0;i<n;i++)
607     {
608       if (lhs[i]!=0.0)
609 	{rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
610     }
611   grop_hc(rsum,rw,2,op,level);
612   rsum[0]+=0.1;
613   rsum[1]+=0.1;
614   /*  if (!my_id)
615     {
616       printf("xxt n unique = %d (%g)\n",(int) rsum[0], rsum[0]);
617       printf("xxt n shared = %d (%g)\n",(int) rsum[1], rsum[1]);
618       }*/
619 
620   if (fabs(rsum[0]-rsum[1])>EPS)
621     {shared=TRUE;}
622 
623   xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0];
624   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0];
625 
626   /* determine separator sets top down */
627   if (shared)
628     {
629       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
630 	{
631 	  /* set rsh of hc, fire, and collect lhs responses */
632 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
633 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
634 
635 	  /* set lsh of hc, fire, and collect rhs responses */
636 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
637 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
638 
639 	  for (i=0;i<n;i++)
640 	    {
641 	      if (id< mask)
642 		{
643 		  if (lhs[i]!=0.0)
644 		    {lhs[i]=1.0;}
645 		}
646 	      if (id>=mask)
647 		{
648 		  if (rhs[i]!=0.0)
649 		    {rhs[i]=1.0;}
650 		}
651 	    }
652 
653 	  if (id< mask)
654 	    {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
655 	  else
656 	    {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}
657 
658 	  /* count number of dofs I own that have signal and not in sep set */
659 	  rvec_zero(rsum,4);
660 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
661 	    {
662 	      if (!used[i])
663 		{
664 		  /* number of unmarked dofs on node */
665 		  ct++;
666 		  /* number of dofs to be marked on lhs hc */
667 		  if (id< mask)
668 		    {
669 		      if (lhs[i]!=0.0)
670 			{sum[0]++; rsum[0]+=1.0/lhs[i];}
671 		    }
672 		  /* number of dofs to be marked on rhs hc */
673 		  if (id>=mask)
674 		    {
675 		      if (rhs[i]!=0.0)
676 			{sum[1]++; rsum[1]+=1.0/rhs[i];}
677 		    }
678 		}
679 	    }
680 
681 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
682 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
683 	  (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
684 	  giop_hc(sum,w,4,op,edge);
685 	  grop_hc(rsum,rw,4,op,edge);
686 	  rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
687 
688 	  if (id<mask)
689 	    {
690 	      /* mark dofs I own that have signal and not in sep set */
691 	      for (ct=i=0;i<n;i++)
692 		{
693 		  if ((!used[i])&&(lhs[i]!=0.0))
694 		    {
695 		      ct++; nfo++;
696 
697 		      if (nfo>n)
698 			{SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}
699 
700 		      *--iptr = local2global[i];
701 		      used[i]=edge;
702 		    }
703 		}
704 	      if (ct>1) {ivec_sort(iptr,ct);}
705 
706 	      lnsep[edge]=ct;
707 	      nsep[edge]=(PetscInt) rsum[0];
708 	      dir [edge]=LEFT;
709 	    }
710 
711 	  if (id>=mask)
712 	    {
713 	      /* mark dofs I own that have signal and not in sep set */
714 	      for (ct=i=0;i<n;i++)
715 		{
716 		  if ((!used[i])&&(rhs[i]!=0.0))
717 		    {
718 		      ct++; nfo++;
719 
720 		      if (nfo>n)
721 			{SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}
722 
723 		      *--iptr = local2global[i];
724 		      used[i]=edge;
725 		    }
726 		}
727 	      if (ct>1) {ivec_sort(iptr,ct);}
728 
729 	      lnsep[edge]=ct;
730 	      nsep[edge]= (PetscInt) rsum[1];
731 	      dir [edge]=RIGHT;
732 	    }
733 
734 	  /* LATER or we can recur on these to order seps at this level */
735 	  /* do we need full set of separators for this?                */
736 
737 	  /* fold rhs hc into lower */
738 	  if (id>=mask)
739 	    {id-=mask;}
740 	}
741     }
742   else
743     {
744       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
745 	{
746 	  /* set rsh of hc, fire, and collect lhs responses */
747 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
748 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
749 
750 	  /* set lsh of hc, fire, and collect rhs responses */
751 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
752 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
753 
754 	  /* count number of dofs I own that have signal and not in sep set */
755 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
756 	    {
757 	      if (!used[i])
758 		{
759 		  /* number of unmarked dofs on node */
760 		  ct++;
761 		  /* number of dofs to be marked on lhs hc */
762 		  if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
763 		  /* number of dofs to be marked on rhs hc */
764 		  if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
765 		}
766 	    }
767 
768 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
769 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
770 	  giop_hc(sum,w,4,op,edge);
771 
772 	  /* lhs hc wins */
773 	  if (sum[2]>=sum[3])
774 	    {
775 	      if (id<mask)
776 		{
777 		  /* mark dofs I own that have signal and not in sep set */
778 		  for (ct=i=0;i<n;i++)
779 		    {
780 		      if ((!used[i])&&(lhs[i]!=0.0))
781 			{
782 			  ct++; nfo++;
783 			  *--iptr = local2global[i];
784 			  used[i]=edge;
785 			}
786 		    }
787 		  if (ct>1) {ivec_sort(iptr,ct);}
788 		  lnsep[edge]=ct;
789 		}
790 	      nsep[edge]=sum[0];
791 	      dir [edge]=LEFT;
792 	    }
793 	  /* rhs hc wins */
794 	  else
795 	    {
796 	      if (id>=mask)
797 		{
798 		  /* mark dofs I own that have signal and not in sep set */
799 		  for (ct=i=0;i<n;i++)
800 		    {
801 		      if ((!used[i])&&(rhs[i]!=0.0))
802 			{
803 			  ct++; nfo++;
804 			  *--iptr = local2global[i];
805 			  used[i]=edge;
806 			}
807 		    }
808 		  if (ct>1) {ivec_sort(iptr,ct);}
809 		  lnsep[edge]=ct;
810 		}
811 	      nsep[edge]=sum[1];
812 	      dir [edge]=RIGHT;
813 	    }
814 	  /* LATER or we can recur on these to order seps at this level */
815 	  /* do we need full set of separators for this?                */
816 
817 	  /* fold rhs hc into lower */
818 	  if (id>=mask)
819 	    {id-=mask;}
820 	}
821     }
822 
823 
824   /* level 0 is on processor case - so mark the remainder */
825   for (ct=i=0;i<n;i++)
826     {
827       if (!used[i])
828 	{
829 	  ct++; nfo++;
830 	  *--iptr = local2global[i];
831 	  used[i]=edge;
832 	}
833     }
834   if (ct>1) {ivec_sort(iptr,ct);}
835   lnsep[edge]=ct;
836   nsep [edge]=ct;
837   dir  [edge]=LEFT;
838 
839   xxt_handle->info->nsep=nsep;
840   xxt_handle->info->lnsep=lnsep;
841   xxt_handle->info->fo=fo;
842   xxt_handle->info->nfo=nfo;
843 
844   free(dir);
845   free(lhs);
846   free(rhs);
847   free(used);
848   PetscFunctionReturn(0);
849 }
850 
851 /**************************************xxt.c***********************************/
852 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
853 {
854   mv_info *mvi;
855 
856 
857   mvi = (mv_info*)malloc(sizeof(mv_info));
858   mvi->n=n;
859   mvi->m=m;
860   mvi->n_global=-1;
861   mvi->m_global=-1;
862   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
863   ivec_copy(mvi->local2global,local2global,m);
864   mvi->local2global[m] = INT_MAX;
865   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
866   mvi->grid_data=grid_data;
867 
868   /* set xxt communication handle to perform restricted matvec */
869   mvi->gs_handle = gs_init(local2global, m, num_nodes);
870 
871   return(mvi);
872 }
873 
874 /**************************************xxt.c***********************************/
875 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
876 {
877   PetscFunctionBegin;
878   A->matvec((mv_info*)A->grid_data,v,u);
879   PetscFunctionReturn(0);
880 }
881 
882 
883 
884