xref: /petsc/src/ksp/pc/impls/tfs/xyt.c (revision 7d0a6c19129e7069c8a40e210b34ed62989173db)
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***********************************/
197c4f633dSBarry Smith #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;
41827bd09bSSatish Balay   gs_ADT 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);
6452f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, 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              */
86827bd09bSSatish Balay 	   void *matvec,       /* b_loc=A_local.x_loc         */
87827bd09bSSatish Balay 	   void *grid_data     /* grid data for matvec        */
88827bd09bSSatish Balay 	   )
89827bd09bSSatish Balay {
90827bd09bSSatish Balay 
91a40341efSBarry Smith   comm_init();
92827bd09bSSatish Balay   check_handle(xyt_handle);
93827bd09bSSatish Balay 
94827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
95c1235816SBarry Smith   if ((1<<(xyt_handle->level=i_log2_num_nodes))!=num_nodes) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %D != %D\n",1<<i_log2_num_nodes,num_nodes);
96827bd09bSSatish Balay 
97827bd09bSSatish Balay   /* space for X info */
98a501084fSBarry Smith   xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));
99827bd09bSSatish Balay 
100827bd09bSSatish Balay   /* set up matvec handles */
101827bd09bSSatish Balay   xyt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);
102827bd09bSSatish Balay 
103827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
104827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
105827bd09bSSatish Balay   xyt_handle->ns=0;
106827bd09bSSatish Balay 
107827bd09bSSatish Balay   /* determine separators and generate firing order - NB xyt info set here */
108827bd09bSSatish Balay   det_separators(xyt_handle);
109827bd09bSSatish Balay 
110827bd09bSSatish Balay   return(do_xyt_factor(xyt_handle));
111827bd09bSSatish Balay }
112827bd09bSSatish Balay 
1137b1ae94cSBarry Smith /**************************************xyt.c***********************************/
1148cda6cd7SBarry Smith PetscInt XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b)
115827bd09bSSatish Balay {
116a40341efSBarry Smith   comm_init();
117827bd09bSSatish Balay   check_handle(xyt_handle);
118827bd09bSSatish Balay 
119827bd09bSSatish Balay   /* need to copy b into x? */
120827bd09bSSatish Balay   if (b)
121827bd09bSSatish Balay     {rvec_copy(x,b,xyt_handle->mvi->n);}
122827bd09bSSatish Balay   do_xyt_solve(xyt_handle,x);
123827bd09bSSatish Balay 
124827bd09bSSatish Balay   return(0);
125827bd09bSSatish Balay }
126827bd09bSSatish Balay 
1277b1ae94cSBarry Smith /**************************************xyt.c***********************************/
12852f87cdaSBarry Smith PetscInt XYT_free(xyt_ADT xyt_handle)
129827bd09bSSatish Balay {
130a40341efSBarry Smith   comm_init();
131827bd09bSSatish Balay   check_handle(xyt_handle);
132827bd09bSSatish Balay   n_xyt_handles--;
133827bd09bSSatish Balay 
134a501084fSBarry Smith   free(xyt_handle->info->nsep);
135a501084fSBarry Smith   free(xyt_handle->info->lnsep);
136a501084fSBarry Smith   free(xyt_handle->info->fo);
137a501084fSBarry Smith   free(xyt_handle->info->stages);
138a501084fSBarry Smith   free(xyt_handle->info->solve_uu);
139a501084fSBarry Smith   free(xyt_handle->info->solve_w);
140a501084fSBarry Smith   free(xyt_handle->info->x);
141a501084fSBarry Smith   free(xyt_handle->info->xcol_vals);
142a501084fSBarry Smith   free(xyt_handle->info->xcol_sz);
143a501084fSBarry Smith   free(xyt_handle->info->xcol_indices);
144a501084fSBarry Smith   free(xyt_handle->info->y);
145a501084fSBarry Smith   free(xyt_handle->info->ycol_vals);
146a501084fSBarry Smith   free(xyt_handle->info->ycol_sz);
147a501084fSBarry Smith   free(xyt_handle->info->ycol_indices);
148a501084fSBarry Smith   free(xyt_handle->info);
149a501084fSBarry Smith   free(xyt_handle->mvi->local2global);
150827bd09bSSatish Balay   gs_free(xyt_handle->mvi->gs_handle);
151a501084fSBarry Smith   free(xyt_handle->mvi);
152a501084fSBarry Smith   free(xyt_handle);
153827bd09bSSatish Balay 
154827bd09bSSatish Balay 
155827bd09bSSatish Balay   /* if the check fails we nuke */
156a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
157827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
158827bd09bSSatish Balay   return(0);
159827bd09bSSatish Balay }
160827bd09bSSatish Balay 
1617b1ae94cSBarry Smith /**************************************xyt.c***********************************/
16252f87cdaSBarry Smith PetscInt XYT_stats(xyt_ADT xyt_handle)
163827bd09bSSatish Balay {
16452f87cdaSBarry Smith   PetscInt    op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
16552f87cdaSBarry Smith   PetscInt   fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
16652f87cdaSBarry Smith   PetscInt    vals[9],  work[9];
167a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
168827bd09bSSatish Balay 
169a40341efSBarry Smith   comm_init();
170827bd09bSSatish Balay   check_handle(xyt_handle);
171827bd09bSSatish Balay 
172827bd09bSSatish Balay   /* if factorization not done there are no stats */
173*7d0a6c19SBarry Smith   if (!xyt_handle->info||!xyt_handle->mvi) {
174*7d0a6c19SBarry Smith     if (!my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n");
175827bd09bSSatish Balay     return 1;
176827bd09bSSatish Balay   }
177827bd09bSSatish Balay 
178827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
179827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
180827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
181827bd09bSSatish Balay   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
182827bd09bSSatish Balay 
183827bd09bSSatish Balay   fvals[0]=fvals[1]=fvals[2]
184827bd09bSSatish Balay     =xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
185827bd09bSSatish Balay   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
186827bd09bSSatish Balay 
187*7d0a6c19SBarry Smith   if (!my_id) {
188*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_nnz=%D\n",my_id,vals[0]);
189*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_nnz=%D\n",my_id,vals[1]);
190*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
191*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_nnz=%D\n",my_id,vals[2]);
192*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
193*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
194*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_n  =%D\n",my_id,vals[3]);
195*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_n  =%D\n",my_id,vals[4]);
196*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
197*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_n  =%D\n",my_id,vals[5]);
198*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_buf=%D\n",my_id,vals[6]);
199*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_buf=%D\n",my_id,vals[7]);
200*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
201*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_slv=%g\n",my_id,fvals[0]);
202*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_slv=%g\n",my_id,fvals[1]);
203*7d0a6c19SBarry Smith     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_slv=%g\n",my_id,fvals[2]/num_nodes);
204827bd09bSSatish Balay   }
205827bd09bSSatish Balay 
206827bd09bSSatish Balay   return(0);
207827bd09bSSatish Balay }
208827bd09bSSatish Balay 
209827bd09bSSatish Balay 
210827bd09bSSatish Balay /*************************************xyt.c************************************
211827bd09bSSatish Balay 
212827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
213827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
214827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
215827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
216827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
217827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
218827bd09bSSatish Balay    gs_init/gop).
219827bd09bSSatish Balay 
220827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
221827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
222827bd09bSSatish Balay **************************************xyt.c***********************************/
22352f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle)
224827bd09bSSatish Balay {
2257b1ae94cSBarry Smith   return xyt_generate(xyt_handle);
226827bd09bSSatish Balay }
227827bd09bSSatish Balay 
2287b1ae94cSBarry Smith /**************************************xyt.c***********************************/
22952f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle)
230827bd09bSSatish Balay {
23152f87cdaSBarry Smith   PetscInt i,j,k,idx;
23252f87cdaSBarry Smith   PetscInt dim, col;
233a501084fSBarry Smith   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
23452f87cdaSBarry Smith   PetscInt *segs;
23552f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
23652f87cdaSBarry Smith   PetscInt off, len;
237a501084fSBarry Smith   PetscScalar *x_ptr, *y_ptr;
23852f87cdaSBarry Smith   PetscInt *iptr, flag;
23952f87cdaSBarry Smith   PetscInt start=0, end, work;
24052f87cdaSBarry Smith   PetscInt op2[] = {GL_MIN,0};
241827bd09bSSatish Balay   gs_ADT gs_handle;
24252f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo;
24352f87cdaSBarry Smith   PetscInt a_n=xyt_handle->mvi->n;
24452f87cdaSBarry Smith   PetscInt a_m=xyt_handle->mvi->m;
24552f87cdaSBarry Smith   PetscInt *a_local2global=xyt_handle->mvi->local2global;
24652f87cdaSBarry Smith   PetscInt level;
24752f87cdaSBarry Smith   PetscInt n, m;
24852f87cdaSBarry Smith   PetscInt *xcol_sz, *xcol_indices, *stages;
249a501084fSBarry Smith   PetscScalar **xcol_vals, *x;
25052f87cdaSBarry Smith   PetscInt *ycol_sz, *ycol_indices;
251a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
25252f87cdaSBarry Smith   PetscInt n_global;
25352f87cdaSBarry Smith   PetscInt xt_nnz=0, xt_max_nnz=0;
25452f87cdaSBarry Smith   PetscInt yt_nnz=0, yt_max_nnz=0;
25552f87cdaSBarry Smith   PetscInt xt_zero_nnz  =0;
25652f87cdaSBarry Smith   PetscInt xt_zero_nnz_0=0;
25752f87cdaSBarry Smith   PetscInt yt_zero_nnz  =0;
25852f87cdaSBarry Smith   PetscInt yt_zero_nnz_0=0;
2594a0de3f6SBarry Smith   PetscBLASInt i1 = 1,dlen;
260a501084fSBarry Smith   PetscScalar dm1 = -1.0;
261d1528f56SBarry Smith   PetscErrorCode ierr;
262827bd09bSSatish Balay 
263827bd09bSSatish Balay   n=xyt_handle->mvi->n;
264827bd09bSSatish Balay   nsep=xyt_handle->info->nsep;
265827bd09bSSatish Balay   lnsep=xyt_handle->info->lnsep;
266827bd09bSSatish Balay   fo=xyt_handle->info->fo;
267827bd09bSSatish Balay   end=lnsep[0];
268827bd09bSSatish Balay   level=xyt_handle->level;
269827bd09bSSatish Balay   gs_handle=xyt_handle->mvi->gs_handle;
270827bd09bSSatish Balay 
271827bd09bSSatish Balay   /* is there a null space? */
272827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
273827bd09bSSatish Balay   for (i=0, j=0; i<=level; i++)
274827bd09bSSatish Balay     {j+=nsep[i];}
275827bd09bSSatish Balay 
276827bd09bSSatish Balay   m = j-xyt_handle->ns;
277827bd09bSSatish Balay   if (m!=j)
2784a0de3f6SBarry Smith     {ierr = PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);}
279827bd09bSSatish Balay 
2804a0de3f6SBarry Smith   ierr = PetscInfo2(0,"xyt_generate() :: X(%D,%D)\n",n,m);CHKERRQ(ierr);
281827bd09bSSatish Balay 
282827bd09bSSatish Balay   /* get and initialize storage for x local         */
283827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
28452f87cdaSBarry Smith   xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
28552f87cdaSBarry Smith   xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
286a501084fSBarry Smith   xcol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
287827bd09bSSatish Balay   for (i=j=0; i<m; i++, j+=2)
288827bd09bSSatish Balay     {
289827bd09bSSatish Balay       xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
290827bd09bSSatish Balay       xcol_vals[i] = NULL;
291827bd09bSSatish Balay     }
292827bd09bSSatish Balay   xcol_indices[j]=-1;
293827bd09bSSatish Balay 
294827bd09bSSatish Balay   /* get and initialize storage for y local         */
295827bd09bSSatish Balay   /* note that y local is nxm and stored by columns */
29652f87cdaSBarry Smith   ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
29752f87cdaSBarry Smith   ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
298a501084fSBarry Smith   ycol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
299827bd09bSSatish Balay   for (i=j=0; i<m; i++, j+=2)
300827bd09bSSatish Balay     {
301827bd09bSSatish Balay       ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
302827bd09bSSatish Balay       ycol_vals[i] = NULL;
303827bd09bSSatish Balay     }
304827bd09bSSatish Balay   ycol_indices[j]=-1;
305827bd09bSSatish Balay 
306827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
307827bd09bSSatish Balay   /* this looks like nsep[]=segments */
30852f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
30952f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
310827bd09bSSatish Balay   ivec_zero(stages,level+1);
311827bd09bSSatish Balay   ivec_copy(segs,nsep,level+1);
312827bd09bSSatish Balay   for (i=0; i<level; i++)
313827bd09bSSatish Balay     {segs[i+1] += segs[i];}
314827bd09bSSatish Balay   stages[0] = segs[0];
315827bd09bSSatish Balay 
316827bd09bSSatish Balay   /* temporary vectors  */
317a501084fSBarry Smith   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
318a501084fSBarry Smith   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
319a501084fSBarry Smith   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
320a501084fSBarry Smith   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
321a501084fSBarry Smith   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
322827bd09bSSatish Balay 
323827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
324827bd09bSSatish Balay   for (i=1, j=0; i<=level; i++)
325827bd09bSSatish Balay     {j+=nsep[i];}
326827bd09bSSatish Balay 
327827bd09bSSatish Balay   /* storage for sparse x values */
328827bd09bSSatish Balay   n_global = xyt_handle->info->n_global;
32952f87cdaSBarry Smith   xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
330a501084fSBarry Smith   x = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
331a501084fSBarry Smith   y = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
332827bd09bSSatish Balay 
333827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
334827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
335827bd09bSSatish Balay   for (dim=i=j=0;i<m;i++)
336827bd09bSSatish Balay     {
337827bd09bSSatish Balay       /* time to move to the next level? */
338d4af0d30SBarry Smith       while (i==segs[dim]){
339e32f2f54SBarry Smith         if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n");
340827bd09bSSatish Balay 	stages[dim++]=i;
341827bd09bSSatish Balay 	end+=lnsep[dim];
342827bd09bSSatish Balay       }
343827bd09bSSatish Balay       stages[dim]=i;
344827bd09bSSatish Balay 
345827bd09bSSatish Balay       /* which column are we firing? */
346827bd09bSSatish Balay       /* i.e. set v_l */
347827bd09bSSatish Balay       /* use new seps and do global min across hc to determine which one to fire */
348827bd09bSSatish Balay       (start<end) ? (col=fo[start]) : (col=INT_MAX);
349827bd09bSSatish Balay       giop_hc(&col,&work,1,op2,dim);
350827bd09bSSatish Balay 
351827bd09bSSatish Balay       /* shouldn't need this */
352827bd09bSSatish Balay       if (col==INT_MAX)
353827bd09bSSatish Balay 	{
354f1ed62a8SBarry Smith 	  ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
355827bd09bSSatish Balay 	  continue;
356827bd09bSSatish Balay 	}
357827bd09bSSatish Balay 
358827bd09bSSatish Balay       /* do I own it? I should */
359827bd09bSSatish Balay       rvec_zero(v ,a_m);
360827bd09bSSatish Balay       if (col==fo[start])
361827bd09bSSatish Balay 	{
362827bd09bSSatish Balay 	  start++;
363827bd09bSSatish Balay 	  idx=ivec_linear_search(col, a_local2global, a_n);
364827bd09bSSatish Balay 	  if (idx!=-1)
365827bd09bSSatish Balay 	    {v[idx] = 1.0; j++;}
366c1235816SBarry Smith 	  else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");
367827bd09bSSatish Balay 	}
368827bd09bSSatish Balay       else
369827bd09bSSatish Balay 	{
370827bd09bSSatish Balay 	  idx=ivec_linear_search(col, a_local2global, a_m);
371827bd09bSSatish Balay 	  if (idx!=-1)
372827bd09bSSatish Balay 	    {v[idx] = 1.0;}
373827bd09bSSatish Balay 	}
374827bd09bSSatish Balay 
375827bd09bSSatish Balay       /* perform u = A.v_l */
376827bd09bSSatish Balay       rvec_zero(u,n);
377827bd09bSSatish Balay       do_matvec(xyt_handle->mvi,v,u);
378827bd09bSSatish Balay 
379827bd09bSSatish Balay       /* uu =  X^T.u_l (local portion) */
380827bd09bSSatish Balay       /* technically only need to zero out first i entries */
381827bd09bSSatish Balay       /* later turn this into an XYT_solve call ? */
382827bd09bSSatish Balay       rvec_zero(uu,m);
383827bd09bSSatish Balay       y_ptr=y;
384827bd09bSSatish Balay       iptr = ycol_indices;
385827bd09bSSatish Balay       for (k=0; k<i; k++)
386827bd09bSSatish Balay 	{
387827bd09bSSatish Balay 	  off = *iptr++;
388827bd09bSSatish Balay 	  len = *iptr++;
3890805154bSBarry Smith           dlen = PetscBLASIntCast(len);
3904a0de3f6SBarry Smith 	  uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1);
391827bd09bSSatish Balay 	  y_ptr+=len;
392827bd09bSSatish Balay 	}
393827bd09bSSatish Balay 
394827bd09bSSatish Balay       /* uu = X^T.u_l (comm portion) */
395827bd09bSSatish Balay       ssgl_radd  (uu, w, dim, stages);
396827bd09bSSatish Balay 
397827bd09bSSatish Balay       /* z = X.uu */
398827bd09bSSatish Balay       rvec_zero(z,n);
399827bd09bSSatish Balay       x_ptr=x;
400827bd09bSSatish Balay       iptr = xcol_indices;
401827bd09bSSatish Balay       for (k=0; k<i; k++)
402827bd09bSSatish Balay 	{
403827bd09bSSatish Balay 	  off = *iptr++;
404827bd09bSSatish Balay 	  len = *iptr++;
4050805154bSBarry Smith           dlen = PetscBLASIntCast(len);
4064a0de3f6SBarry Smith 	  BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
407827bd09bSSatish Balay 	  x_ptr+=len;
408827bd09bSSatish Balay 	}
409827bd09bSSatish Balay 
410827bd09bSSatish Balay       /* compute v_l = v_l - z */
411827bd09bSSatish Balay       rvec_zero(v+a_n,a_m-a_n);
4120805154bSBarry Smith       dlen = PetscBLASIntCast(n);
4134a0de3f6SBarry Smith       BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);
414827bd09bSSatish Balay 
415827bd09bSSatish Balay       /* compute u_l = A.v_l */
416827bd09bSSatish Balay       if (a_n!=a_m)
417827bd09bSSatish Balay 	{gs_gop_hc(gs_handle,v,"+\0",dim);}
418827bd09bSSatish Balay       rvec_zero(u,n);
419827bd09bSSatish Balay      do_matvec(xyt_handle->mvi,v,u);
420827bd09bSSatish Balay 
421827bd09bSSatish Balay       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
4220805154bSBarry Smith      dlen = PetscBLASIntCast(n);
4234a0de3f6SBarry Smith       alpha = BLASdot_(&dlen,u,&i1,u,&i1);
424827bd09bSSatish Balay       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
425827bd09bSSatish Balay       grop_hc(&alpha, &alpha_w, 1, op, dim);
426827bd09bSSatish Balay 
427a501084fSBarry Smith       alpha = (PetscScalar) sqrt((double)alpha);
428827bd09bSSatish Balay 
429827bd09bSSatish Balay       /* check for small alpha                             */
430827bd09bSSatish Balay       /* LATER use this to detect and determine null space */
431e32f2f54SBarry Smith       if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);
432827bd09bSSatish Balay 
433827bd09bSSatish Balay       /* compute v_l = v_l/sqrt(alpha) */
434827bd09bSSatish Balay       rvec_scale(v,1.0/alpha,n);
435827bd09bSSatish Balay       rvec_scale(u,1.0/alpha,n);
436827bd09bSSatish Balay 
437827bd09bSSatish Balay       /* add newly generated column, v_l, to X */
438827bd09bSSatish Balay       flag = 1;
439827bd09bSSatish Balay       off=len=0;
440827bd09bSSatish Balay       for (k=0; k<n; k++)
441827bd09bSSatish Balay 	{
442827bd09bSSatish Balay 	  if (v[k]!=0.0)
443827bd09bSSatish Balay 	    {
444827bd09bSSatish Balay 	      len=k;
445827bd09bSSatish Balay 	      if (flag)
446827bd09bSSatish Balay 		{off=k; flag=0;}
447827bd09bSSatish Balay 	    }
448827bd09bSSatish Balay 	}
449827bd09bSSatish Balay 
450827bd09bSSatish Balay       len -= (off-1);
451827bd09bSSatish Balay 
452827bd09bSSatish Balay       if (len>0)
453827bd09bSSatish Balay 	{
454827bd09bSSatish Balay 	  if ((xt_nnz+len)>xt_max_nnz)
455827bd09bSSatish Balay 	    {
456f1ed62a8SBarry Smith 	      ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
457827bd09bSSatish Balay 	      xt_max_nnz *= 2;
458a501084fSBarry Smith 	      x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
459827bd09bSSatish Balay 	      rvec_copy(x_ptr,x,xt_nnz);
460a501084fSBarry Smith 	      free(x);
461827bd09bSSatish Balay 	      x = x_ptr;
462827bd09bSSatish Balay 	      x_ptr+=xt_nnz;
463827bd09bSSatish Balay 	    }
464827bd09bSSatish Balay 	  xt_nnz += len;
465827bd09bSSatish Balay 	  rvec_copy(x_ptr,v+off,len);
466827bd09bSSatish Balay 
467827bd09bSSatish Balay           /* keep track of number of zeros */
468827bd09bSSatish Balay 	  if (dim)
469827bd09bSSatish Balay 	    {
470827bd09bSSatish Balay 	      for (k=0; k<len; k++)
471827bd09bSSatish Balay 		{
472827bd09bSSatish Balay 		  if (x_ptr[k]==0.0)
473827bd09bSSatish Balay 		    {xt_zero_nnz++;}
474827bd09bSSatish Balay 		}
475827bd09bSSatish Balay 	    }
476827bd09bSSatish Balay 	  else
477827bd09bSSatish Balay 	    {
478827bd09bSSatish Balay 	      for (k=0; k<len; k++)
479827bd09bSSatish Balay 		{
480827bd09bSSatish Balay 		  if (x_ptr[k]==0.0)
481827bd09bSSatish Balay 		    {xt_zero_nnz_0++;}
482827bd09bSSatish Balay 		}
483827bd09bSSatish Balay 	    }
484827bd09bSSatish Balay 	  xcol_indices[2*i] = off;
485827bd09bSSatish Balay 	  xcol_sz[i] = xcol_indices[2*i+1] = len;
486827bd09bSSatish Balay 	  xcol_vals[i] = x_ptr;
487827bd09bSSatish Balay 	}
488827bd09bSSatish Balay       else
489827bd09bSSatish Balay 	{
490827bd09bSSatish Balay 	  xcol_indices[2*i] = 0;
491827bd09bSSatish Balay 	  xcol_sz[i] = xcol_indices[2*i+1] = 0;
492827bd09bSSatish Balay 	  xcol_vals[i] = x_ptr;
493827bd09bSSatish Balay 	}
494827bd09bSSatish Balay 
495827bd09bSSatish Balay 
496827bd09bSSatish Balay       /* add newly generated column, u_l, to Y */
497827bd09bSSatish Balay       flag = 1;
498827bd09bSSatish Balay       off=len=0;
499827bd09bSSatish Balay       for (k=0; k<n; k++)
500827bd09bSSatish Balay 	{
501827bd09bSSatish Balay 	  if (u[k]!=0.0)
502827bd09bSSatish Balay 	    {
503827bd09bSSatish Balay 	      len=k;
504827bd09bSSatish Balay 	      if (flag)
505827bd09bSSatish Balay 		{off=k; flag=0;}
506827bd09bSSatish Balay 	    }
507827bd09bSSatish Balay 	}
508827bd09bSSatish Balay 
509827bd09bSSatish Balay       len -= (off-1);
510827bd09bSSatish Balay 
511827bd09bSSatish Balay       if (len>0)
512827bd09bSSatish Balay 	{
513827bd09bSSatish Balay 	  if ((yt_nnz+len)>yt_max_nnz)
514827bd09bSSatish Balay 	    {
515f1ed62a8SBarry Smith 	      ierr = PetscInfo(0,"increasing space for Y by 2x!\n");CHKERRQ(ierr);
516827bd09bSSatish Balay 	      yt_max_nnz *= 2;
517a501084fSBarry Smith 	      y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
518827bd09bSSatish Balay 	      rvec_copy(y_ptr,y,yt_nnz);
519a501084fSBarry Smith 	      free(y);
520827bd09bSSatish Balay 	      y = y_ptr;
521827bd09bSSatish Balay 	      y_ptr+=yt_nnz;
522827bd09bSSatish Balay 	    }
523827bd09bSSatish Balay 	  yt_nnz += len;
524827bd09bSSatish Balay 	  rvec_copy(y_ptr,u+off,len);
525827bd09bSSatish Balay 
526827bd09bSSatish Balay           /* keep track of number of zeros */
527827bd09bSSatish Balay 	  if (dim)
528827bd09bSSatish Balay 	    {
529827bd09bSSatish Balay 	      for (k=0; k<len; k++)
530827bd09bSSatish Balay 		{
531827bd09bSSatish Balay 		  if (y_ptr[k]==0.0)
532827bd09bSSatish Balay 		    {yt_zero_nnz++;}
533827bd09bSSatish Balay 		}
534827bd09bSSatish Balay 	    }
535827bd09bSSatish Balay 	  else
536827bd09bSSatish Balay 	    {
537827bd09bSSatish Balay 	      for (k=0; k<len; k++)
538827bd09bSSatish Balay 		{
539827bd09bSSatish Balay 		  if (y_ptr[k]==0.0)
540827bd09bSSatish Balay 		    {yt_zero_nnz_0++;}
541827bd09bSSatish Balay 		}
542827bd09bSSatish Balay 	    }
543827bd09bSSatish Balay 	  ycol_indices[2*i] = off;
544827bd09bSSatish Balay 	  ycol_sz[i] = ycol_indices[2*i+1] = len;
545827bd09bSSatish Balay 	  ycol_vals[i] = y_ptr;
546827bd09bSSatish Balay 	}
547827bd09bSSatish Balay       else
548827bd09bSSatish Balay 	{
549827bd09bSSatish Balay 	  ycol_indices[2*i] = 0;
550827bd09bSSatish Balay 	  ycol_sz[i] = ycol_indices[2*i+1] = 0;
551827bd09bSSatish Balay 	  ycol_vals[i] = y_ptr;
552827bd09bSSatish Balay 	}
553827bd09bSSatish Balay     }
554827bd09bSSatish Balay 
555827bd09bSSatish Balay   /* close off stages for execution phase */
556827bd09bSSatish Balay   while (dim!=level)
557827bd09bSSatish Balay     {
558827bd09bSSatish Balay       stages[dim++]=i;
5594a0de3f6SBarry Smith       ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
560827bd09bSSatish Balay     }
561827bd09bSSatish Balay   stages[dim]=i;
562827bd09bSSatish Balay 
563827bd09bSSatish Balay   xyt_handle->info->n=xyt_handle->mvi->n;
564827bd09bSSatish Balay   xyt_handle->info->m=m;
565827bd09bSSatish Balay   xyt_handle->info->nnz=xt_nnz + yt_nnz;
566827bd09bSSatish Balay   xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz;
567827bd09bSSatish Balay   xyt_handle->info->msg_buf_sz=stages[level]-stages[0];
568a501084fSBarry Smith   xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
569a501084fSBarry Smith   xyt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
570827bd09bSSatish Balay   xyt_handle->info->x=x;
571827bd09bSSatish Balay   xyt_handle->info->xcol_vals=xcol_vals;
572827bd09bSSatish Balay   xyt_handle->info->xcol_sz=xcol_sz;
573827bd09bSSatish Balay   xyt_handle->info->xcol_indices=xcol_indices;
574827bd09bSSatish Balay   xyt_handle->info->stages=stages;
575827bd09bSSatish Balay   xyt_handle->info->y=y;
576827bd09bSSatish Balay   xyt_handle->info->ycol_vals=ycol_vals;
577827bd09bSSatish Balay   xyt_handle->info->ycol_sz=ycol_sz;
578827bd09bSSatish Balay   xyt_handle->info->ycol_indices=ycol_indices;
579827bd09bSSatish Balay 
580a501084fSBarry Smith   free(segs);
581a501084fSBarry Smith   free(u);
582a501084fSBarry Smith   free(v);
583a501084fSBarry Smith   free(uu);
584a501084fSBarry Smith   free(z);
585a501084fSBarry Smith   free(w);
586827bd09bSSatish Balay 
587827bd09bSSatish Balay   return(0);
588827bd09bSSatish Balay }
589827bd09bSSatish Balay 
5907b1ae94cSBarry Smith /**************************************xyt.c***********************************/
5910924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle,  PetscScalar *uc)
592827bd09bSSatish Balay {
59352f87cdaSBarry Smith   PetscInt off, len, *iptr;
59452f87cdaSBarry Smith   PetscInt level       =xyt_handle->level;
59552f87cdaSBarry Smith   PetscInt n           =xyt_handle->info->n;
59652f87cdaSBarry Smith   PetscInt m           =xyt_handle->info->m;
59752f87cdaSBarry Smith   PetscInt *stages     =xyt_handle->info->stages;
59852f87cdaSBarry Smith   PetscInt *xcol_indices=xyt_handle->info->xcol_indices;
59952f87cdaSBarry Smith   PetscInt *ycol_indices=xyt_handle->info->ycol_indices;
600a501084fSBarry Smith    PetscScalar *x_ptr, *y_ptr, *uu_ptr;
601a501084fSBarry Smith   PetscScalar *solve_uu=xyt_handle->info->solve_uu;
602a501084fSBarry Smith   PetscScalar *solve_w =xyt_handle->info->solve_w;
603a501084fSBarry Smith   PetscScalar *x       =xyt_handle->info->x;
604a501084fSBarry Smith   PetscScalar *y       =xyt_handle->info->y;
6054a0de3f6SBarry Smith   PetscBLASInt i1 = 1,dlen;
606827bd09bSSatish Balay 
6070924e98cSBarry Smith   PetscFunctionBegin;
608827bd09bSSatish Balay   uu_ptr=solve_uu;
609827bd09bSSatish Balay   rvec_zero(uu_ptr,m);
610827bd09bSSatish Balay 
611827bd09bSSatish Balay   /* x  = X.Y^T.b */
612827bd09bSSatish Balay   /* uu = Y^T.b */
613827bd09bSSatish Balay   for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len)
614827bd09bSSatish Balay     {
615827bd09bSSatish Balay       off=*iptr++; len=*iptr++;
6160805154bSBarry Smith       dlen = PetscBLASIntCast(len);
6174a0de3f6SBarry Smith       *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1);
618827bd09bSSatish Balay     }
619827bd09bSSatish Balay 
620827bd09bSSatish Balay   /* comunication of beta */
621827bd09bSSatish Balay   uu_ptr=solve_uu;
622827bd09bSSatish Balay   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}
623827bd09bSSatish Balay 
624827bd09bSSatish Balay   rvec_zero(uc,n);
625827bd09bSSatish Balay 
626827bd09bSSatish Balay   /* x = X.uu */
627827bd09bSSatish Balay   for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len)
628827bd09bSSatish Balay     {
629827bd09bSSatish Balay       off=*iptr++; len=*iptr++;
6300805154bSBarry Smith       dlen = PetscBLASIntCast(len);
6314a0de3f6SBarry Smith       BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
632827bd09bSSatish Balay     }
6330924e98cSBarry Smith   PetscFunctionReturn(0);
634827bd09bSSatish Balay }
635827bd09bSSatish Balay 
6367b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6370924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle)
638827bd09bSSatish Balay {
63952f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
640827bd09bSSatish Balay 
6410924e98cSBarry Smith   PetscFunctionBegin;
642e32f2f54SBarry Smith   if (!xyt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle);
643827bd09bSSatish Balay 
644827bd09bSSatish Balay   vals[0]=vals[1]=xyt_handle->id;
645827bd09bSSatish Balay   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
646e32f2f54SBarry 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);
6470924e98cSBarry Smith   PetscFunctionReturn(0);
648827bd09bSSatish Balay }
649827bd09bSSatish Balay 
6507b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6510924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle)
652827bd09bSSatish Balay {
65352f87cdaSBarry Smith   PetscInt i, ct, id;
65452f87cdaSBarry Smith   PetscInt mask, edge, *iptr;
65552f87cdaSBarry Smith   PetscInt *dir, *used;
65652f87cdaSBarry Smith   PetscInt sum[4], w[4];
657a501084fSBarry Smith   PetscScalar rsum[4], rw[4];
65852f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
659a501084fSBarry Smith   PetscScalar *lhs, *rhs;
66052f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo=0;
661827bd09bSSatish Balay   gs_ADT gs_handle=xyt_handle->mvi->gs_handle;
66252f87cdaSBarry Smith   PetscInt *local2global=xyt_handle->mvi->local2global;
66352f87cdaSBarry Smith   PetscInt  n=xyt_handle->mvi->n;
66452f87cdaSBarry Smith   PetscInt  m=xyt_handle->mvi->m;
66552f87cdaSBarry Smith   PetscInt level=xyt_handle->level;
66652f87cdaSBarry Smith   PetscInt shared=FALSE;
667d1528f56SBarry Smith   PetscErrorCode ierr;
668827bd09bSSatish Balay 
6690924e98cSBarry Smith   PetscFunctionBegin;
67052f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
67152f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
67252f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
67352f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
67452f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
675827bd09bSSatish Balay 
676827bd09bSSatish Balay   ivec_zero(dir  ,level+1);
677827bd09bSSatish Balay   ivec_zero(nsep ,level+1);
678827bd09bSSatish Balay   ivec_zero(lnsep,level+1);
679827bd09bSSatish Balay   ivec_set (fo   ,-1,n+1);
680827bd09bSSatish Balay   ivec_zero(used,n);
681827bd09bSSatish Balay 
6828cda6cd7SBarry Smith   lhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
6838cda6cd7SBarry Smith   rhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
684827bd09bSSatish Balay 
685827bd09bSSatish Balay   /* determine the # of unique dof */
686827bd09bSSatish Balay   rvec_zero(lhs,m);
687827bd09bSSatish Balay   rvec_set(lhs,1.0,n);
688827bd09bSSatish Balay   gs_gop_hc(gs_handle,lhs,"+\0",level);
689d1528f56SBarry Smith   ierr = PetscInfo(0,"done first gs_gop_hc\n");CHKERRQ(ierr);
690827bd09bSSatish Balay   rvec_zero(rsum,2);
691827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
692827bd09bSSatish Balay     {
693827bd09bSSatish Balay       if (lhs[i]!=0.0)
694827bd09bSSatish Balay 	{rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
695827bd09bSSatish Balay 
696827bd09bSSatish Balay       if (lhs[i]!=1.0)
697827bd09bSSatish Balay 	{
698827bd09bSSatish Balay           shared=TRUE;
699827bd09bSSatish Balay         }
700827bd09bSSatish Balay     }
701827bd09bSSatish Balay 
702827bd09bSSatish Balay   grop_hc(rsum,rw,2,op,level);
703827bd09bSSatish Balay   rsum[0]+=0.1;
704827bd09bSSatish Balay   rsum[1]+=0.1;
705827bd09bSSatish Balay 
70652f87cdaSBarry Smith   xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0];
70752f87cdaSBarry Smith   xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0];
708827bd09bSSatish Balay 
709827bd09bSSatish Balay   /* determine separator sets top down */
710827bd09bSSatish Balay   if (shared)
711827bd09bSSatish Balay     {
712827bd09bSSatish Balay       /* solution is to do as in the symmetric shared case but then */
713827bd09bSSatish Balay       /* pick the sub-hc with the most free dofs and do a mat-vec   */
714827bd09bSSatish Balay       /* and pick up the responses on the other sub-hc from the     */
715827bd09bSSatish Balay       /* initial separator set obtained from the symm. shared case  */
716e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n");
717827bd09bSSatish Balay       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
718827bd09bSSatish Balay 	{
719827bd09bSSatish Balay 	  /* set rsh of hc, fire, and collect lhs responses */
720827bd09bSSatish Balay 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
721827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
722827bd09bSSatish Balay 
723827bd09bSSatish Balay 	  /* set lsh of hc, fire, and collect rhs responses */
724827bd09bSSatish Balay 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
725827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
726827bd09bSSatish Balay 
727827bd09bSSatish Balay 	  for (i=0;i<n;i++)
728827bd09bSSatish Balay 	    {
729827bd09bSSatish Balay 	      if (id< mask)
730827bd09bSSatish Balay 		{
731827bd09bSSatish Balay 		  if (lhs[i]!=0.0)
732827bd09bSSatish Balay 		    {lhs[i]=1.0;}
733827bd09bSSatish Balay 		}
734827bd09bSSatish Balay 	      if (id>=mask)
735827bd09bSSatish Balay 		{
736827bd09bSSatish Balay 		  if (rhs[i]!=0.0)
737827bd09bSSatish Balay 		    {rhs[i]=1.0;}
738827bd09bSSatish Balay 		}
739827bd09bSSatish Balay 	    }
740827bd09bSSatish Balay 
741827bd09bSSatish Balay 	  if (id< mask)
742827bd09bSSatish Balay 	    {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
743827bd09bSSatish Balay 	  else
744827bd09bSSatish Balay 	    {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}
745827bd09bSSatish Balay 
746827bd09bSSatish Balay 	  /* count number of dofs I own that have signal and not in sep set */
747827bd09bSSatish Balay 	  rvec_zero(rsum,4);
748827bd09bSSatish Balay 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
749827bd09bSSatish Balay 	    {
750827bd09bSSatish Balay 	      if (!used[i])
751827bd09bSSatish Balay 		{
752827bd09bSSatish Balay 		  /* number of unmarked dofs on node */
753827bd09bSSatish Balay 		  ct++;
754827bd09bSSatish Balay 		  /* number of dofs to be marked on lhs hc */
755827bd09bSSatish Balay 		  if (id< mask)
756827bd09bSSatish Balay 		    {
757827bd09bSSatish Balay 		      if (lhs[i]!=0.0)
758827bd09bSSatish Balay 			{sum[0]++; rsum[0]+=1.0/lhs[i];}
759827bd09bSSatish Balay 		    }
760827bd09bSSatish Balay 		  /* number of dofs to be marked on rhs hc */
761827bd09bSSatish Balay 		  if (id>=mask)
762827bd09bSSatish Balay 		    {
763827bd09bSSatish Balay 		      if (rhs[i]!=0.0)
764827bd09bSSatish Balay 			{sum[1]++; rsum[1]+=1.0/rhs[i];}
765827bd09bSSatish Balay 		    }
766827bd09bSSatish Balay 		}
767827bd09bSSatish Balay 	    }
768827bd09bSSatish Balay 
769827bd09bSSatish Balay 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
770827bd09bSSatish Balay 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
771827bd09bSSatish Balay 	  (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
772827bd09bSSatish Balay 	  giop_hc(sum,w,4,op,edge);
773827bd09bSSatish Balay 	  grop_hc(rsum,rw,4,op,edge);
774827bd09bSSatish Balay 	  rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
775827bd09bSSatish Balay 
776827bd09bSSatish Balay 	  if (id<mask)
777827bd09bSSatish Balay 	    {
778827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
779827bd09bSSatish Balay 	      for (ct=i=0;i<n;i++)
780827bd09bSSatish Balay 		{
781827bd09bSSatish Balay 		  if ((!used[i])&&(lhs[i]!=0.0))
782827bd09bSSatish Balay 		    {
783827bd09bSSatish Balay 		      ct++; nfo++;
784827bd09bSSatish Balay 
785c1235816SBarry Smith 		      if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
786827bd09bSSatish Balay 
787827bd09bSSatish Balay 		      *--iptr = local2global[i];
788827bd09bSSatish Balay 		      used[i]=edge;
789827bd09bSSatish Balay 		    }
790827bd09bSSatish Balay 		}
791827bd09bSSatish Balay 	      if (ct>1) {ivec_sort(iptr,ct);}
792827bd09bSSatish Balay 
793827bd09bSSatish Balay 	      lnsep[edge]=ct;
79452f87cdaSBarry Smith 	      nsep[edge]=(PetscInt) rsum[0];
795827bd09bSSatish Balay 	      dir [edge]=LEFT;
796827bd09bSSatish Balay 	    }
797827bd09bSSatish Balay 
798827bd09bSSatish Balay 	  if (id>=mask)
799827bd09bSSatish Balay 	    {
800827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
801827bd09bSSatish Balay 	      for (ct=i=0;i<n;i++)
802827bd09bSSatish Balay 		{
803827bd09bSSatish Balay 		  if ((!used[i])&&(rhs[i]!=0.0))
804827bd09bSSatish Balay 		    {
805827bd09bSSatish Balay 		      ct++; nfo++;
806827bd09bSSatish Balay 
807c1235816SBarry Smith 		      if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
808827bd09bSSatish Balay 
809827bd09bSSatish Balay 		      *--iptr = local2global[i];
810827bd09bSSatish Balay 		      used[i]=edge;
811827bd09bSSatish Balay 		    }
812827bd09bSSatish Balay 		}
813827bd09bSSatish Balay 	      if (ct>1) {ivec_sort(iptr,ct);}
814827bd09bSSatish Balay 
815827bd09bSSatish Balay 	      lnsep[edge]=ct;
81652f87cdaSBarry Smith 	      nsep[edge]= (PetscInt) rsum[1];
817827bd09bSSatish Balay 	      dir [edge]=RIGHT;
818827bd09bSSatish Balay 	    }
819827bd09bSSatish Balay 
820827bd09bSSatish Balay 	  /* LATER or we can recur on these to order seps at this level */
821827bd09bSSatish Balay 	  /* do we need full set of separators for this?                */
822827bd09bSSatish Balay 
823827bd09bSSatish Balay 	  /* fold rhs hc into lower */
824827bd09bSSatish Balay 	  if (id>=mask)
825827bd09bSSatish Balay 	    {id-=mask;}
826827bd09bSSatish Balay 	}
827827bd09bSSatish Balay     }
828827bd09bSSatish Balay   else
829827bd09bSSatish Balay     {
830827bd09bSSatish Balay       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
831827bd09bSSatish Balay 	{
832827bd09bSSatish Balay 	  /* set rsh of hc, fire, and collect lhs responses */
833827bd09bSSatish Balay 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
834827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
835827bd09bSSatish Balay 
836827bd09bSSatish Balay 	  /* set lsh of hc, fire, and collect rhs responses */
837827bd09bSSatish Balay 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
838827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
839827bd09bSSatish Balay 
840827bd09bSSatish Balay 	  /* count number of dofs I own that have signal and not in sep set */
841827bd09bSSatish Balay 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
842827bd09bSSatish Balay 	    {
843827bd09bSSatish Balay 	      if (!used[i])
844827bd09bSSatish Balay 		{
845827bd09bSSatish Balay 		  /* number of unmarked dofs on node */
846827bd09bSSatish Balay 		  ct++;
847827bd09bSSatish Balay 		  /* number of dofs to be marked on lhs hc */
848827bd09bSSatish Balay 		  if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
849827bd09bSSatish Balay 		  /* number of dofs to be marked on rhs hc */
850827bd09bSSatish Balay 		  if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
851827bd09bSSatish Balay 		}
852827bd09bSSatish Balay 	    }
853827bd09bSSatish Balay 
854827bd09bSSatish Balay 	  /* for the non-symmetric case we need separators of width 2 */
855827bd09bSSatish Balay 	  /* so take both sides */
856827bd09bSSatish Balay 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
857827bd09bSSatish Balay 	  giop_hc(sum,w,4,op,edge);
858827bd09bSSatish Balay 
859827bd09bSSatish Balay 	  ct=0;
860827bd09bSSatish Balay 	  if (id<mask)
861827bd09bSSatish Balay 	    {
862827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
863827bd09bSSatish Balay 	      for (i=0;i<n;i++)
864827bd09bSSatish Balay 		{
865827bd09bSSatish Balay 		  if ((!used[i])&&(lhs[i]!=0.0))
866827bd09bSSatish Balay 		    {
867827bd09bSSatish Balay 		      ct++; nfo++;
868827bd09bSSatish Balay 		      *--iptr = local2global[i];
869827bd09bSSatish Balay 		      used[i]=edge;
870827bd09bSSatish Balay 		    }
871827bd09bSSatish Balay 		}
872827bd09bSSatish Balay 	      /* LSH hc summation of ct should be sum[0] */
873827bd09bSSatish Balay 	    }
874827bd09bSSatish Balay 	  else
875827bd09bSSatish Balay 	    {
876827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
877827bd09bSSatish Balay 	      for (i=0;i<n;i++)
878827bd09bSSatish Balay 		{
879827bd09bSSatish Balay 		  if ((!used[i])&&(rhs[i]!=0.0))
880827bd09bSSatish Balay 		    {
881827bd09bSSatish Balay 		      ct++; nfo++;
882827bd09bSSatish Balay 		      *--iptr = local2global[i];
883827bd09bSSatish Balay 		      used[i]=edge;
884827bd09bSSatish Balay 		    }
885827bd09bSSatish Balay 		}
886827bd09bSSatish Balay 	      /* RSH hc summation of ct should be sum[1] */
887827bd09bSSatish Balay 	    }
888827bd09bSSatish Balay 
889827bd09bSSatish Balay 	  if (ct>1) {ivec_sort(iptr,ct);}
890827bd09bSSatish Balay 	  lnsep[edge]=ct;
891827bd09bSSatish Balay 	  nsep[edge]=sum[0]+sum[1];
892827bd09bSSatish Balay 	  dir [edge]=BOTH;
893827bd09bSSatish Balay 
894827bd09bSSatish Balay 	  /* LATER or we can recur on these to order seps at this level */
895827bd09bSSatish Balay 	  /* do we need full set of separators for this?                */
896827bd09bSSatish Balay 
897827bd09bSSatish Balay 	  /* fold rhs hc into lower */
898827bd09bSSatish Balay 	  if (id>=mask)
899827bd09bSSatish Balay 	    {id-=mask;}
900827bd09bSSatish Balay 	}
901827bd09bSSatish Balay     }
902827bd09bSSatish Balay 
903827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
904827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
905827bd09bSSatish Balay     {
906827bd09bSSatish Balay       if (!used[i])
907827bd09bSSatish Balay 	{
908827bd09bSSatish Balay 	  ct++; nfo++;
909827bd09bSSatish Balay 	  *--iptr = local2global[i];
910827bd09bSSatish Balay 	  used[i]=edge;
911827bd09bSSatish Balay 	}
912827bd09bSSatish Balay     }
913827bd09bSSatish Balay   if (ct>1) {ivec_sort(iptr,ct);}
914827bd09bSSatish Balay   lnsep[edge]=ct;
915827bd09bSSatish Balay   nsep [edge]=ct;
916827bd09bSSatish Balay   dir  [edge]=BOTH;
917827bd09bSSatish Balay 
918827bd09bSSatish Balay   xyt_handle->info->nsep=nsep;
919827bd09bSSatish Balay   xyt_handle->info->lnsep=lnsep;
920827bd09bSSatish Balay   xyt_handle->info->fo=fo;
921827bd09bSSatish Balay   xyt_handle->info->nfo=nfo;
922827bd09bSSatish Balay 
923a501084fSBarry Smith   free(dir);
924a501084fSBarry Smith   free(lhs);
925a501084fSBarry Smith   free(rhs);
926a501084fSBarry Smith   free(used);
9270924e98cSBarry Smith   PetscFunctionReturn(0);
928827bd09bSSatish Balay }
929827bd09bSSatish Balay 
9307b1ae94cSBarry Smith /**************************************xyt.c***********************************/
93152f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
932827bd09bSSatish Balay {
933827bd09bSSatish Balay   mv_info *mvi;
934827bd09bSSatish Balay 
935827bd09bSSatish Balay 
936a501084fSBarry Smith   mvi = (mv_info*)malloc(sizeof(mv_info));
937827bd09bSSatish Balay   mvi->n=n;
938827bd09bSSatish Balay   mvi->m=m;
939827bd09bSSatish Balay   mvi->n_global=-1;
940827bd09bSSatish Balay   mvi->m_global=-1;
94152f87cdaSBarry Smith   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
942827bd09bSSatish Balay   ivec_copy(mvi->local2global,local2global,m);
943827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
944a501084fSBarry Smith   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
945827bd09bSSatish Balay   mvi->grid_data=grid_data;
946827bd09bSSatish Balay 
947827bd09bSSatish Balay   /* set xyt communication handle to perform restricted matvec */
948827bd09bSSatish Balay   mvi->gs_handle = gs_init(local2global, m, num_nodes);
949827bd09bSSatish Balay 
950827bd09bSSatish Balay   return(mvi);
951827bd09bSSatish Balay }
952827bd09bSSatish Balay 
9537b1ae94cSBarry Smith /**************************************xyt.c***********************************/
9543fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
955827bd09bSSatish Balay {
9560924e98cSBarry Smith   PetscFunctionBegin;
957827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
9580924e98cSBarry Smith   PetscFunctionReturn(0);
959827bd09bSSatish Balay }
960827bd09bSSatish Balay 
961827bd09bSSatish Balay 
962827bd09bSSatish Balay 
963