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