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