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