xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision e7e72b3d0edcd0d15e7f68c03be08666507fc872)
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_COMM_SELF,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, PetscScalar *x, PetscScalar *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   PetscErrorCode ierr;
163 
164   comm_init();
165   check_handle(xxt_handle);
166 
167   /* if factorization not done there are no stats */
168   if (!xxt_handle->info||!xxt_handle->mvi)
169     {
170       if (!my_id)
171 	{ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");}
172       return 1;
173     }
174 
175   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
176   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
177   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
178   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
179 
180   fvals[0]=fvals[1]=fvals[2]
181     =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
182   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
183 
184   if (!my_id)
185     {
186       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_nnz=%D\n",my_id,vals[0]);
187       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_nnz=%D\n",my_id,vals[1]);
188       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
189       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_nnz=%D\n",my_id,vals[2]);
190       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
191       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
192       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_n  =%D\n",my_id,vals[3]);
193       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_n  =%D\n",my_id,vals[4]);
194       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
195       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_n  =%D\n",my_id,vals[5]);
196       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_buf=%D\n",my_id,vals[6]);
197       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_buf=%D\n",my_id,vals[7]);
198       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
199       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_slv=%g\n",my_id,fvals[0]);
200       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_slv=%g\n",my_id,fvals[1]);
201       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_slv=%g\n",my_id,fvals[2]/num_nodes);
202     }
203 
204   return(0);
205 }
206 
207 /*************************************xxt.c************************************
208 
209 Description: get A_local, local portion of global coarse matrix which
210 is a row dist. nxm matrix w/ n<m.
211    o my_ml holds address of ML struct associated w/A_local and coarse grid
212    o local2global holds global number of column i (i=0,...,m-1)
213    o local2global holds global number of row    i (i=0,...,n-1)
214    o mylocmatvec performs A_local . vec_local (note that gs is performed using
215    gs_init/gop).
216 
217 mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
218 mylocmatvec (void :: void *data, double *in, double *out)
219 **************************************xxt.c***********************************/
220 static PetscInt do_xxt_factor(xxt_ADT xxt_handle)
221 {
222   return xxt_generate(xxt_handle);
223 }
224 
225 /**************************************xxt.c***********************************/
226 static PetscInt xxt_generate(xxt_ADT xxt_handle)
227 {
228   PetscInt i,j,k,idex;
229   PetscInt dim, col;
230   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
231   PetscInt *segs;
232   PetscInt op[] = {GL_ADD,0};
233   PetscInt off, len;
234   PetscScalar *x_ptr;
235   PetscInt *iptr, flag;
236   PetscInt start=0, end, work;
237   PetscInt op2[] = {GL_MIN,0};
238   gs_ADT gs_handle;
239   PetscInt *nsep, *lnsep, *fo;
240   PetscInt a_n=xxt_handle->mvi->n;
241   PetscInt a_m=xxt_handle->mvi->m;
242   PetscInt *a_local2global=xxt_handle->mvi->local2global;
243   PetscInt level;
244   PetscInt xxt_nnz=0, xxt_max_nnz=0;
245   PetscInt n, m;
246   PetscInt *col_sz, *col_indices, *stages;
247   PetscScalar **col_vals, *x;
248   PetscInt n_global;
249   PetscInt xxt_zero_nnz=0;
250   PetscInt xxt_zero_nnz_0=0;
251   PetscBLASInt i1 = 1,dlen;
252   PetscScalar dm1 = -1.0;
253   PetscErrorCode ierr;
254 
255   n=xxt_handle->mvi->n;
256   nsep=xxt_handle->info->nsep;
257   lnsep=xxt_handle->info->lnsep;
258   fo=xxt_handle->info->fo;
259   end=lnsep[0];
260   level=xxt_handle->level;
261   gs_handle=xxt_handle->mvi->gs_handle;
262 
263   /* is there a null space? */
264   /* LATER add in ability to detect null space by checking alpha */
265   for (i=0, j=0; i<=level; i++)
266     {j+=nsep[i];}
267 
268   m = j-xxt_handle->ns;
269   if (m!=j)
270     {ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);}
271 
272   /* get and initialize storage for x local         */
273   /* note that x local is nxm and stored by columns */
274   col_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
275   col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
276   col_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
277   for (i=j=0; i<m; i++, j+=2)
278     {
279       col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
280       col_vals[i] = NULL;
281     }
282   col_indices[j]=-1;
283 
284   /* size of separators for each sub-hc working from bottom of tree to top */
285   /* this looks like nsep[]=segments */
286   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
287   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
288   ivec_zero(stages,level+1);
289   ivec_copy(segs,nsep,level+1);
290   for (i=0; i<level; i++)
291     {segs[i+1] += segs[i];}
292   stages[0] = segs[0];
293 
294   /* temporary vectors  */
295   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
296   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
297   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
298   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
299   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
300 
301   /* extra nnz due to replication of vertices across separators */
302   for (i=1, j=0; i<=level; i++)
303     {j+=nsep[i];}
304 
305   /* storage for sparse x values */
306   n_global = xxt_handle->info->n_global;
307   xxt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
308   x = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
309   xxt_nnz = 0;
310 
311   /* LATER - can embed next sep to fire in gs */
312   /* time to make the donuts - generate X factor */
313   for (dim=i=j=0;i<m;i++)
314     {
315       /* time to move to the next level? */
316       while (i==segs[dim]) {
317         if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n");
318         stages[dim++]=i;
319         end+=lnsep[dim];
320       }
321       stages[dim]=i;
322 
323       /* which column are we firing? */
324       /* i.e. set v_l */
325       /* use new seps and do global min across hc to determine which one to fire */
326       (start<end) ? (col=fo[start]) : (col=INT_MAX);
327       giop_hc(&col,&work,1,op2,dim);
328 
329       /* shouldn't need this */
330       if (col==INT_MAX)
331 	{
332 	  ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
333 	  continue;
334 	}
335 
336       /* do I own it? I should */
337       rvec_zero(v ,a_m);
338       if (col==fo[start])
339 	{
340 	  start++;
341 	  idex=ivec_linear_search(col, a_local2global, a_n);
342 	  if (idex!=-1){
343             v[idex] = 1.0; j++;
344           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");
345 	}
346       else
347 	{
348 	  idex=ivec_linear_search(col, a_local2global, a_m);
349 	  if (idex!=-1)
350 	    {v[idex] = 1.0;}
351 	}
352 
353       /* perform u = A.v_l */
354       rvec_zero(u,n);
355       do_matvec(xxt_handle->mvi,v,u);
356 
357       /* uu =  X^T.u_l (local portion) */
358       /* technically only need to zero out first i entries */
359       /* later turn this into an XXT_solve call ? */
360       rvec_zero(uu,m);
361       x_ptr=x;
362       iptr = col_indices;
363       for (k=0; k<i; k++)
364 	{
365 	  off = *iptr++;
366 	  len = *iptr++;
367           dlen = PetscBLASIntCast(len);
368 	  uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1);
369 	  x_ptr+=len;
370 	}
371 
372 
373       /* uu = X^T.u_l (comm portion) */
374       ssgl_radd  (uu, w, dim, stages);
375 
376       /* z = X.uu */
377       rvec_zero(z,n);
378       x_ptr=x;
379       iptr = col_indices;
380       for (k=0; k<i; k++)
381 	{
382 	  off = *iptr++;
383 	  len = *iptr++;
384           dlen = PetscBLASIntCast(len);
385 	  BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
386 	  x_ptr+=len;
387 	}
388 
389       /* compute v_l = v_l - z */
390       rvec_zero(v+a_n,a_m-a_n);
391       dlen = PetscBLASIntCast(n);
392       BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);
393 
394       /* compute u_l = A.v_l */
395       if (a_n!=a_m)
396 	{gs_gop_hc(gs_handle,v,"+\0",dim);}
397       rvec_zero(u,n);
398       do_matvec(xxt_handle->mvi,v,u);
399 
400       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
401       dlen = PetscBLASIntCast(n);
402       alpha = BLASdot_(&dlen,u,&i1,v,&i1);
403       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
404       grop_hc(&alpha, &alpha_w, 1, op, dim);
405 
406       alpha = (PetscScalar) sqrt((double)alpha);
407 
408       /* check for small alpha                             */
409       /* LATER use this to detect and determine null space */
410       if (fabs(alpha)<1.0e-14)
411 	{SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);}
412 
413       /* compute v_l = v_l/sqrt(alpha) */
414       rvec_scale(v,1.0/alpha,n);
415 
416       /* add newly generated column, v_l, to X */
417       flag = 1;
418       off=len=0;
419       for (k=0; k<n; k++)
420 	{
421 	  if (v[k]!=0.0)
422 	    {
423 	      len=k;
424 	      if (flag)
425 		{off=k; flag=0;}
426 	    }
427 	}
428 
429       len -= (off-1);
430 
431       if (len>0)
432 	{
433 	  if ((xxt_nnz+len)>xxt_max_nnz)
434 	    {
435 	      ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
436 	      xxt_max_nnz *= 2;
437 	      x_ptr = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
438 	      rvec_copy(x_ptr,x,xxt_nnz);
439 	      free(x);
440 	      x = x_ptr;
441 	      x_ptr+=xxt_nnz;
442 	    }
443 	  xxt_nnz += len;
444 	  rvec_copy(x_ptr,v+off,len);
445 
446           /* keep track of number of zeros */
447 	  if (dim)
448 	    {
449 	      for (k=0; k<len; k++)
450 		{
451 		  if (x_ptr[k]==0.0)
452 		    {xxt_zero_nnz++;}
453 		}
454 	    }
455 	  else
456 	    {
457 	      for (k=0; k<len; k++)
458 		{
459 		  if (x_ptr[k]==0.0)
460 		    {xxt_zero_nnz_0++;}
461 		}
462 	    }
463 	  col_indices[2*i] = off;
464 	  col_sz[i] = col_indices[2*i+1] = len;
465 	  col_vals[i] = x_ptr;
466 	}
467       else
468 	{
469 	  col_indices[2*i] = 0;
470 	  col_sz[i] = col_indices[2*i+1] = 0;
471 	  col_vals[i] = x_ptr;
472 	}
473     }
474 
475   /* close off stages for execution phase */
476   while (dim!=level)
477     {
478       stages[dim++]=i;
479       ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
480     }
481   stages[dim]=i;
482 
483   xxt_handle->info->n=xxt_handle->mvi->n;
484   xxt_handle->info->m=m;
485   xxt_handle->info->nnz=xxt_nnz;
486   xxt_handle->info->max_nnz=xxt_max_nnz;
487   xxt_handle->info->msg_buf_sz=stages[level]-stages[0];
488   xxt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
489   xxt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
490   xxt_handle->info->x=x;
491   xxt_handle->info->col_vals=col_vals;
492   xxt_handle->info->col_sz=col_sz;
493   xxt_handle->info->col_indices=col_indices;
494   xxt_handle->info->stages=stages;
495   xxt_handle->info->nsolves=0;
496   xxt_handle->info->tot_solve_time=0.0;
497 
498   free(segs);
499   free(u);
500   free(v);
501   free(uu);
502   free(z);
503   free(w);
504 
505   return(0);
506 }
507 
508 /**************************************xxt.c***********************************/
509 static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle,  PetscScalar *uc)
510 {
511    PetscInt off, len, *iptr;
512   PetscInt level       =xxt_handle->level;
513   PetscInt n           =xxt_handle->info->n;
514   PetscInt m           =xxt_handle->info->m;
515   PetscInt *stages     =xxt_handle->info->stages;
516   PetscInt *col_indices=xxt_handle->info->col_indices;
517   PetscScalar *x_ptr, *uu_ptr;
518   PetscScalar *solve_uu=xxt_handle->info->solve_uu;
519   PetscScalar *solve_w =xxt_handle->info->solve_w;
520   PetscScalar *x       =xxt_handle->info->x;
521   PetscBLASInt i1 = 1,dlen;
522 
523   PetscFunctionBegin;
524   uu_ptr=solve_uu;
525   rvec_zero(uu_ptr,m);
526 
527   /* x  = X.Y^T.b */
528   /* uu = Y^T.b */
529   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
530     {
531       off=*iptr++; len=*iptr++;
532       dlen = PetscBLASIntCast(len);
533       *uu_ptr++ = BLASdot_(&dlen,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       dlen = PetscBLASIntCast(len);
547       BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
548     }
549   PetscFunctionReturn(0);
550 }
551 
552 /**************************************xxt.c***********************************/
553 static PetscErrorCode check_handle(xxt_ADT xxt_handle)
554 {
555   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
556 
557   PetscFunctionBegin;
558   if (!xxt_handle) SETERRQ1(PETSC_COMM_SELF,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)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %D/%D %D\n",vals[0],vals[1], xxt_handle->id);
563   PetscFunctionReturn(0);
564 }
565 
566 /**************************************xxt.c***********************************/
567 static  PetscErrorCode det_separators(xxt_ADT xxt_handle)
568 {
569   PetscInt i, ct, id;
570   PetscInt mask, edge, *iptr;
571   PetscInt *dir, *used;
572   PetscInt sum[4], w[4];
573   PetscScalar rsum[4], rw[4];
574   PetscInt op[] = {GL_ADD,0};
575   PetscScalar *lhs, *rhs;
576   PetscInt *nsep, *lnsep, *fo, nfo=0;
577   gs_ADT gs_handle=xxt_handle->mvi->gs_handle;
578   PetscInt *local2global=xxt_handle->mvi->local2global;
579   PetscInt  n=xxt_handle->mvi->n;
580   PetscInt  m=xxt_handle->mvi->m;
581   PetscInt level=xxt_handle->level;
582   PetscInt shared=FALSE;
583 
584   PetscFunctionBegin;
585   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
586   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
587   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
588   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
589   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
590 
591   ivec_zero(dir  ,level+1);
592   ivec_zero(nsep ,level+1);
593   ivec_zero(lnsep,level+1);
594   ivec_set (fo   ,-1,n+1);
595   ivec_zero(used,n);
596 
597   lhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
598   rhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
599 
600   /* determine the # of unique dof */
601   rvec_zero(lhs,m);
602   rvec_set(lhs,1.0,n);
603   gs_gop_hc(gs_handle,lhs,"+\0",level);
604   rvec_zero(rsum,2);
605   for (ct=i=0;i<n;i++)
606     {
607       if (lhs[i]!=0.0)
608 	{rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
609     }
610   grop_hc(rsum,rw,2,op,level);
611   rsum[0]+=0.1;
612   rsum[1]+=0.1;
613 
614   if (fabs(rsum[0]-rsum[1])>EPS)
615     {shared=TRUE;}
616 
617   xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0];
618   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0];
619 
620   /* determine separator sets top down */
621   if (shared)
622     {
623       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
624 	{
625 	  /* set rsh of hc, fire, and collect lhs responses */
626 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
627 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
628 
629 	  /* set lsh of hc, fire, and collect rhs responses */
630 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
631 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
632 
633 	  for (i=0;i<n;i++)
634 	    {
635 	      if (id< mask)
636 		{
637 		  if (lhs[i]!=0.0)
638 		    {lhs[i]=1.0;}
639 		}
640 	      if (id>=mask)
641 		{
642 		  if (rhs[i]!=0.0)
643 		    {rhs[i]=1.0;}
644 		}
645 	    }
646 
647 	  if (id< mask)
648 	    {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
649 	  else
650 	    {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}
651 
652 	  /* count number of dofs I own that have signal and not in sep set */
653 	  rvec_zero(rsum,4);
654 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
655 	    {
656 	      if (!used[i])
657 		{
658 		  /* number of unmarked dofs on node */
659 		  ct++;
660 		  /* number of dofs to be marked on lhs hc */
661 		  if (id< mask)
662 		    {
663 		      if (lhs[i]!=0.0)
664 			{sum[0]++; rsum[0]+=1.0/lhs[i];}
665 		    }
666 		  /* number of dofs to be marked on rhs hc */
667 		  if (id>=mask)
668 		    {
669 		      if (rhs[i]!=0.0)
670 			{sum[1]++; rsum[1]+=1.0/rhs[i];}
671 		    }
672 		}
673 	    }
674 
675 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
676 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
677 	  (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
678 	  giop_hc(sum,w,4,op,edge);
679 	  grop_hc(rsum,rw,4,op,edge);
680 	  rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
681 
682 	  if (id<mask)
683 	    {
684 	      /* mark dofs I own that have signal and not in sep set */
685 	      for (ct=i=0;i<n;i++)
686 		{
687 		  if ((!used[i])&&(lhs[i]!=0.0))
688 		    {
689 		      ct++; nfo++;
690 
691 		      if (nfo>n)
692 			{SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");}
693 
694 		      *--iptr = local2global[i];
695 		      used[i]=edge;
696 		    }
697 		}
698 	      if (ct>1) {ivec_sort(iptr,ct);}
699 
700 	      lnsep[edge]=ct;
701 	      nsep[edge]=(PetscInt) rsum[0];
702 	      dir [edge]=LEFT;
703 	    }
704 
705 	  if (id>=mask)
706 	    {
707 	      /* mark dofs I own that have signal and not in sep set */
708 	      for (ct=i=0;i<n;i++)
709 		{
710 		  if ((!used[i])&&(rhs[i]!=0.0))
711 		    {
712 		      ct++; nfo++;
713 
714 		      if (nfo>n)
715 			{SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");}
716 
717 		      *--iptr = local2global[i];
718 		      used[i]=edge;
719 		    }
720 		}
721 	      if (ct>1) {ivec_sort(iptr,ct);}
722 
723 	      lnsep[edge]=ct;
724 	      nsep[edge]= (PetscInt) rsum[1];
725 	      dir [edge]=RIGHT;
726 	    }
727 
728 	  /* LATER or we can recur on these to order seps at this level */
729 	  /* do we need full set of separators for this?                */
730 
731 	  /* fold rhs hc into lower */
732 	  if (id>=mask)
733 	    {id-=mask;}
734 	}
735     }
736   else
737     {
738       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
739 	{
740 	  /* set rsh of hc, fire, and collect lhs responses */
741 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
742 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
743 
744 	  /* set lsh of hc, fire, and collect rhs responses */
745 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
746 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
747 
748 	  /* count number of dofs I own that have signal and not in sep set */
749 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
750 	    {
751 	      if (!used[i])
752 		{
753 		  /* number of unmarked dofs on node */
754 		  ct++;
755 		  /* number of dofs to be marked on lhs hc */
756 		  if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
757 		  /* number of dofs to be marked on rhs hc */
758 		  if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
759 		}
760 	    }
761 
762 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
763 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
764 	  giop_hc(sum,w,4,op,edge);
765 
766 	  /* lhs hc wins */
767 	  if (sum[2]>=sum[3])
768 	    {
769 	      if (id<mask)
770 		{
771 		  /* mark dofs I own that have signal and not in sep set */
772 		  for (ct=i=0;i<n;i++)
773 		    {
774 		      if ((!used[i])&&(lhs[i]!=0.0))
775 			{
776 			  ct++; nfo++;
777 			  *--iptr = local2global[i];
778 			  used[i]=edge;
779 			}
780 		    }
781 		  if (ct>1) {ivec_sort(iptr,ct);}
782 		  lnsep[edge]=ct;
783 		}
784 	      nsep[edge]=sum[0];
785 	      dir [edge]=LEFT;
786 	    }
787 	  /* rhs hc wins */
788 	  else
789 	    {
790 	      if (id>=mask)
791 		{
792 		  /* mark dofs I own that have signal and not in sep set */
793 		  for (ct=i=0;i<n;i++)
794 		    {
795 		      if ((!used[i])&&(rhs[i]!=0.0))
796 			{
797 			  ct++; nfo++;
798 			  *--iptr = local2global[i];
799 			  used[i]=edge;
800 			}
801 		    }
802 		  if (ct>1) {ivec_sort(iptr,ct);}
803 		  lnsep[edge]=ct;
804 		}
805 	      nsep[edge]=sum[1];
806 	      dir [edge]=RIGHT;
807 	    }
808 	  /* LATER or we can recur on these to order seps at this level */
809 	  /* do we need full set of separators for this?                */
810 
811 	  /* fold rhs hc into lower */
812 	  if (id>=mask)
813 	    {id-=mask;}
814 	}
815     }
816 
817 
818   /* level 0 is on processor case - so mark the remainder */
819   for (ct=i=0;i<n;i++)
820     {
821       if (!used[i])
822 	{
823 	  ct++; nfo++;
824 	  *--iptr = local2global[i];
825 	  used[i]=edge;
826 	}
827     }
828   if (ct>1) {ivec_sort(iptr,ct);}
829   lnsep[edge]=ct;
830   nsep [edge]=ct;
831   dir  [edge]=LEFT;
832 
833   xxt_handle->info->nsep=nsep;
834   xxt_handle->info->lnsep=lnsep;
835   xxt_handle->info->fo=fo;
836   xxt_handle->info->nfo=nfo;
837 
838   free(dir);
839   free(lhs);
840   free(rhs);
841   free(used);
842   PetscFunctionReturn(0);
843 }
844 
845 /**************************************xxt.c***********************************/
846 static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
847 {
848   mv_info *mvi;
849 
850 
851   mvi = (mv_info*)malloc(sizeof(mv_info));
852   mvi->n=n;
853   mvi->m=m;
854   mvi->n_global=-1;
855   mvi->m_global=-1;
856   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
857   ivec_copy(mvi->local2global,local2global,m);
858   mvi->local2global[m] = INT_MAX;
859   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
860   mvi->grid_data=grid_data;
861 
862   /* set xxt communication handle to perform restricted matvec */
863   mvi->gs_handle = gs_init(local2global, m, num_nodes);
864 
865   return(mvi);
866 }
867 
868 /**************************************xxt.c***********************************/
869 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
870 {
871   PetscFunctionBegin;
872   A->matvec((mv_info*)A->grid_data,v,u);
873   PetscFunctionReturn(0);
874 }
875 
876 
877 
878