xref: /petsc/src/ksp/pc/impls/tfs/xyt.c (revision 5c8f6a953e7ed1c81f507d64aebddb11080b60e9)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /*************************************xyt.c************************************
3827bd09bSSatish Balay Module Name: xyt
4827bd09bSSatish Balay Module Info:
5827bd09bSSatish Balay 
6827bd09bSSatish Balay author:  Henry M. Tufo III
7827bd09bSSatish Balay e-mail:  hmt@asci.uchicago.edu
8827bd09bSSatish Balay contact:
9827bd09bSSatish Balay +--------------------------------+--------------------------------+
10827bd09bSSatish Balay |MCS Division - Building 221     |Department of Computer Science  |
11827bd09bSSatish Balay |Argonne National Laboratory     |Ryerson 152                     |
12827bd09bSSatish Balay |9700 S. Cass Avenue             |The University of Chicago       |
13827bd09bSSatish Balay |Argonne, IL  60439              |Chicago, IL  60637              |
14827bd09bSSatish Balay |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
15827bd09bSSatish Balay +--------------------------------+--------------------------------+
16827bd09bSSatish Balay 
17827bd09bSSatish Balay Last Modification: 3.20.01
18827bd09bSSatish Balay **************************************xyt.c***********************************/
19c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
20827bd09bSSatish Balay 
21827bd09bSSatish Balay #define LEFT  -1
22827bd09bSSatish Balay #define RIGHT  1
23827bd09bSSatish Balay #define BOTH   0
24827bd09bSSatish Balay 
25827bd09bSSatish Balay typedef struct xyt_solver_info {
2652f87cdaSBarry Smith   PetscInt n, m, n_global, m_global;
2752f87cdaSBarry Smith   PetscInt nnz, max_nnz, msg_buf_sz;
2852f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo, *stages;
2952f87cdaSBarry Smith   PetscInt *xcol_sz, *xcol_indices;
30a501084fSBarry Smith   PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
3152f87cdaSBarry Smith   PetscInt *ycol_sz, *ycol_indices;
32a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
3352f87cdaSBarry Smith   PetscInt nsolves;
34a501084fSBarry Smith   PetscScalar tot_solve_time;
35827bd09bSSatish Balay } xyt_info;
36827bd09bSSatish Balay 
37827bd09bSSatish Balay 
38827bd09bSSatish Balay typedef struct matvec_info {
3952f87cdaSBarry Smith   PetscInt n, m, n_global, m_global;
4052f87cdaSBarry Smith   PetscInt *local2global;
41ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle;
42a501084fSBarry Smith   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
43827bd09bSSatish Balay   void *grid_data;
44827bd09bSSatish Balay } mv_info;
45827bd09bSSatish Balay 
46827bd09bSSatish Balay struct xyt_CDT{
4752f87cdaSBarry Smith   PetscInt id;
4852f87cdaSBarry Smith   PetscInt ns;
4952f87cdaSBarry Smith   PetscInt level;
50827bd09bSSatish Balay   xyt_info *info;
51827bd09bSSatish Balay   mv_info  *mvi;
52827bd09bSSatish Balay };
53827bd09bSSatish Balay 
5452f87cdaSBarry Smith static PetscInt n_xyt=0;
5552f87cdaSBarry Smith static PetscInt n_xyt_handles=0;
56827bd09bSSatish Balay 
57827bd09bSSatish Balay /* prototypes */
583fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
593fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle);
603fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle);
613fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
6252f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle);
6352f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle);
64*5c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data);
65827bd09bSSatish Balay 
667b1ae94cSBarry Smith /**************************************xyt.c***********************************/
677b1ae94cSBarry Smith xyt_ADT XYT_new(void)
68827bd09bSSatish Balay {
69827bd09bSSatish Balay   xyt_ADT xyt_handle;
70827bd09bSSatish Balay 
71827bd09bSSatish Balay   /* rolling count on n_xyt ... pot. problem here */
72827bd09bSSatish Balay   n_xyt_handles++;
73a501084fSBarry Smith   xyt_handle       = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
74827bd09bSSatish Balay   xyt_handle->id   = ++n_xyt;
75827bd09bSSatish Balay   xyt_handle->info = NULL;
76827bd09bSSatish Balay   xyt_handle->mvi  = NULL;
77827bd09bSSatish Balay 
78827bd09bSSatish Balay   return(xyt_handle);
79827bd09bSSatish Balay }
80827bd09bSSatish Balay 
817b1ae94cSBarry Smith /**************************************xyt.c***********************************/
8252f87cdaSBarry Smith PetscInt XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt  handle */
8352f87cdaSBarry Smith            PetscInt *local2global,      /* global column mapping       */
8452f87cdaSBarry Smith            PetscInt n,                  /* local num rows              */
8552f87cdaSBarry Smith            PetscInt m,                  /* local num cols              */
86*5c8f6a95SKarl Rupp            PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*),                /* b_loc=A_local.x_loc         */
871147fc2aSKarl Rupp            void *grid_data)             /* grid data for matvec        */
88827bd09bSSatish Balay {
89827bd09bSSatish Balay 
90b1c944f5SJed Brown   PCTFS_comm_init();
91827bd09bSSatish Balay   check_handle(xyt_handle);
92827bd09bSSatish Balay 
93827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
94b1c944f5SJed Brown   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);
95827bd09bSSatish Balay 
96827bd09bSSatish Balay   /* space for X info */
97a501084fSBarry Smith   xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));
98827bd09bSSatish Balay 
99827bd09bSSatish Balay   /* set up matvec handles */
100*5c8f6a95SKarl Rupp   xyt_handle->mvi  = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data);
101827bd09bSSatish Balay 
102827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
103827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
104827bd09bSSatish Balay   xyt_handle->ns=0;
105827bd09bSSatish Balay 
106827bd09bSSatish Balay   /* determine separators and generate firing order - NB xyt info set here */
107827bd09bSSatish Balay   det_separators(xyt_handle);
108827bd09bSSatish Balay 
109827bd09bSSatish Balay   return(do_xyt_factor(xyt_handle));
110827bd09bSSatish Balay }
111827bd09bSSatish Balay 
1127b1ae94cSBarry Smith /**************************************xyt.c***********************************/
1138cda6cd7SBarry Smith PetscInt XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b)
114827bd09bSSatish Balay {
115b1c944f5SJed Brown   PCTFS_comm_init();
116827bd09bSSatish Balay   check_handle(xyt_handle);
117827bd09bSSatish Balay 
118827bd09bSSatish Balay   /* need to copy b into x? */
119827bd09bSSatish Balay   if (b)
120ca8e9878SJed Brown     {PCTFS_rvec_copy(x,b,xyt_handle->mvi->n);}
121827bd09bSSatish Balay   do_xyt_solve(xyt_handle,x);
122827bd09bSSatish Balay 
123827bd09bSSatish Balay   return(0);
124827bd09bSSatish Balay }
125827bd09bSSatish Balay 
1267b1ae94cSBarry Smith /**************************************xyt.c***********************************/
12752f87cdaSBarry Smith PetscInt XYT_free(xyt_ADT xyt_handle)
128827bd09bSSatish Balay {
129b1c944f5SJed Brown   PCTFS_comm_init();
130827bd09bSSatish Balay   check_handle(xyt_handle);
131827bd09bSSatish Balay   n_xyt_handles--;
132827bd09bSSatish Balay 
133a501084fSBarry Smith   free(xyt_handle->info->nsep);
134a501084fSBarry Smith   free(xyt_handle->info->lnsep);
135a501084fSBarry Smith   free(xyt_handle->info->fo);
136a501084fSBarry Smith   free(xyt_handle->info->stages);
137a501084fSBarry Smith   free(xyt_handle->info->solve_uu);
138a501084fSBarry Smith   free(xyt_handle->info->solve_w);
139a501084fSBarry Smith   free(xyt_handle->info->x);
140a501084fSBarry Smith   free(xyt_handle->info->xcol_vals);
141a501084fSBarry Smith   free(xyt_handle->info->xcol_sz);
142a501084fSBarry Smith   free(xyt_handle->info->xcol_indices);
143a501084fSBarry Smith   free(xyt_handle->info->y);
144a501084fSBarry Smith   free(xyt_handle->info->ycol_vals);
145a501084fSBarry Smith   free(xyt_handle->info->ycol_sz);
146a501084fSBarry Smith   free(xyt_handle->info->ycol_indices);
147a501084fSBarry Smith   free(xyt_handle->info);
148a501084fSBarry Smith   free(xyt_handle->mvi->local2global);
149ca8e9878SJed Brown   PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle);
150a501084fSBarry Smith   free(xyt_handle->mvi);
151a501084fSBarry Smith   free(xyt_handle);
152827bd09bSSatish Balay 
153827bd09bSSatish Balay 
154827bd09bSSatish Balay   /* if the check fails we nuke */
155a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
156827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
157827bd09bSSatish Balay   return(0);
158827bd09bSSatish Balay }
159827bd09bSSatish Balay 
1607b1ae94cSBarry Smith /**************************************xyt.c***********************************/
16152f87cdaSBarry Smith PetscInt XYT_stats(xyt_ADT xyt_handle)
162827bd09bSSatish Balay {
16352f87cdaSBarry Smith   PetscInt    op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
16452f87cdaSBarry Smith   PetscInt   fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
16552f87cdaSBarry Smith   PetscInt    vals[9],  work[9];
166a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
167827bd09bSSatish Balay 
168b1c944f5SJed Brown   PCTFS_comm_init();
169827bd09bSSatish Balay   check_handle(xyt_handle);
170827bd09bSSatish Balay 
171827bd09bSSatish Balay   /* if factorization not done there are no stats */
1727d0a6c19SBarry Smith   if (!xyt_handle->info||!xyt_handle->mvi) {
173b1c944f5SJed Brown     if (!PCTFS_my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n");
174827bd09bSSatish Balay     return 1;
175827bd09bSSatish Balay   }
176827bd09bSSatish Balay 
177827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
178827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
179827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
180b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
181827bd09bSSatish Balay 
182827bd09bSSatish Balay   fvals[0]=fvals[1]=fvals[2]
183827bd09bSSatish Balay     =xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
184b1c944f5SJed Brown   PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
185827bd09bSSatish Balay 
186b1c944f5SJed Brown   if (!PCTFS_my_id) {
187b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_nnz=%D\n",PCTFS_my_id,vals[0]);
188b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_nnz=%D\n",PCTFS_my_id,vals[1]);
189b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);
190b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_nnz=%D\n",PCTFS_my_id,vals[2]);
191b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(2d)  =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.5)));
192b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(3d)  =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
193b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_n  =%D\n",PCTFS_my_id,vals[3]);
194b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_n  =%D\n",PCTFS_my_id,vals[4]);
195b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_n  =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);
196b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_n  =%D\n",PCTFS_my_id,vals[5]);
197b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_buf=%D\n",PCTFS_my_id,vals[6]);
198b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_buf=%D\n",PCTFS_my_id,vals[7]);
199b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);
200b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_slv=%g\n",PCTFS_my_id,fvals[0]);
201b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_slv=%g\n",PCTFS_my_id,fvals[1]);
202b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);
203827bd09bSSatish Balay   }
204827bd09bSSatish Balay 
205827bd09bSSatish Balay   return(0);
206827bd09bSSatish Balay }
207827bd09bSSatish Balay 
208827bd09bSSatish Balay 
209827bd09bSSatish Balay /*************************************xyt.c************************************
210827bd09bSSatish Balay 
211827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
212827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
213827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
214827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
215827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
216827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
217ca8e9878SJed Brown    PCTFS_gs_init/gop).
218827bd09bSSatish Balay 
219827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
220827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
221827bd09bSSatish Balay **************************************xyt.c***********************************/
22252f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle)
223827bd09bSSatish Balay {
2247b1ae94cSBarry Smith   return xyt_generate(xyt_handle);
225827bd09bSSatish Balay }
226827bd09bSSatish Balay 
2277b1ae94cSBarry Smith /**************************************xyt.c***********************************/
22852f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle)
229827bd09bSSatish Balay {
23052f87cdaSBarry Smith   PetscInt i,j,k,idx;
23152f87cdaSBarry Smith   PetscInt dim, col;
232a501084fSBarry Smith   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
23352f87cdaSBarry Smith   PetscInt *segs;
23452f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
23552f87cdaSBarry Smith   PetscInt off, len;
236a501084fSBarry Smith   PetscScalar *x_ptr, *y_ptr;
23752f87cdaSBarry Smith   PetscInt *iptr, flag;
23852f87cdaSBarry Smith   PetscInt start=0, end, work;
23952f87cdaSBarry Smith   PetscInt op2[] = {GL_MIN,0};
240ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle;
24152f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo;
24252f87cdaSBarry Smith   PetscInt a_n=xyt_handle->mvi->n;
24352f87cdaSBarry Smith   PetscInt a_m=xyt_handle->mvi->m;
24452f87cdaSBarry Smith   PetscInt *a_local2global=xyt_handle->mvi->local2global;
24552f87cdaSBarry Smith   PetscInt level;
24652f87cdaSBarry Smith   PetscInt n, m;
24752f87cdaSBarry Smith   PetscInt *xcol_sz, *xcol_indices, *stages;
248a501084fSBarry Smith   PetscScalar **xcol_vals, *x;
24952f87cdaSBarry Smith   PetscInt *ycol_sz, *ycol_indices;
250a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
25152f87cdaSBarry Smith   PetscInt n_global;
25252f87cdaSBarry Smith   PetscInt xt_nnz=0, xt_max_nnz=0;
25352f87cdaSBarry Smith   PetscInt yt_nnz=0, yt_max_nnz=0;
25452f87cdaSBarry Smith   PetscInt xt_zero_nnz  =0;
25552f87cdaSBarry Smith   PetscInt xt_zero_nnz_0=0;
25652f87cdaSBarry Smith   PetscInt yt_zero_nnz  =0;
25752f87cdaSBarry Smith   PetscInt yt_zero_nnz_0=0;
2584a0de3f6SBarry Smith   PetscBLASInt i1 = 1,dlen;
259a501084fSBarry Smith   PetscScalar dm1 = -1.0;
260d1528f56SBarry Smith   PetscErrorCode ierr;
261827bd09bSSatish Balay 
262827bd09bSSatish Balay   n=xyt_handle->mvi->n;
263827bd09bSSatish Balay   nsep=xyt_handle->info->nsep;
264827bd09bSSatish Balay   lnsep=xyt_handle->info->lnsep;
265827bd09bSSatish Balay   fo=xyt_handle->info->fo;
266827bd09bSSatish Balay   end=lnsep[0];
267827bd09bSSatish Balay   level=xyt_handle->level;
268ca8e9878SJed Brown   PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle;
269827bd09bSSatish Balay 
270827bd09bSSatish Balay   /* is there a null space? */
271827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
272db4deed7SKarl Rupp   for (i=0, j=0; i<=level; i++) { j+=nsep[i]; }
273827bd09bSSatish Balay 
274827bd09bSSatish Balay   m = j-xyt_handle->ns;
27522d28d08SBarry Smith   if (m!=j) {
27622d28d08SBarry Smith     ierr = PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);CHKERRQ(ierr);
27722d28d08SBarry Smith   }
278827bd09bSSatish Balay 
2794a0de3f6SBarry Smith   ierr = PetscInfo2(0,"xyt_generate() :: X(%D,%D)\n",n,m);CHKERRQ(ierr);
280827bd09bSSatish Balay 
281827bd09bSSatish Balay   /* get and initialize storage for x local         */
282827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
28352f87cdaSBarry Smith   xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
28452f87cdaSBarry Smith   xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
285a501084fSBarry Smith   xcol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
286db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
287827bd09bSSatish Balay     xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
288827bd09bSSatish Balay     xcol_vals[i] = NULL;
289827bd09bSSatish Balay   }
290827bd09bSSatish Balay   xcol_indices[j]=-1;
291827bd09bSSatish Balay 
292827bd09bSSatish Balay   /* get and initialize storage for y local         */
293827bd09bSSatish Balay   /* note that y local is nxm and stored by columns */
29452f87cdaSBarry Smith   ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
29552f87cdaSBarry Smith   ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
296a501084fSBarry Smith   ycol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
297db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
298827bd09bSSatish Balay     ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
299827bd09bSSatish Balay     ycol_vals[i] = NULL;
300827bd09bSSatish Balay   }
301827bd09bSSatish Balay   ycol_indices[j]=-1;
302827bd09bSSatish Balay 
303827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
304827bd09bSSatish Balay   /* this looks like nsep[]=segments */
30552f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
30652f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
307ca8e9878SJed Brown   PCTFS_ivec_zero(stages,level+1);
308ca8e9878SJed Brown   PCTFS_ivec_copy(segs,nsep,level+1);
309db4deed7SKarl Rupp   for (i=0; i<level; i++) { segs[i+1] += segs[i]; }
310827bd09bSSatish Balay   stages[0] = segs[0];
311827bd09bSSatish Balay 
312827bd09bSSatish Balay   /* temporary vectors  */
313a501084fSBarry Smith   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
314a501084fSBarry Smith   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
315a501084fSBarry Smith   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
316a501084fSBarry Smith   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
317a501084fSBarry Smith   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
318827bd09bSSatish Balay 
319827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
320db4deed7SKarl Rupp   for (i=1, j=0; i<=level; i++) { j+=nsep[i]; }
321827bd09bSSatish Balay 
322827bd09bSSatish Balay   /* storage for sparse x values */
323827bd09bSSatish Balay   n_global = xyt_handle->info->n_global;
324b1c944f5SJed Brown   xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes;
325a501084fSBarry Smith   x = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
326a501084fSBarry Smith   y = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
327827bd09bSSatish Balay 
328827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
329827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
330db4deed7SKarl Rupp   for (dim=i=j=0;i<m;i++) {
331827bd09bSSatish Balay     /* time to move to the next level? */
332d4af0d30SBarry Smith     while (i==segs[dim]) {
333e32f2f54SBarry Smith       if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n");
334827bd09bSSatish Balay       stages[dim++]=i;
335827bd09bSSatish Balay       end+=lnsep[dim];
336827bd09bSSatish Balay     }
337827bd09bSSatish Balay     stages[dim]=i;
338827bd09bSSatish Balay 
339827bd09bSSatish Balay     /* which column are we firing? */
340827bd09bSSatish Balay     /* i.e. set v_l */
341827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
342827bd09bSSatish Balay     (start<end) ? (col=fo[start]) : (col=INT_MAX);
343b1c944f5SJed Brown     PCTFS_giop_hc(&col,&work,1,op2,dim);
344827bd09bSSatish Balay 
345827bd09bSSatish Balay     /* shouldn't need this */
346db4deed7SKarl Rupp     if (col==INT_MAX) {
347f1ed62a8SBarry Smith       ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
348827bd09bSSatish Balay       continue;
349827bd09bSSatish Balay     }
350827bd09bSSatish Balay 
351827bd09bSSatish Balay     /* do I own it? I should */
352ca8e9878SJed Brown     PCTFS_rvec_zero(v ,a_m);
353827bd09bSSatish Balay     if (col==fo[start])
354827bd09bSSatish Balay     {
355827bd09bSSatish Balay       start++;
356ca8e9878SJed Brown       idx=PCTFS_ivec_linear_search(col, a_local2global, a_n);
357db4deed7SKarl Rupp       if (idx!=-1) { v[idx] = 1.0; j++; }
358c1235816SBarry Smith       else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");
359db4deed7SKarl Rupp     } else {
360ca8e9878SJed Brown       idx=PCTFS_ivec_linear_search(col, a_local2global, a_m);
361db4deed7SKarl Rupp       if (idx!=-1) {v[idx] = 1.0;}
362827bd09bSSatish Balay     }
363827bd09bSSatish Balay 
364827bd09bSSatish Balay     /* perform u = A.v_l */
365ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
366827bd09bSSatish Balay     do_matvec(xyt_handle->mvi,v,u);
367827bd09bSSatish Balay 
368827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
369827bd09bSSatish Balay     /* technically only need to zero out first i entries */
370827bd09bSSatish Balay     /* later turn this into an XYT_solve call ? */
371ca8e9878SJed Brown     PCTFS_rvec_zero(uu,m);
372827bd09bSSatish Balay     y_ptr=y;
373827bd09bSSatish Balay     iptr = ycol_indices;
374db4deed7SKarl Rupp     for (k=0; k<i; k++) {
375827bd09bSSatish Balay       off = *iptr++;
376827bd09bSSatish Balay       len = *iptr++;
3770805154bSBarry Smith       dlen = PetscBLASIntCast(len);
3784a0de3f6SBarry Smith       uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1);
379827bd09bSSatish Balay       y_ptr+=len;
380827bd09bSSatish Balay     }
381827bd09bSSatish Balay 
382827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
383b1c944f5SJed Brown     PCTFS_ssgl_radd  (uu, w, dim, stages);
384827bd09bSSatish Balay 
385827bd09bSSatish Balay     /* z = X.uu */
386ca8e9878SJed Brown     PCTFS_rvec_zero(z,n);
387827bd09bSSatish Balay     x_ptr=x;
388827bd09bSSatish Balay     iptr = xcol_indices;
389db4deed7SKarl Rupp     for (k=0; k<i; k++) {
390827bd09bSSatish Balay       off = *iptr++;
391827bd09bSSatish Balay       len = *iptr++;
3920805154bSBarry Smith       dlen = PetscBLASIntCast(len);
3934a0de3f6SBarry Smith       BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
394827bd09bSSatish Balay       x_ptr+=len;
395827bd09bSSatish Balay     }
396827bd09bSSatish Balay 
397827bd09bSSatish Balay     /* compute v_l = v_l - z */
398ca8e9878SJed Brown     PCTFS_rvec_zero(v+a_n,a_m-a_n);
3990805154bSBarry Smith     dlen = PetscBLASIntCast(n);
4004a0de3f6SBarry Smith     BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);
401827bd09bSSatish Balay 
402827bd09bSSatish Balay     /* compute u_l = A.v_l */
403db4deed7SKarl Rupp     if (a_n!=a_m) { PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); }
404ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
405827bd09bSSatish Balay     do_matvec(xyt_handle->mvi,v,u);
406827bd09bSSatish Balay 
407827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
4080805154bSBarry Smith     dlen = PetscBLASIntCast(n);
4094a0de3f6SBarry Smith     alpha = BLASdot_(&dlen,u,&i1,u,&i1);
410827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
411b1c944f5SJed Brown     PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
412827bd09bSSatish Balay 
4138f1a2a5eSBarry Smith     alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha);
414827bd09bSSatish Balay 
415827bd09bSSatish Balay     /* check for small alpha                             */
416827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
417e32f2f54SBarry Smith     if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);
418827bd09bSSatish Balay 
419827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
420ca8e9878SJed Brown     PCTFS_rvec_scale(v,1.0/alpha,n);
421ca8e9878SJed Brown     PCTFS_rvec_scale(u,1.0/alpha,n);
422827bd09bSSatish Balay 
423827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
424827bd09bSSatish Balay     flag = 1;
425827bd09bSSatish Balay     off=len=0;
426db4deed7SKarl Rupp     for (k=0; k<n; k++) {
427db4deed7SKarl Rupp       if (v[k]!=0.0) {
428827bd09bSSatish Balay         len=k;
429db4deed7SKarl Rupp         if (flag) {off=k; flag=0;}
430827bd09bSSatish Balay       }
431827bd09bSSatish Balay     }
432827bd09bSSatish Balay 
433827bd09bSSatish Balay     len -= (off-1);
434827bd09bSSatish Balay 
435db4deed7SKarl Rupp     if (len>0) {
436db4deed7SKarl Rupp       if ((xt_nnz+len)>xt_max_nnz) {
437f1ed62a8SBarry Smith         ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
438827bd09bSSatish Balay         xt_max_nnz *= 2;
439a501084fSBarry Smith         x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
440ca8e9878SJed Brown         PCTFS_rvec_copy(x_ptr,x,xt_nnz);
441a501084fSBarry Smith         free(x);
442827bd09bSSatish Balay         x = x_ptr;
443827bd09bSSatish Balay         x_ptr+=xt_nnz;
444827bd09bSSatish Balay       }
445827bd09bSSatish Balay       xt_nnz += len;
446ca8e9878SJed Brown       PCTFS_rvec_copy(x_ptr,v+off,len);
447827bd09bSSatish Balay 
448827bd09bSSatish Balay       /* keep track of number of zeros */
449db4deed7SKarl Rupp       if (dim) {
450db4deed7SKarl Rupp         for (k=0; k<len; k++) {
451db4deed7SKarl Rupp           if (x_ptr[k]==0.0) {xt_zero_nnz++;}
452827bd09bSSatish Balay         }
453db4deed7SKarl Rupp       } else {
454db4deed7SKarl Rupp         for (k=0; k<len; k++) {
455db4deed7SKarl Rupp           if (x_ptr[k]==0.0) { xt_zero_nnz_0++; }
456827bd09bSSatish Balay         }
457827bd09bSSatish Balay       }
458827bd09bSSatish Balay       xcol_indices[2*i] = off;
459827bd09bSSatish Balay       xcol_sz[i] = xcol_indices[2*i+1] = len;
460827bd09bSSatish Balay       xcol_vals[i] = x_ptr;
461db4deed7SKarl Rupp     } else {
462827bd09bSSatish Balay       xcol_indices[2*i] = 0;
463827bd09bSSatish Balay       xcol_sz[i] = xcol_indices[2*i+1] = 0;
464827bd09bSSatish Balay       xcol_vals[i] = x_ptr;
465827bd09bSSatish Balay     }
466827bd09bSSatish Balay 
467827bd09bSSatish Balay 
468827bd09bSSatish Balay     /* add newly generated column, u_l, to Y */
469827bd09bSSatish Balay     flag = 1;
470827bd09bSSatish Balay     off=len=0;
471db4deed7SKarl Rupp     for (k=0; k<n; k++) {
472db4deed7SKarl Rupp       if (u[k]!=0.0) {
473827bd09bSSatish Balay         len=k;
474db4deed7SKarl Rupp         if (flag) { off=k; flag=0; }
475827bd09bSSatish Balay       }
476827bd09bSSatish Balay     }
477827bd09bSSatish Balay 
478827bd09bSSatish Balay     len -= (off-1);
479827bd09bSSatish Balay 
480db4deed7SKarl Rupp     if (len>0) {
481db4deed7SKarl Rupp       if ((yt_nnz+len)>yt_max_nnz) {
482f1ed62a8SBarry Smith         ierr = PetscInfo(0,"increasing space for Y by 2x!\n");CHKERRQ(ierr);
483827bd09bSSatish Balay         yt_max_nnz *= 2;
484a501084fSBarry Smith         y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
485ca8e9878SJed Brown         PCTFS_rvec_copy(y_ptr,y,yt_nnz);
486a501084fSBarry Smith         free(y);
487827bd09bSSatish Balay         y = y_ptr;
488827bd09bSSatish Balay         y_ptr+=yt_nnz;
489827bd09bSSatish Balay       }
490827bd09bSSatish Balay       yt_nnz += len;
491ca8e9878SJed Brown       PCTFS_rvec_copy(y_ptr,u+off,len);
492827bd09bSSatish Balay 
493827bd09bSSatish Balay       /* keep track of number of zeros */
494db4deed7SKarl Rupp       if (dim) {
495db4deed7SKarl Rupp         for (k=0; k<len; k++) {
496db4deed7SKarl Rupp           if (y_ptr[k]==0.0) { yt_zero_nnz++; }
497827bd09bSSatish Balay         }
498db4deed7SKarl Rupp       } else {
499db4deed7SKarl Rupp         for (k=0; k<len; k++) {
500db4deed7SKarl Rupp           if (y_ptr[k]==0.0) { yt_zero_nnz_0++; }
501827bd09bSSatish Balay         }
502827bd09bSSatish Balay       }
503827bd09bSSatish Balay       ycol_indices[2*i] = off;
504827bd09bSSatish Balay       ycol_sz[i] = ycol_indices[2*i+1] = len;
505827bd09bSSatish Balay       ycol_vals[i] = y_ptr;
506db4deed7SKarl Rupp     } else {
507827bd09bSSatish Balay       ycol_indices[2*i] = 0;
508827bd09bSSatish Balay       ycol_sz[i] = ycol_indices[2*i+1] = 0;
509827bd09bSSatish Balay       ycol_vals[i] = y_ptr;
510827bd09bSSatish Balay     }
511827bd09bSSatish Balay   }
512827bd09bSSatish Balay 
513827bd09bSSatish Balay   /* close off stages for execution phase */
514db4deed7SKarl Rupp   while (dim!=level) {
515827bd09bSSatish Balay     stages[dim++]=i;
5164a0de3f6SBarry Smith     ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
517827bd09bSSatish Balay   }
518827bd09bSSatish Balay   stages[dim]=i;
519827bd09bSSatish Balay 
520827bd09bSSatish Balay   xyt_handle->info->n=xyt_handle->mvi->n;
521827bd09bSSatish Balay   xyt_handle->info->m=m;
522827bd09bSSatish Balay   xyt_handle->info->nnz=xt_nnz + yt_nnz;
523827bd09bSSatish Balay   xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz;
524827bd09bSSatish Balay   xyt_handle->info->msg_buf_sz=stages[level]-stages[0];
525a501084fSBarry Smith   xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
526a501084fSBarry Smith   xyt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
527827bd09bSSatish Balay   xyt_handle->info->x=x;
528827bd09bSSatish Balay   xyt_handle->info->xcol_vals=xcol_vals;
529827bd09bSSatish Balay   xyt_handle->info->xcol_sz=xcol_sz;
530827bd09bSSatish Balay   xyt_handle->info->xcol_indices=xcol_indices;
531827bd09bSSatish Balay   xyt_handle->info->stages=stages;
532827bd09bSSatish Balay   xyt_handle->info->y=y;
533827bd09bSSatish Balay   xyt_handle->info->ycol_vals=ycol_vals;
534827bd09bSSatish Balay   xyt_handle->info->ycol_sz=ycol_sz;
535827bd09bSSatish Balay   xyt_handle->info->ycol_indices=ycol_indices;
536827bd09bSSatish Balay 
537a501084fSBarry Smith   free(segs);
538a501084fSBarry Smith   free(u);
539a501084fSBarry Smith   free(v);
540a501084fSBarry Smith   free(uu);
541a501084fSBarry Smith   free(z);
542a501084fSBarry Smith   free(w);
543827bd09bSSatish Balay 
544827bd09bSSatish Balay   return(0);
545827bd09bSSatish Balay }
546827bd09bSSatish Balay 
5477b1ae94cSBarry Smith /**************************************xyt.c***********************************/
5480924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle,  PetscScalar *uc)
549827bd09bSSatish Balay {
55052f87cdaSBarry Smith   PetscInt off, len, *iptr;
55152f87cdaSBarry Smith   PetscInt level       =xyt_handle->level;
55252f87cdaSBarry Smith   PetscInt n           =xyt_handle->info->n;
55352f87cdaSBarry Smith   PetscInt m           =xyt_handle->info->m;
55452f87cdaSBarry Smith   PetscInt *stages     =xyt_handle->info->stages;
55552f87cdaSBarry Smith   PetscInt *xcol_indices=xyt_handle->info->xcol_indices;
55652f87cdaSBarry Smith   PetscInt *ycol_indices=xyt_handle->info->ycol_indices;
557a501084fSBarry Smith   PetscScalar *x_ptr, *y_ptr, *uu_ptr;
558a501084fSBarry Smith   PetscScalar *solve_uu=xyt_handle->info->solve_uu;
559a501084fSBarry Smith   PetscScalar *solve_w =xyt_handle->info->solve_w;
560a501084fSBarry Smith   PetscScalar *x       =xyt_handle->info->x;
561a501084fSBarry Smith   PetscScalar *y       =xyt_handle->info->y;
5624a0de3f6SBarry Smith   PetscBLASInt i1 = 1,dlen;
563827bd09bSSatish Balay 
5640924e98cSBarry Smith   PetscFunctionBegin;
565827bd09bSSatish Balay   uu_ptr=solve_uu;
566ca8e9878SJed Brown   PCTFS_rvec_zero(uu_ptr,m);
567827bd09bSSatish Balay 
568827bd09bSSatish Balay   /* x  = X.Y^T.b */
569827bd09bSSatish Balay   /* uu = Y^T.b */
570827bd09bSSatish Balay   for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len)
571827bd09bSSatish Balay   {
572827bd09bSSatish Balay     off=*iptr++; len=*iptr++;
5730805154bSBarry Smith     dlen = PetscBLASIntCast(len);
5744a0de3f6SBarry Smith     *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1);
575827bd09bSSatish Balay   }
576827bd09bSSatish Balay 
577827bd09bSSatish Balay   /* comunication of beta */
578827bd09bSSatish Balay   uu_ptr=solve_uu;
579b1c944f5SJed Brown   if (level) { PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages); }
580827bd09bSSatish Balay 
581ca8e9878SJed Brown   PCTFS_rvec_zero(uc,n);
582827bd09bSSatish Balay 
583827bd09bSSatish Balay   /* x = X.uu */
584db4deed7SKarl Rupp   for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) {
585827bd09bSSatish Balay     off=*iptr++; len=*iptr++;
5860805154bSBarry Smith     dlen = PetscBLASIntCast(len);
5874a0de3f6SBarry Smith     BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
588827bd09bSSatish Balay   }
5890924e98cSBarry Smith   PetscFunctionReturn(0);
590827bd09bSSatish Balay }
591827bd09bSSatish Balay 
5927b1ae94cSBarry Smith /**************************************xyt.c***********************************/
5930924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle)
594827bd09bSSatish Balay {
59552f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
596827bd09bSSatish Balay 
5970924e98cSBarry Smith   PetscFunctionBegin;
598e32f2f54SBarry Smith   if (!xyt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle);
599827bd09bSSatish Balay 
600827bd09bSSatish Balay   vals[0]=vals[1]=xyt_handle->id;
601b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
602e32f2f54SBarry Smith   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);
6030924e98cSBarry Smith   PetscFunctionReturn(0);
604827bd09bSSatish Balay }
605827bd09bSSatish Balay 
6067b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6070924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle)
608827bd09bSSatish Balay {
60952f87cdaSBarry Smith   PetscInt i, ct, id;
61052f87cdaSBarry Smith   PetscInt mask, edge, *iptr;
61152f87cdaSBarry Smith   PetscInt *dir, *used;
61252f87cdaSBarry Smith   PetscInt sum[4], w[4];
613a501084fSBarry Smith   PetscScalar rsum[4], rw[4];
61452f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
615a501084fSBarry Smith   PetscScalar *lhs, *rhs;
61652f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo=0;
617ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle;
61852f87cdaSBarry Smith   PetscInt *local2global=xyt_handle->mvi->local2global;
61952f87cdaSBarry Smith   PetscInt  n=xyt_handle->mvi->n;
62052f87cdaSBarry Smith   PetscInt  m=xyt_handle->mvi->m;
62152f87cdaSBarry Smith   PetscInt level=xyt_handle->level;
622ab824b78SBarry Smith   PetscInt shared=0;
623d1528f56SBarry Smith   PetscErrorCode ierr;
624827bd09bSSatish Balay 
6250924e98cSBarry Smith   PetscFunctionBegin;
62652f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
62752f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
62852f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
62952f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
63052f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
631827bd09bSSatish Balay 
632ca8e9878SJed Brown   PCTFS_ivec_zero(dir  ,level+1);
633ca8e9878SJed Brown   PCTFS_ivec_zero(nsep ,level+1);
634ca8e9878SJed Brown   PCTFS_ivec_zero(lnsep,level+1);
635ca8e9878SJed Brown   PCTFS_ivec_set (fo   ,-1,n+1);
636ca8e9878SJed Brown   PCTFS_ivec_zero(used,n);
637827bd09bSSatish Balay 
6388cda6cd7SBarry Smith   lhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
6398cda6cd7SBarry Smith   rhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
640827bd09bSSatish Balay 
641827bd09bSSatish Balay   /* determine the # of unique dof */
642ca8e9878SJed Brown   PCTFS_rvec_zero(lhs,m);
643ca8e9878SJed Brown   PCTFS_rvec_set(lhs,1.0,n);
644ca8e9878SJed Brown   PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level);
645ca8e9878SJed Brown   ierr = PetscInfo(0,"done first PCTFS_gs_gop_hc\n");CHKERRQ(ierr);
646ca8e9878SJed Brown   PCTFS_rvec_zero(rsum,2);
647827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
648827bd09bSSatish Balay   {
649db4deed7SKarl Rupp     if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; }
650827bd09bSSatish Balay 
651db4deed7SKarl Rupp     if (lhs[i]!=1.0) {
652ab824b78SBarry Smith       shared=1;
653827bd09bSSatish Balay     }
654827bd09bSSatish Balay   }
655827bd09bSSatish Balay 
656b1c944f5SJed Brown   PCTFS_grop_hc(rsum,rw,2,op,level);
657827bd09bSSatish Balay   rsum[0]+=0.1;
658827bd09bSSatish Balay   rsum[1]+=0.1;
659827bd09bSSatish Balay 
66052f87cdaSBarry Smith   xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0];
66152f87cdaSBarry Smith   xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0];
662827bd09bSSatish Balay 
663827bd09bSSatish Balay   /* determine separator sets top down */
664827bd09bSSatish Balay   if (shared)
665827bd09bSSatish Balay   {
666827bd09bSSatish Balay     /* solution is to do as in the symmetric shared case but then */
667827bd09bSSatish Balay     /* pick the sub-hc with the most free dofs and do a mat-vec   */
668827bd09bSSatish Balay     /* and pick up the responses on the other sub-hc from the     */
669827bd09bSSatish Balay     /* initial separator set obtained from the symm. shared case  */
670e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n");
671db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
672db4deed7SKarl Rupp 
673827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
674ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
675ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
676827bd09bSSatish Balay 
677827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
678ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
679ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
680827bd09bSSatish Balay 
681db4deed7SKarl Rupp       for (i=0;i<n;i++) {
682db4deed7SKarl Rupp         if (id< mask) {
683db4deed7SKarl Rupp           if (lhs[i]!=0.0) { lhs[i]=1.0; }
684827bd09bSSatish Balay         }
685db4deed7SKarl Rupp         if (id>=mask) {
686db4deed7SKarl Rupp           if (rhs[i]!=0.0) { rhs[i]=1.0; }
687827bd09bSSatish Balay         }
688827bd09bSSatish Balay       }
689827bd09bSSatish Balay 
690db4deed7SKarl Rupp       if (id< mask) { PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1); }
691db4deed7SKarl Rupp       else { PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1); }
692827bd09bSSatish Balay 
693827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
694ca8e9878SJed Brown       PCTFS_rvec_zero(rsum,4);
695db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
696db4deed7SKarl Rupp         if (!used[i]) {
697827bd09bSSatish Balay           /* number of unmarked dofs on node */
698827bd09bSSatish Balay           ct++;
699827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
700db4deed7SKarl Rupp           if (id< mask) {
701db4deed7SKarl Rupp             if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; }
702827bd09bSSatish Balay           }
703827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
704db4deed7SKarl Rupp           if (id>=mask) {
705db4deed7SKarl Rupp             if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; }
706827bd09bSSatish Balay           }
707827bd09bSSatish Balay         }
708827bd09bSSatish Balay       }
709827bd09bSSatish Balay 
710827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
711827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
712827bd09bSSatish Balay       (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
713b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
714b1c944f5SJed Brown       PCTFS_grop_hc(rsum,rw,4,op,edge);
715827bd09bSSatish Balay       rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
716827bd09bSSatish Balay 
717db4deed7SKarl Rupp       if (id<mask) {
718827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
719db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
720db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
721827bd09bSSatish Balay             ct++; nfo++;
722827bd09bSSatish Balay 
723c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
724827bd09bSSatish Balay 
725827bd09bSSatish Balay             *--iptr = local2global[i];
726827bd09bSSatish Balay             used[i]=edge;
727827bd09bSSatish Balay           }
728827bd09bSSatish Balay         }
729ca8e9878SJed Brown         if (ct>1) {PCTFS_ivec_sort(iptr,ct);}
730827bd09bSSatish Balay 
731827bd09bSSatish Balay         lnsep[edge]=ct;
73252f87cdaSBarry Smith         nsep[edge]=(PetscInt) rsum[0];
733827bd09bSSatish Balay         dir [edge]=LEFT;
734827bd09bSSatish Balay       }
735827bd09bSSatish Balay 
736db4deed7SKarl Rupp       if (id>=mask) {
737827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
738db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
739db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
740827bd09bSSatish Balay             ct++; nfo++;
741827bd09bSSatish Balay 
742c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
743827bd09bSSatish Balay 
744827bd09bSSatish Balay             *--iptr = local2global[i];
745827bd09bSSatish Balay             used[i]=edge;
746827bd09bSSatish Balay           }
747827bd09bSSatish Balay         }
748ca8e9878SJed Brown         if (ct>1) { PCTFS_ivec_sort(iptr,ct); }
749827bd09bSSatish Balay 
750827bd09bSSatish Balay         lnsep[edge]=ct;
75152f87cdaSBarry Smith         nsep[edge]= (PetscInt) rsum[1];
752827bd09bSSatish Balay         dir [edge]=RIGHT;
753827bd09bSSatish Balay       }
754827bd09bSSatish Balay 
755827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
756827bd09bSSatish Balay       /* do we need full set of separators for this?                */
757827bd09bSSatish Balay 
758827bd09bSSatish Balay       /* fold rhs hc into lower */
759db4deed7SKarl Rupp       if (id>=mask) { id-=mask; }
760827bd09bSSatish Balay     }
761db4deed7SKarl Rupp   } else {
762db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
763827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
764ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
765ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
766827bd09bSSatish Balay 
767827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
768ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
769ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
770827bd09bSSatish Balay 
771827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
772db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
773db4deed7SKarl Rupp         if (!used[i]) {
774827bd09bSSatish Balay           /* number of unmarked dofs on node */
775827bd09bSSatish Balay           ct++;
776827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
777827bd09bSSatish Balay           if ((id< mask)&&(lhs[i]!=0.0)) { sum[0]++; }
778827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
779827bd09bSSatish Balay           if ((id>=mask)&&(rhs[i]!=0.0)) { sum[1]++; }
780827bd09bSSatish Balay         }
781827bd09bSSatish Balay       }
782827bd09bSSatish Balay 
783827bd09bSSatish Balay       /* for the non-symmetric case we need separators of width 2 */
784827bd09bSSatish Balay       /* so take both sides */
785827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
786b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
787827bd09bSSatish Balay 
788827bd09bSSatish Balay       ct=0;
789db4deed7SKarl Rupp       if (id<mask) {
790827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
791db4deed7SKarl Rupp         for (i=0;i<n;i++) {
792db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
793827bd09bSSatish Balay             ct++; nfo++;
794827bd09bSSatish Balay             *--iptr = local2global[i];
795827bd09bSSatish Balay             used[i]=edge;
796827bd09bSSatish Balay           }
797827bd09bSSatish Balay         }
798827bd09bSSatish Balay         /* LSH hc summation of ct should be sum[0] */
799db4deed7SKarl Rupp       } else {
800827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
801db4deed7SKarl Rupp         for (i=0;i<n;i++) {
802db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
803827bd09bSSatish Balay             ct++; nfo++;
804827bd09bSSatish Balay             *--iptr = local2global[i];
805827bd09bSSatish Balay             used[i]=edge;
806827bd09bSSatish Balay           }
807827bd09bSSatish Balay         }
808827bd09bSSatish Balay         /* RSH hc summation of ct should be sum[1] */
809827bd09bSSatish Balay       }
810827bd09bSSatish Balay 
811ca8e9878SJed Brown       if (ct>1) { PCTFS_ivec_sort(iptr,ct); }
812827bd09bSSatish Balay       lnsep[edge]=ct;
813827bd09bSSatish Balay       nsep[edge]=sum[0]+sum[1];
814827bd09bSSatish Balay       dir [edge]=BOTH;
815827bd09bSSatish Balay 
816827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
817827bd09bSSatish Balay       /* do we need full set of separators for this?                */
818827bd09bSSatish Balay 
819827bd09bSSatish Balay       /* fold rhs hc into lower */
820827bd09bSSatish Balay       if (id>=mask)
821827bd09bSSatish Balay         {id-=mask;}
822827bd09bSSatish Balay     }
823827bd09bSSatish Balay   }
824827bd09bSSatish Balay 
825827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
826db4deed7SKarl Rupp   for (ct=i=0;i<n;i++) {
827db4deed7SKarl Rupp     if (!used[i]) {
828827bd09bSSatish Balay       ct++; nfo++;
829827bd09bSSatish Balay       *--iptr = local2global[i];
830827bd09bSSatish Balay       used[i]=edge;
831827bd09bSSatish Balay     }
832827bd09bSSatish Balay   }
833ca8e9878SJed Brown   if (ct>1) {PCTFS_ivec_sort(iptr,ct);}
834827bd09bSSatish Balay   lnsep[edge]=ct;
835827bd09bSSatish Balay   nsep [edge]=ct;
836827bd09bSSatish Balay   dir  [edge]=BOTH;
837827bd09bSSatish Balay 
838827bd09bSSatish Balay   xyt_handle->info->nsep=nsep;
839827bd09bSSatish Balay   xyt_handle->info->lnsep=lnsep;
840827bd09bSSatish Balay   xyt_handle->info->fo=fo;
841827bd09bSSatish Balay   xyt_handle->info->nfo=nfo;
842827bd09bSSatish Balay 
843a501084fSBarry Smith   free(dir);
844a501084fSBarry Smith   free(lhs);
845a501084fSBarry Smith   free(rhs);
846a501084fSBarry Smith   free(used);
8470924e98cSBarry Smith   PetscFunctionReturn(0);
848827bd09bSSatish Balay }
849827bd09bSSatish Balay 
8507b1ae94cSBarry Smith /**************************************xyt.c***********************************/
851*5c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data)
852827bd09bSSatish Balay {
853827bd09bSSatish Balay   mv_info *mvi;
854827bd09bSSatish Balay 
855827bd09bSSatish Balay 
856a501084fSBarry Smith   mvi = (mv_info*)malloc(sizeof(mv_info));
857827bd09bSSatish Balay   mvi->n=n;
858827bd09bSSatish Balay   mvi->m=m;
859827bd09bSSatish Balay   mvi->n_global=-1;
860827bd09bSSatish Balay   mvi->m_global=-1;
86152f87cdaSBarry Smith   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
862ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global,local2global,m);
863827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
864*5c8f6a95SKarl Rupp   mvi->matvec=matvec;
865827bd09bSSatish Balay   mvi->grid_data=grid_data;
866827bd09bSSatish Balay 
867827bd09bSSatish Balay   /* set xyt communication handle to perform restricted matvec */
868ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
869827bd09bSSatish Balay 
870827bd09bSSatish Balay   return(mvi);
871827bd09bSSatish Balay }
872827bd09bSSatish Balay 
8737b1ae94cSBarry Smith /**************************************xyt.c***********************************/
8743fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
875827bd09bSSatish Balay {
8760924e98cSBarry Smith   PetscFunctionBegin;
877827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
8780924e98cSBarry Smith   PetscFunctionReturn(0);
879827bd09bSSatish Balay }
880827bd09bSSatish Balay 
881827bd09bSSatish Balay 
882827bd09bSSatish Balay 
883