xref: /petsc/src/ksp/pc/impls/tfs/xyt.c (revision d9fead3d4a4aa17f1e85c6884ab255f226095a30)
1 #define PETSCKSP_DLL
2 
3 /*************************************xyt.c************************************
4 Module Name: xyt
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 **************************************xyt.c***********************************/
20 
21 
22 /*************************************xyt.c************************************
23 NOTES ON USAGE:
24 
25 **************************************xyt.c***********************************/
26 #include "src/ksp/pc/impls/tfs/tfs.h"
27 
28 #define LEFT  -1
29 #define RIGHT  1
30 #define BOTH   0
31 
32 typedef struct xyt_solver_info {
33   int n, m, n_global, m_global;
34   int nnz, max_nnz, msg_buf_sz;
35   int *nsep, *lnsep, *fo, nfo, *stages;
36   int *xcol_sz, *xcol_indices;
37   PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
38   int *ycol_sz, *ycol_indices;
39   PetscScalar **ycol_vals, *y;
40   int nsolves;
41   PetscScalar tot_solve_time;
42 } xyt_info;
43 
44 
45 typedef struct matvec_info {
46   int n, m, n_global, m_global;
47   int *local2global;
48   gs_ADT gs_handle;
49   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
50   void *grid_data;
51 } mv_info;
52 
53 struct xyt_CDT{
54   int id;
55   int ns;
56   int level;
57   xyt_info *info;
58   mv_info  *mvi;
59 };
60 
61 static int n_xyt=0;
62 static int n_xyt_handles=0;
63 
64 /* prototypes */
65 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
66 static PetscErrorCode check_handle(xyt_ADT xyt_handle);
67 static PetscErrorCode det_separators(xyt_ADT xyt_handle);
68 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
69 static int xyt_generate(xyt_ADT xyt_handle);
70 static int do_xyt_factor(xyt_ADT xyt_handle);
71 static mv_info *set_mvi(int *local2global, int n, int m, void *matvec, void *grid_data);
72 
73 
74 /*************************************xyt.c************************************
75 Function: XYT_new()
76 
77 Input :
78 Output:
79 Return:
80 Description:
81 **************************************xyt.c***********************************/
82 xyt_ADT
83 XYT_new(void)
84 {
85   xyt_ADT xyt_handle;
86 
87 
88 
89   /* rolling count on n_xyt ... pot. problem here */
90   n_xyt_handles++;
91   xyt_handle       = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
92   xyt_handle->id   = ++n_xyt;
93   xyt_handle->info = NULL;
94   xyt_handle->mvi  = NULL;
95 
96   return(xyt_handle);
97 }
98 
99 
100 /*************************************xyt.c************************************
101 Function: XYT_factor()
102 
103 Input :
104 Output:
105 Return:
106 Description:
107 **************************************xyt.c***********************************/
108 int
109 XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt  handle */
110 	   int *local2global,  /* global column mapping       */
111 	   int n,              /* local num rows              */
112 	   int m,              /* local num cols              */
113 	   void *matvec,       /* b_loc=A_local.x_loc         */
114 	   void *grid_data     /* grid data for matvec        */
115 	   )
116 {
117 
118   comm_init();
119   check_handle(xyt_handle);
120 
121   /* only 2^k for now and all nodes participating */
122   if ((1<<(xyt_handle->level=i_log2_num_nodes))!=num_nodes)
123     {error_msg_fatal("only 2^k for now and MPI_COMM_WORLD!!! %d != %d\n",1<<i_log2_num_nodes,num_nodes);}
124 
125   /* space for X info */
126   xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));
127 
128   /* set up matvec handles */
129   xyt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);
130 
131   /* matrix is assumed to be of full rank */
132   /* LATER we can reset to indicate rank def. */
133   xyt_handle->ns=0;
134 
135   /* determine separators and generate firing order - NB xyt info set here */
136   det_separators(xyt_handle);
137 
138   return(do_xyt_factor(xyt_handle));
139 }
140 
141 
142 /*************************************xyt.c************************************
143 Function: XYT_solve
144 
145 Input :
146 Output:
147 Return:
148 Description:
149 **************************************xyt.c***********************************/
150 int
151 XYT_solve(xyt_ADT xyt_handle, double *x, double *b)
152 {
153   comm_init();
154   check_handle(xyt_handle);
155 
156   /* need to copy b into x? */
157   if (b)
158     {rvec_copy(x,b,xyt_handle->mvi->n);}
159   do_xyt_solve(xyt_handle,x);
160 
161   return(0);
162 }
163 
164 
165 /*************************************xyt.c************************************
166 Function: XYT_free()
167 
168 Input :
169 Output:
170 Return:
171 Description:
172 **************************************xyt.c***********************************/
173 int
174 XYT_free(xyt_ADT xyt_handle)
175 {
176   comm_init();
177   check_handle(xyt_handle);
178   n_xyt_handles--;
179 
180   free(xyt_handle->info->nsep);
181   free(xyt_handle->info->lnsep);
182   free(xyt_handle->info->fo);
183   free(xyt_handle->info->stages);
184   free(xyt_handle->info->solve_uu);
185   free(xyt_handle->info->solve_w);
186   free(xyt_handle->info->x);
187   free(xyt_handle->info->xcol_vals);
188   free(xyt_handle->info->xcol_sz);
189   free(xyt_handle->info->xcol_indices);
190   free(xyt_handle->info->y);
191   free(xyt_handle->info->ycol_vals);
192   free(xyt_handle->info->ycol_sz);
193   free(xyt_handle->info->ycol_indices);
194   free(xyt_handle->info);
195   free(xyt_handle->mvi->local2global);
196    gs_free(xyt_handle->mvi->gs_handle);
197   free(xyt_handle->mvi);
198   free(xyt_handle);
199 
200 
201   /* if the check fails we nuke */
202   /* if NULL pointer passed to free we nuke */
203   /* if the calls to free fail that's not my problem */
204   return(0);
205 }
206 
207 
208 
209 /*************************************xyt.c************************************
210 Function:
211 
212 Input :
213 Output:
214 Return:
215 Description:
216 **************************************xyt.c***********************************/
217 int
218 XYT_stats(xyt_ADT xyt_handle)
219 {
220   int  op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
221   int fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
222   int   vals[9],  work[9];
223   PetscScalar fvals[3], fwork[3];
224 
225 
226   comm_init();
227   check_handle(xyt_handle);
228 
229   /* if factorization not done there are no stats */
230   if (!xyt_handle->info||!xyt_handle->mvi)
231     {
232       if (!my_id)
233 	{printf("XYT_stats() :: no stats available!\n");}
234       return 1;
235     }
236 
237   vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
238   vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
239   vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
240   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
241 
242   fvals[0]=fvals[1]=fvals[2]
243     =xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
244   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
245 
246   if (!my_id)
247     {
248       printf("%d :: min   xyt_nnz=%d\n",my_id,vals[0]);
249       printf("%d :: max   xyt_nnz=%d\n",my_id,vals[1]);
250       printf("%d :: avg   xyt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
251       printf("%d :: tot   xyt_nnz=%d\n",my_id,vals[2]);
252       printf("%d :: xyt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
253       printf("%d :: xyt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
254       printf("%d :: min   xyt_n  =%d\n",my_id,vals[3]);
255       printf("%d :: max   xyt_n  =%d\n",my_id,vals[4]);
256       printf("%d :: avg   xyt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
257       printf("%d :: tot   xyt_n  =%d\n",my_id,vals[5]);
258       printf("%d :: min   xyt_buf=%d\n",my_id,vals[6]);
259       printf("%d :: max   xyt_buf=%d\n",my_id,vals[7]);
260       printf("%d :: avg   xyt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
261       printf("%d :: min   xyt_slv=%g\n",my_id,fvals[0]);
262       printf("%d :: max   xyt_slv=%g\n",my_id,fvals[1]);
263       printf("%d :: avg   xyt_slv=%g\n",my_id,fvals[2]/num_nodes);
264     }
265 
266   return(0);
267 }
268 
269 
270 /*************************************xyt.c************************************
271 Function: do_xyt_factor
272 
273 Input :
274 Output:
275 Return:
276 Description: get A_local, local portion of global coarse matrix which
277 is a row dist. nxm matrix w/ n<m.
278    o my_ml holds address of ML struct associated w/A_local and coarse grid
279    o local2global holds global number of column i (i=0,...,m-1)
280    o local2global holds global number of row    i (i=0,...,n-1)
281    o mylocmatvec performs A_local . vec_local (note that gs is performed using
282    gs_init/gop).
283 
284 mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
285 mylocmatvec (void :: void *data, double *in, double *out)
286 **************************************xyt.c***********************************/
287 static
288 int
289 do_xyt_factor(xyt_ADT xyt_handle)
290 {
291   int flag;
292 
293 
294   flag=xyt_generate(xyt_handle);
295   return(flag);
296 }
297 
298 
299 /*************************************xyt.c************************************
300 Function:
301 
302 Input :
303 Output:
304 Return:
305 Description:
306 **************************************xyt.c***********************************/
307 static
308 int
309 xyt_generate(xyt_ADT xyt_handle)
310 {
311   int i,j,k,idx;
312   int dim, col;
313   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
314   int *segs;
315   int op[] = {GL_ADD,0};
316   int off, len;
317   PetscScalar *x_ptr, *y_ptr;
318   int *iptr, flag;
319   int start=0, end, work;
320   int op2[] = {GL_MIN,0};
321   gs_ADT gs_handle;
322   int *nsep, *lnsep, *fo;
323   int a_n=xyt_handle->mvi->n;
324   int a_m=xyt_handle->mvi->m;
325   int *a_local2global=xyt_handle->mvi->local2global;
326   int level;
327   int n, m;
328   int *xcol_sz, *xcol_indices, *stages;
329   PetscScalar **xcol_vals, *x;
330   int *ycol_sz, *ycol_indices;
331   PetscScalar **ycol_vals, *y;
332   int n_global;
333   int xt_nnz=0, xt_max_nnz=0;
334   int yt_nnz=0, yt_max_nnz=0;
335   int xt_zero_nnz  =0;
336   int xt_zero_nnz_0=0;
337   int yt_zero_nnz  =0;
338   int yt_zero_nnz_0=0;
339   PetscBLASInt i1 = 1;
340   PetscScalar dm1 = -1.0;
341 
342   n=xyt_handle->mvi->n;
343   nsep=xyt_handle->info->nsep;
344   lnsep=xyt_handle->info->lnsep;
345   fo=xyt_handle->info->fo;
346   end=lnsep[0];
347   level=xyt_handle->level;
348   gs_handle=xyt_handle->mvi->gs_handle;
349 
350   /* is there a null space? */
351   /* LATER add in ability to detect null space by checking alpha */
352   for (i=0, j=0; i<=level; i++)
353     {j+=nsep[i];}
354 
355   m = j-xyt_handle->ns;
356   if (m!=j)
357     {printf("xyt_generate() :: null space exists %d %d %d\n",m,j,xyt_handle->ns);}
358 
359   error_msg_warning("xyt_generate() :: X(%d,%d)\n",n,m);
360 
361   /* get and initialize storage for x local         */
362   /* note that x local is nxm and stored by columns */
363   xcol_sz = (int*) malloc(m*sizeof(PetscInt));
364   xcol_indices = (int*) malloc((2*m+1)*sizeof(int));
365   xcol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
366   for (i=j=0; i<m; i++, j+=2)
367     {
368       xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
369       xcol_vals[i] = NULL;
370     }
371   xcol_indices[j]=-1;
372 
373   /* get and initialize storage for y local         */
374   /* note that y local is nxm and stored by columns */
375   ycol_sz = (int*) malloc(m*sizeof(PetscInt));
376   ycol_indices = (int*) malloc((2*m+1)*sizeof(int));
377   ycol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
378   for (i=j=0; i<m; i++, j+=2)
379     {
380       ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
381       ycol_vals[i] = NULL;
382     }
383   ycol_indices[j]=-1;
384 
385   /* size of separators for each sub-hc working from bottom of tree to top */
386   /* this looks like nsep[]=segments */
387   stages = (int*) malloc((level+1)*sizeof(PetscInt));
388   segs   = (int*) malloc((level+1)*sizeof(PetscInt));
389   ivec_zero(stages,level+1);
390   ivec_copy(segs,nsep,level+1);
391   for (i=0; i<level; i++)
392     {segs[i+1] += segs[i];}
393   stages[0] = segs[0];
394 
395   /* temporary vectors  */
396   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
397   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
398   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
399   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
400   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
401 
402   /* extra nnz due to replication of vertices across separators */
403   for (i=1, j=0; i<=level; i++)
404     {j+=nsep[i];}
405 
406   /* storage for sparse x values */
407   n_global = xyt_handle->info->n_global;
408   xt_max_nnz = yt_max_nnz = (int)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
409   x = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
410   y = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
411 
412   /* LATER - can embed next sep to fire in gs */
413   /* time to make the donuts - generate X factor */
414   for (dim=i=j=0;i<m;i++)
415     {
416       /* time to move to the next level? */
417       while (i==segs[dim])
418 	{
419 	  if (dim==level)
420 	    {error_msg_fatal("dim about to exceed level\n"); break;}
421 
422 	  stages[dim++]=i;
423 	  end+=lnsep[dim];
424 	}
425       stages[dim]=i;
426 
427       /* which column are we firing? */
428       /* i.e. set v_l */
429       /* use new seps and do global min across hc to determine which one to fire */
430       (start<end) ? (col=fo[start]) : (col=INT_MAX);
431       giop_hc(&col,&work,1,op2,dim);
432 
433       /* shouldn't need this */
434       if (col==INT_MAX)
435 	{
436 	  error_msg_warning("hey ... col==INT_MAX??\n");
437 	  continue;
438 	}
439 
440       /* do I own it? I should */
441       rvec_zero(v ,a_m);
442       if (col==fo[start])
443 	{
444 	  start++;
445 	  idx=ivec_linear_search(col, a_local2global, a_n);
446 	  if (idx!=-1)
447 	    {v[idx] = 1.0; j++;}
448 	  else
449 	    {error_msg_fatal("NOT FOUND!\n");}
450 	}
451       else
452 	{
453 	  idx=ivec_linear_search(col, a_local2global, a_m);
454 	  if (idx!=-1)
455 	    {v[idx] = 1.0;}
456 	}
457 
458       /* perform u = A.v_l */
459       rvec_zero(u,n);
460       do_matvec(xyt_handle->mvi,v,u);
461 
462       /* uu =  X^T.u_l (local portion) */
463       /* technically only need to zero out first i entries */
464       /* later turn this into an XYT_solve call ? */
465       rvec_zero(uu,m);
466       y_ptr=y;
467       iptr = ycol_indices;
468       for (k=0; k<i; k++)
469 	{
470 	  off = *iptr++;
471 	  len = *iptr++;
472 
473 	  uu[k] = BLASdot_(&len,u+off,&i1,y_ptr,&i1);
474 	  y_ptr+=len;
475 	}
476 
477       /* uu = X^T.u_l (comm portion) */
478       ssgl_radd  (uu, w, dim, stages);
479 
480       /* z = X.uu */
481       rvec_zero(z,n);
482       x_ptr=x;
483       iptr = xcol_indices;
484       for (k=0; k<i; k++)
485 	{
486 	  off = *iptr++;
487 	  len = *iptr++;
488 
489 	  BLASaxpy_(&len,&uu[k],x_ptr,&i1,z+off,&i1);
490 	  x_ptr+=len;
491 	}
492 
493       /* compute v_l = v_l - z */
494       rvec_zero(v+a_n,a_m-a_n);
495       BLASaxpy_(&n,&dm1,z,&i1,v,&i1);
496 
497       /* compute u_l = A.v_l */
498       if (a_n!=a_m)
499 	{gs_gop_hc(gs_handle,v,"+\0",dim);}
500       rvec_zero(u,n);
501      do_matvec(xyt_handle->mvi,v,u);
502 
503       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
504       alpha = BLASdot_(&n,u,&i1,u,&i1);
505       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
506       grop_hc(&alpha, &alpha_w, 1, op, dim);
507 
508       alpha = (PetscScalar) sqrt((double)alpha);
509 
510       /* check for small alpha                             */
511       /* LATER use this to detect and determine null space */
512       if (fabs(alpha)<1.0e-14)
513 	{error_msg_fatal("bad alpha! %g\n",alpha);}
514 
515       /* compute v_l = v_l/sqrt(alpha) */
516       rvec_scale(v,1.0/alpha,n);
517       rvec_scale(u,1.0/alpha,n);
518 
519       /* add newly generated column, v_l, to X */
520       flag = 1;
521       off=len=0;
522       for (k=0; k<n; k++)
523 	{
524 	  if (v[k]!=0.0)
525 	    {
526 	      len=k;
527 	      if (flag)
528 		{off=k; flag=0;}
529 	    }
530 	}
531 
532       len -= (off-1);
533 
534       if (len>0)
535 	{
536 	  if ((xt_nnz+len)>xt_max_nnz)
537 	    {
538 	      error_msg_warning("increasing space for X by 2x!\n");
539 	      xt_max_nnz *= 2;
540 	      x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
541 	      rvec_copy(x_ptr,x,xt_nnz);
542 	      free(x);
543 	      x = x_ptr;
544 	      x_ptr+=xt_nnz;
545 	    }
546 	  xt_nnz += len;
547 	  rvec_copy(x_ptr,v+off,len);
548 
549           /* keep track of number of zeros */
550 	  if (dim)
551 	    {
552 	      for (k=0; k<len; k++)
553 		{
554 		  if (x_ptr[k]==0.0)
555 		    {xt_zero_nnz++;}
556 		}
557 	    }
558 	  else
559 	    {
560 	      for (k=0; k<len; k++)
561 		{
562 		  if (x_ptr[k]==0.0)
563 		    {xt_zero_nnz_0++;}
564 		}
565 	    }
566 	  xcol_indices[2*i] = off;
567 	  xcol_sz[i] = xcol_indices[2*i+1] = len;
568 	  xcol_vals[i] = x_ptr;
569 	}
570       else
571 	{
572 	  xcol_indices[2*i] = 0;
573 	  xcol_sz[i] = xcol_indices[2*i+1] = 0;
574 	  xcol_vals[i] = x_ptr;
575 	}
576 
577 
578       /* add newly generated column, u_l, to Y */
579       flag = 1;
580       off=len=0;
581       for (k=0; k<n; k++)
582 	{
583 	  if (u[k]!=0.0)
584 	    {
585 	      len=k;
586 	      if (flag)
587 		{off=k; flag=0;}
588 	    }
589 	}
590 
591       len -= (off-1);
592 
593       if (len>0)
594 	{
595 	  if ((yt_nnz+len)>yt_max_nnz)
596 	    {
597 	      error_msg_warning("increasing space for Y by 2x!\n");
598 	      yt_max_nnz *= 2;
599 	      y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
600 	      rvec_copy(y_ptr,y,yt_nnz);
601 	      free(y);
602 	      y = y_ptr;
603 	      y_ptr+=yt_nnz;
604 	    }
605 	  yt_nnz += len;
606 	  rvec_copy(y_ptr,u+off,len);
607 
608           /* keep track of number of zeros */
609 	  if (dim)
610 	    {
611 	      for (k=0; k<len; k++)
612 		{
613 		  if (y_ptr[k]==0.0)
614 		    {yt_zero_nnz++;}
615 		}
616 	    }
617 	  else
618 	    {
619 	      for (k=0; k<len; k++)
620 		{
621 		  if (y_ptr[k]==0.0)
622 		    {yt_zero_nnz_0++;}
623 		}
624 	    }
625 	  ycol_indices[2*i] = off;
626 	  ycol_sz[i] = ycol_indices[2*i+1] = len;
627 	  ycol_vals[i] = y_ptr;
628 	}
629       else
630 	{
631 	  ycol_indices[2*i] = 0;
632 	  ycol_sz[i] = ycol_indices[2*i+1] = 0;
633 	  ycol_vals[i] = y_ptr;
634 	}
635     }
636 
637   /* close off stages for execution phase */
638   while (dim!=level)
639     {
640       stages[dim++]=i;
641       error_msg_warning("disconnected!!! dim(%d)!=level(%d)\n",dim,level);
642     }
643   stages[dim]=i;
644 
645   xyt_handle->info->n=xyt_handle->mvi->n;
646   xyt_handle->info->m=m;
647   xyt_handle->info->nnz=xt_nnz + yt_nnz;
648   xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz;
649   xyt_handle->info->msg_buf_sz=stages[level]-stages[0];
650   xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
651   xyt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
652   xyt_handle->info->x=x;
653   xyt_handle->info->xcol_vals=xcol_vals;
654   xyt_handle->info->xcol_sz=xcol_sz;
655   xyt_handle->info->xcol_indices=xcol_indices;
656   xyt_handle->info->stages=stages;
657   xyt_handle->info->y=y;
658   xyt_handle->info->ycol_vals=ycol_vals;
659   xyt_handle->info->ycol_sz=ycol_sz;
660   xyt_handle->info->ycol_indices=ycol_indices;
661 
662   free(segs);
663   free(u);
664   free(v);
665   free(uu);
666   free(z);
667   free(w);
668 
669   return(0);
670 }
671 
672 
673 /*************************************xyt.c************************************
674 Function:
675 
676 Input :
677 Output:
678 Return:
679 Description:
680 **************************************xyt.c***********************************/
681 static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle,  PetscScalar *uc)
682 {
683   int off, len, *iptr;
684   int level       =xyt_handle->level;
685   int n           =xyt_handle->info->n;
686   int m           =xyt_handle->info->m;
687   int *stages     =xyt_handle->info->stages;
688   int *xcol_indices=xyt_handle->info->xcol_indices;
689   int *ycol_indices=xyt_handle->info->ycol_indices;
690    PetscScalar *x_ptr, *y_ptr, *uu_ptr;
691   PetscScalar *solve_uu=xyt_handle->info->solve_uu;
692   PetscScalar *solve_w =xyt_handle->info->solve_w;
693   PetscScalar *x       =xyt_handle->info->x;
694   PetscScalar *y       =xyt_handle->info->y;
695   PetscBLASInt i1 = 1;
696 
697   PetscFunctionBegin;
698   uu_ptr=solve_uu;
699   rvec_zero(uu_ptr,m);
700 
701   /* x  = X.Y^T.b */
702   /* uu = Y^T.b */
703   for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len)
704     {
705       off=*iptr++; len=*iptr++;
706       *uu_ptr++ = BLASdot_(&len,uc+off,&i1,y_ptr,&i1);
707     }
708 
709   /* comunication of beta */
710   uu_ptr=solve_uu;
711   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}
712 
713   rvec_zero(uc,n);
714 
715   /* x = X.uu */
716   for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len)
717     {
718       off=*iptr++; len=*iptr++;
719       BLASaxpy_(&len,uu_ptr++,x_ptr,&i1,uc+off,&i1);
720     }
721   PetscFunctionReturn(0);
722 }
723 
724 
725 /*************************************xyt.c************************************
726 Function: check_handle()
727 
728 Input :
729 Output:
730 Return:
731 Description:
732 **************************************xyt.c***********************************/
733 static PetscErrorCode check_handle(xyt_ADT xyt_handle)
734 {
735   int vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
736 
737   PetscFunctionBegin;
738   if (xyt_handle==NULL)
739     {error_msg_fatal("check_handle() :: bad handle :: NULL %d\n",xyt_handle);}
740 
741   vals[0]=vals[1]=xyt_handle->id;
742   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
743   if ((vals[0]!=vals[1])||(xyt_handle->id<=0))
744     {error_msg_fatal("check_handle() :: bad handle :: id mismatch min/max %d/%d %d\n",
745 		     vals[0],vals[1], xyt_handle->id);}
746   PetscFunctionReturn(0);
747 }
748 
749 
750 /*************************************xyt.c************************************
751 Function: det_separators
752 
753 Input :
754 Output:
755 Return:
756 Description:
757   det_separators(xyt_handle, local2global, n, m, mylocmatvec, grid_data);
758 **************************************xyt.c***********************************/
759 static PetscErrorCode det_separators(xyt_ADT xyt_handle)
760 {
761   int i, ct, id;
762   int mask, edge, *iptr;
763   int *dir, *used;
764   int sum[4], w[4];
765   PetscScalar rsum[4], rw[4];
766   int op[] = {GL_ADD,0};
767   PetscScalar *lhs, *rhs;
768   int *nsep, *lnsep, *fo, nfo=0;
769   gs_ADT gs_handle=xyt_handle->mvi->gs_handle;
770   int *local2global=xyt_handle->mvi->local2global;
771   int  n=xyt_handle->mvi->n;
772   int  m=xyt_handle->mvi->m;
773   int level=xyt_handle->level;
774   int shared=FALSE;
775 
776   PetscFunctionBegin;
777   dir  = (int*)malloc(sizeof(PetscInt)*(level+1));
778   nsep = (int*)malloc(sizeof(PetscInt)*(level+1));
779   lnsep= (int*)malloc(sizeof(PetscInt)*(level+1));
780   fo   = (int*)malloc(sizeof(PetscInt)*(n+1));
781   used = (int*)malloc(sizeof(PetscInt)*n);
782 
783   ivec_zero(dir  ,level+1);
784   ivec_zero(nsep ,level+1);
785   ivec_zero(lnsep,level+1);
786   ivec_set (fo   ,-1,n+1);
787   ivec_zero(used,n);
788 
789   lhs  = (double*)malloc(sizeof(PetscScalar)*m);
790   rhs  = (double*)malloc(sizeof(PetscScalar)*m);
791 
792   /* determine the # of unique dof */
793   rvec_zero(lhs,m);
794   rvec_set(lhs,1.0,n);
795   gs_gop_hc(gs_handle,lhs,"+\0",level);
796   error_msg_warning("done first gs_gop_hc\n");
797   rvec_zero(rsum,2);
798   for (ct=i=0;i<n;i++)
799     {
800       if (lhs[i]!=0.0)
801 	{rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
802 
803       if (lhs[i]!=1.0)
804 	{
805           shared=TRUE;
806         }
807     }
808 
809   grop_hc(rsum,rw,2,op,level);
810   rsum[0]+=0.1;
811   rsum[1]+=0.1;
812 
813   /*
814       if (!my_id)
815       {
816       printf("xyt n unique = %d (%g)\n",(int) rsum[0], rsum[0]);
817       printf("xyt n shared = %d (%g)\n",(int) rsum[1], rsum[1]);
818       }
819   */
820 
821   xyt_handle->info->n_global=xyt_handle->info->m_global=(int) rsum[0];
822   xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(int) rsum[0];
823 
824   /* determine separator sets top down */
825   if (shared)
826     {
827       /* solution is to do as in the symmetric shared case but then */
828       /* pick the sub-hc with the most free dofs and do a mat-vec   */
829       /* and pick up the responses on the other sub-hc from the     */
830       /* initial separator set obtained from the symm. shared case  */
831       error_msg_fatal("shared dof separator determination not ready ... see hmt!!!\n");
832       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
833 	{
834 	  /* set rsh of hc, fire, and collect lhs responses */
835 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
836 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
837 
838 	  /* set lsh of hc, fire, and collect rhs responses */
839 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
840 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
841 
842 	  for (i=0;i<n;i++)
843 	    {
844 	      if (id< mask)
845 		{
846 		  if (lhs[i]!=0.0)
847 		    {lhs[i]=1.0;}
848 		}
849 	      if (id>=mask)
850 		{
851 		  if (rhs[i]!=0.0)
852 		    {rhs[i]=1.0;}
853 		}
854 	    }
855 
856 	  if (id< mask)
857 	    {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
858 	  else
859 	    {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}
860 
861 	  /* count number of dofs I own that have signal and not in sep set */
862 	  rvec_zero(rsum,4);
863 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
864 	    {
865 	      if (!used[i])
866 		{
867 		  /* number of unmarked dofs on node */
868 		  ct++;
869 		  /* number of dofs to be marked on lhs hc */
870 		  if (id< mask)
871 		    {
872 		      if (lhs[i]!=0.0)
873 			{sum[0]++; rsum[0]+=1.0/lhs[i];}
874 		    }
875 		  /* number of dofs to be marked on rhs hc */
876 		  if (id>=mask)
877 		    {
878 		      if (rhs[i]!=0.0)
879 			{sum[1]++; rsum[1]+=1.0/rhs[i];}
880 		    }
881 		}
882 	    }
883 
884 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
885 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
886 	  (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
887 	  giop_hc(sum,w,4,op,edge);
888 	  grop_hc(rsum,rw,4,op,edge);
889 	  rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
890 
891 	  if (id<mask)
892 	    {
893 	      /* mark dofs I own that have signal and not in sep set */
894 	      for (ct=i=0;i<n;i++)
895 		{
896 		  if ((!used[i])&&(lhs[i]!=0.0))
897 		    {
898 		      ct++; nfo++;
899 
900 		      if (nfo>n)
901 			{error_msg_fatal("nfo about to exceed n\n");}
902 
903 		      *--iptr = local2global[i];
904 		      used[i]=edge;
905 		    }
906 		}
907 	      if (ct>1) {ivec_sort(iptr,ct);}
908 
909 	      lnsep[edge]=ct;
910 	      nsep[edge]=(int) rsum[0];
911 	      dir [edge]=LEFT;
912 	    }
913 
914 	  if (id>=mask)
915 	    {
916 	      /* mark dofs I own that have signal and not in sep set */
917 	      for (ct=i=0;i<n;i++)
918 		{
919 		  if ((!used[i])&&(rhs[i]!=0.0))
920 		    {
921 		      ct++; nfo++;
922 
923 		      if (nfo>n)
924 			{error_msg_fatal("nfo about to exceed n\n");}
925 
926 		      *--iptr = local2global[i];
927 		      used[i]=edge;
928 		    }
929 		}
930 	      if (ct>1) {ivec_sort(iptr,ct);}
931 
932 	      lnsep[edge]=ct;
933 	      nsep[edge]= (int) rsum[1];
934 	      dir [edge]=RIGHT;
935 	    }
936 
937 	  /* LATER or we can recur on these to order seps at this level */
938 	  /* do we need full set of separators for this?                */
939 
940 	  /* fold rhs hc into lower */
941 	  if (id>=mask)
942 	    {id-=mask;}
943 	}
944     }
945   else
946     {
947       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
948 	{
949 	  /* set rsh of hc, fire, and collect lhs responses */
950 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
951 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
952 
953 	  /* set lsh of hc, fire, and collect rhs responses */
954 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
955 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
956 
957 	  /* count number of dofs I own that have signal and not in sep set */
958 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
959 	    {
960 	      if (!used[i])
961 		{
962 		  /* number of unmarked dofs on node */
963 		  ct++;
964 		  /* number of dofs to be marked on lhs hc */
965 		  if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
966 		  /* number of dofs to be marked on rhs hc */
967 		  if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
968 		}
969 	    }
970 
971 	  /* for the non-symmetric case we need separators of width 2 */
972 	  /* so take both sides */
973 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
974 	  giop_hc(sum,w,4,op,edge);
975 
976 	  ct=0;
977 	  if (id<mask)
978 	    {
979 	      /* mark dofs I own that have signal and not in sep set */
980 	      for (i=0;i<n;i++)
981 		{
982 		  if ((!used[i])&&(lhs[i]!=0.0))
983 		    {
984 		      ct++; nfo++;
985 		      *--iptr = local2global[i];
986 		      used[i]=edge;
987 		    }
988 		}
989 	      /* LSH hc summation of ct should be sum[0] */
990 	    }
991 	  else
992 	    {
993 	      /* mark dofs I own that have signal and not in sep set */
994 	      for (i=0;i<n;i++)
995 		{
996 		  if ((!used[i])&&(rhs[i]!=0.0))
997 		    {
998 		      ct++; nfo++;
999 		      *--iptr = local2global[i];
1000 		      used[i]=edge;
1001 		    }
1002 		}
1003 	      /* RSH hc summation of ct should be sum[1] */
1004 	    }
1005 
1006 	  if (ct>1) {ivec_sort(iptr,ct);}
1007 	  lnsep[edge]=ct;
1008 	  nsep[edge]=sum[0]+sum[1];
1009 	  dir [edge]=BOTH;
1010 
1011 	  /* LATER or we can recur on these to order seps at this level */
1012 	  /* do we need full set of separators for this?                */
1013 
1014 	  /* fold rhs hc into lower */
1015 	  if (id>=mask)
1016 	    {id-=mask;}
1017 	}
1018     }
1019 
1020   /* level 0 is on processor case - so mark the remainder */
1021   for (ct=i=0;i<n;i++)
1022     {
1023       if (!used[i])
1024 	{
1025 	  ct++; nfo++;
1026 	  *--iptr = local2global[i];
1027 	  used[i]=edge;
1028 	}
1029     }
1030   if (ct>1) {ivec_sort(iptr,ct);}
1031   lnsep[edge]=ct;
1032   nsep [edge]=ct;
1033   dir  [edge]=BOTH;
1034 
1035   xyt_handle->info->nsep=nsep;
1036   xyt_handle->info->lnsep=lnsep;
1037   xyt_handle->info->fo=fo;
1038   xyt_handle->info->nfo=nfo;
1039 
1040   free(dir);
1041   free(lhs);
1042   free(rhs);
1043   free(used);
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 
1048 /*************************************xyt.c************************************
1049 Function: set_mvi
1050 
1051 Input :
1052 Output:
1053 Return:
1054 Description:
1055 **************************************xyt.c***********************************/
1056 static
1057 mv_info *set_mvi(int *local2global, int n, int m, void *matvec, void *grid_data)
1058 {
1059   mv_info *mvi;
1060 
1061 
1062   mvi = (mv_info*)malloc(sizeof(mv_info));
1063   mvi->n=n;
1064   mvi->m=m;
1065   mvi->n_global=-1;
1066   mvi->m_global=-1;
1067   mvi->local2global=(int*)malloc((m+1)*sizeof(PetscInt));
1068   ivec_copy(mvi->local2global,local2global,m);
1069   mvi->local2global[m] = INT_MAX;
1070   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
1071   mvi->grid_data=grid_data;
1072 
1073   /* set xyt communication handle to perform restricted matvec */
1074   mvi->gs_handle = gs_init(local2global, m, num_nodes);
1075 
1076   return(mvi);
1077 }
1078 
1079 
1080 /*************************************xyt.c************************************
1081 Function: set_mvi
1082 
1083 Input :
1084 Output:
1085 Return:
1086 Description:
1087 
1088       computes u = A.v
1089       do_matvec(xyt_handle->mvi,v,u);
1090 **************************************xyt.c***********************************/
1091 static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
1092 {
1093   PetscFunctionBegin;
1094   A->matvec((mv_info*)A->grid_data,v,u);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 
1099 
1100