xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision 8cda6cd7047ed8bb2d68b8e371c012553cb590ce)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2827bd09bSSatish Balay 
3827bd09bSSatish Balay /*************************************xxt.c************************************
4827bd09bSSatish Balay Module Name: xxt
5827bd09bSSatish Balay Module Info:
6827bd09bSSatish Balay 
7827bd09bSSatish Balay author:  Henry M. Tufo III
8827bd09bSSatish Balay e-mail:  hmt@asci.uchicago.edu
9827bd09bSSatish Balay contact:
10827bd09bSSatish Balay +--------------------------------+--------------------------------+
11827bd09bSSatish Balay |MCS Division - Building 221     |Department of Computer Science  |
12827bd09bSSatish Balay |Argonne National Laboratory     |Ryerson 152                     |
13827bd09bSSatish Balay |9700 S. Cass Avenue             |The University of Chicago       |
14827bd09bSSatish Balay |Argonne, IL  60439              |Chicago, IL  60637              |
15827bd09bSSatish Balay |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
16827bd09bSSatish Balay +--------------------------------+--------------------------------+
17827bd09bSSatish Balay 
18827bd09bSSatish Balay Last Modification: 3.20.01
19827bd09bSSatish Balay **************************************xxt.c***********************************/
207758a8cdSBarry Smith #include "src/ksp/pc/impls/tfs/tfs.h"
21827bd09bSSatish Balay 
22827bd09bSSatish Balay #define LEFT  -1
23827bd09bSSatish Balay #define RIGHT  1
24827bd09bSSatish Balay #define BOTH   0
25827bd09bSSatish Balay 
26827bd09bSSatish Balay typedef struct xxt_solver_info {
2752f87cdaSBarry Smith   PetscInt n, m, n_global, m_global;
2852f87cdaSBarry Smith   PetscInt nnz, max_nnz, msg_buf_sz;
2952f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo, *stages;
3052f87cdaSBarry Smith   PetscInt *col_sz, *col_indices;
31a501084fSBarry Smith   PetscScalar **col_vals, *x, *solve_uu, *solve_w;
3252f87cdaSBarry Smith   PetscInt nsolves;
33a501084fSBarry Smith   PetscScalar tot_solve_time;
34827bd09bSSatish Balay } xxt_info;
35827bd09bSSatish Balay 
36827bd09bSSatish Balay typedef struct matvec_info {
3752f87cdaSBarry Smith   PetscInt n, m, n_global, m_global;
3852f87cdaSBarry Smith   PetscInt *local2global;
39827bd09bSSatish Balay   gs_ADT gs_handle;
40a501084fSBarry Smith   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
41827bd09bSSatish Balay   void *grid_data;
42827bd09bSSatish Balay } mv_info;
43827bd09bSSatish Balay 
44827bd09bSSatish Balay struct xxt_CDT{
4552f87cdaSBarry Smith   PetscInt id;
4652f87cdaSBarry Smith   PetscInt ns;
4752f87cdaSBarry Smith   PetscInt level;
48827bd09bSSatish Balay   xxt_info *info;
49827bd09bSSatish Balay   mv_info  *mvi;
50827bd09bSSatish Balay };
51827bd09bSSatish Balay 
5252f87cdaSBarry Smith static PetscInt n_xxt=0;
5352f87cdaSBarry Smith static PetscInt n_xxt_handles=0;
54827bd09bSSatish Balay 
55827bd09bSSatish Balay /* prototypes */
563fdc5746SBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *rhs);
573fdc5746SBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle);
583fdc5746SBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle);
593fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
6052f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle);
6152f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle);
6252f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data);
63a501084fSBarry Smith 
647b1ae94cSBarry Smith /**************************************xxt.c***********************************/
657b1ae94cSBarry Smith xxt_ADT XXT_new(void)
66827bd09bSSatish Balay {
67827bd09bSSatish Balay   xxt_ADT xxt_handle;
68827bd09bSSatish Balay 
69827bd09bSSatish Balay   /* rolling count on n_xxt ... pot. problem here */
70827bd09bSSatish Balay   n_xxt_handles++;
71a501084fSBarry Smith   xxt_handle       = (xxt_ADT)malloc(sizeof(struct xxt_CDT));
72827bd09bSSatish Balay   xxt_handle->id   = ++n_xxt;
73827bd09bSSatish Balay   xxt_handle->info = NULL; xxt_handle->mvi  = NULL;
74827bd09bSSatish Balay 
75827bd09bSSatish Balay   return(xxt_handle);
76827bd09bSSatish Balay }
77827bd09bSSatish Balay 
787b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7952f87cdaSBarry Smith PetscInt XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt  handle */
8052f87cdaSBarry Smith 	   PetscInt *local2global,  /* global column mapping       */
8152f87cdaSBarry Smith 	   PetscInt n,              /* local num rows              */
8252f87cdaSBarry Smith 	   PetscInt m,              /* local num cols              */
83827bd09bSSatish Balay 	   void *matvec,       /* b_loc=A_local.x_loc         */
84827bd09bSSatish Balay 	   void *grid_data     /* grid data for matvec        */
85827bd09bSSatish Balay 	   )
86827bd09bSSatish Balay {
87a40341efSBarry Smith   comm_init();
88827bd09bSSatish Balay   check_handle(xxt_handle);
89827bd09bSSatish Balay 
90827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
91827bd09bSSatish Balay   if ((1<<(xxt_handle->level=i_log2_num_nodes))!=num_nodes)
9271a0148aSBarry Smith     {SETERRQ2(PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %D != %D\n",1<<i_log2_num_nodes,num_nodes);}
93827bd09bSSatish Balay 
94827bd09bSSatish Balay   /* space for X info */
95a501084fSBarry Smith   xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info));
96827bd09bSSatish Balay 
97827bd09bSSatish Balay   /* set up matvec handles */
98827bd09bSSatish Balay   xxt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);
99827bd09bSSatish Balay 
100827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
101827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
102827bd09bSSatish Balay   xxt_handle->ns=0;
103827bd09bSSatish Balay 
104827bd09bSSatish Balay   /* determine separators and generate firing order - NB xxt info set here */
105827bd09bSSatish Balay   det_separators(xxt_handle);
106827bd09bSSatish Balay 
107827bd09bSSatish Balay   return(do_xxt_factor(xxt_handle));
108827bd09bSSatish Balay }
109827bd09bSSatish Balay 
1107b1ae94cSBarry Smith /**************************************xxt.c***********************************/
111*8cda6cd7SBarry Smith PetscInt XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b)
112827bd09bSSatish Balay {
113827bd09bSSatish Balay 
114a40341efSBarry Smith   comm_init();
115827bd09bSSatish Balay   check_handle(xxt_handle);
116827bd09bSSatish Balay 
117827bd09bSSatish Balay   /* need to copy b into x? */
118827bd09bSSatish Balay   if (b)
119827bd09bSSatish Balay     {rvec_copy(x,b,xxt_handle->mvi->n);}
120827bd09bSSatish Balay   do_xxt_solve(xxt_handle,x);
121827bd09bSSatish Balay 
122827bd09bSSatish Balay   return(0);
123827bd09bSSatish Balay }
124827bd09bSSatish Balay 
1257b1ae94cSBarry Smith /**************************************xxt.c***********************************/
12652f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle)
127827bd09bSSatish Balay {
128827bd09bSSatish Balay 
129a40341efSBarry Smith   comm_init();
130827bd09bSSatish Balay   check_handle(xxt_handle);
131827bd09bSSatish Balay   n_xxt_handles--;
132827bd09bSSatish Balay 
133a501084fSBarry Smith   free(xxt_handle->info->nsep);
134a501084fSBarry Smith   free(xxt_handle->info->lnsep);
135a501084fSBarry Smith   free(xxt_handle->info->fo);
136a501084fSBarry Smith   free(xxt_handle->info->stages);
137a501084fSBarry Smith   free(xxt_handle->info->solve_uu);
138a501084fSBarry Smith   free(xxt_handle->info->solve_w);
139a501084fSBarry Smith   free(xxt_handle->info->x);
140a501084fSBarry Smith   free(xxt_handle->info->col_vals);
141a501084fSBarry Smith   free(xxt_handle->info->col_sz);
142a501084fSBarry Smith   free(xxt_handle->info->col_indices);
143a501084fSBarry Smith   free(xxt_handle->info);
144a501084fSBarry Smith   free(xxt_handle->mvi->local2global);
145827bd09bSSatish Balay    gs_free(xxt_handle->mvi->gs_handle);
146a501084fSBarry Smith   free(xxt_handle->mvi);
147a501084fSBarry Smith   free(xxt_handle);
148827bd09bSSatish Balay 
149827bd09bSSatish Balay   /* if the check fails we nuke */
150a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
151827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
152827bd09bSSatish Balay   return(0);
153827bd09bSSatish Balay }
154827bd09bSSatish Balay 
1557b1ae94cSBarry Smith /**************************************xxt.c***********************************/
15652f87cdaSBarry Smith PetscInt XXT_stats(xxt_ADT xxt_handle)
157827bd09bSSatish Balay {
15852f87cdaSBarry Smith   PetscInt  op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
15952f87cdaSBarry Smith   PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
16052f87cdaSBarry Smith   PetscInt   vals[9],  work[9];
161a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
16271a0148aSBarry Smith   PetscErrorCode ierr;
163827bd09bSSatish Balay 
164a40341efSBarry Smith   comm_init();
165827bd09bSSatish Balay   check_handle(xxt_handle);
166827bd09bSSatish Balay 
167827bd09bSSatish Balay   /* if factorization not done there are no stats */
168827bd09bSSatish Balay   if (!xxt_handle->info||!xxt_handle->mvi)
169827bd09bSSatish Balay     {
170827bd09bSSatish Balay       if (!my_id)
17171a0148aSBarry Smith 	{ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");}
172827bd09bSSatish Balay       return 1;
173827bd09bSSatish Balay     }
174827bd09bSSatish Balay 
175827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
176827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
177827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
178827bd09bSSatish Balay   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
179827bd09bSSatish Balay 
180827bd09bSSatish Balay   fvals[0]=fvals[1]=fvals[2]
181827bd09bSSatish Balay     =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
182827bd09bSSatish Balay   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
183827bd09bSSatish Balay 
184827bd09bSSatish Balay   if (!my_id)
185827bd09bSSatish Balay     {
18671a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_nnz=%D\n",my_id,vals[0]);
18771a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_nnz=%D\n",my_id,vals[1]);
18871a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
18971a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_nnz=%D\n",my_id,vals[2]);
19071a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
19171a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
19271a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_n  =%D\n",my_id,vals[3]);
19371a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_n  =%D\n",my_id,vals[4]);
19471a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
19571a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_n  =%D\n",my_id,vals[5]);
19671a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_buf=%D\n",my_id,vals[6]);
19771a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_buf=%D\n",my_id,vals[7]);
19871a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
19971a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_slv=%g\n",my_id,fvals[0]);
20071a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_slv=%g\n",my_id,fvals[1]);
20171a0148aSBarry Smith       ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_slv=%g\n",my_id,fvals[2]/num_nodes);
202827bd09bSSatish Balay     }
203827bd09bSSatish Balay 
204827bd09bSSatish Balay   return(0);
205827bd09bSSatish Balay }
206827bd09bSSatish Balay 
207827bd09bSSatish Balay /*************************************xxt.c************************************
208827bd09bSSatish Balay 
209827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
210827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
211827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
212827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
213827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
214827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
215827bd09bSSatish Balay    gs_init/gop).
216827bd09bSSatish Balay 
217827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
218827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
219827bd09bSSatish Balay **************************************xxt.c***********************************/
22052f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle)
221827bd09bSSatish Balay {
2227b1ae94cSBarry Smith   return xxt_generate(xxt_handle);
223827bd09bSSatish Balay }
224827bd09bSSatish Balay 
2257b1ae94cSBarry Smith /**************************************xxt.c***********************************/
22652f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle)
227827bd09bSSatish Balay {
22852f87cdaSBarry Smith   PetscInt i,j,k,idex;
22952f87cdaSBarry Smith   PetscInt dim, col;
230a501084fSBarry Smith   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
23152f87cdaSBarry Smith   PetscInt *segs;
23252f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
23352f87cdaSBarry Smith   PetscInt off, len;
234a501084fSBarry Smith   PetscScalar *x_ptr;
23552f87cdaSBarry Smith   PetscInt *iptr, flag;
23652f87cdaSBarry Smith   PetscInt start=0, end, work;
23752f87cdaSBarry Smith   PetscInt op2[] = {GL_MIN,0};
238827bd09bSSatish Balay   gs_ADT gs_handle;
23952f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo;
24052f87cdaSBarry Smith   PetscInt a_n=xxt_handle->mvi->n;
24152f87cdaSBarry Smith   PetscInt a_m=xxt_handle->mvi->m;
24252f87cdaSBarry Smith   PetscInt *a_local2global=xxt_handle->mvi->local2global;
24352f87cdaSBarry Smith   PetscInt level;
24452f87cdaSBarry Smith   PetscInt xxt_nnz=0, xxt_max_nnz=0;
24552f87cdaSBarry Smith   PetscInt n, m;
24652f87cdaSBarry Smith   PetscInt *col_sz, *col_indices, *stages;
247a501084fSBarry Smith   PetscScalar **col_vals, *x;
24852f87cdaSBarry Smith   PetscInt n_global;
24952f87cdaSBarry Smith   PetscInt xxt_zero_nnz=0;
25052f87cdaSBarry Smith   PetscInt xxt_zero_nnz_0=0;
25171a0148aSBarry Smith   PetscBLASInt i1 = 1,dlen;
252a501084fSBarry Smith   PetscScalar dm1 = -1.0;
253d1528f56SBarry Smith   PetscErrorCode ierr;
254827bd09bSSatish Balay 
255827bd09bSSatish Balay   n=xxt_handle->mvi->n;
256827bd09bSSatish Balay   nsep=xxt_handle->info->nsep;
257827bd09bSSatish Balay   lnsep=xxt_handle->info->lnsep;
258827bd09bSSatish Balay   fo=xxt_handle->info->fo;
259827bd09bSSatish Balay   end=lnsep[0];
260827bd09bSSatish Balay   level=xxt_handle->level;
261827bd09bSSatish Balay   gs_handle=xxt_handle->mvi->gs_handle;
262827bd09bSSatish Balay 
263827bd09bSSatish Balay   /* is there a null space? */
264827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
265827bd09bSSatish Balay   for (i=0, j=0; i<=level; i++)
266827bd09bSSatish Balay     {j+=nsep[i];}
267827bd09bSSatish Balay 
268827bd09bSSatish Balay   m = j-xxt_handle->ns;
269827bd09bSSatish Balay   if (m!=j)
27071a0148aSBarry Smith     {ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);}
271827bd09bSSatish Balay 
272827bd09bSSatish Balay   /* get and initialize storage for x local         */
273827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
27452f87cdaSBarry Smith   col_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
27552f87cdaSBarry Smith   col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
276a501084fSBarry Smith   col_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
277827bd09bSSatish Balay   for (i=j=0; i<m; i++, j+=2)
278827bd09bSSatish Balay     {
279827bd09bSSatish Balay       col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
280827bd09bSSatish Balay       col_vals[i] = NULL;
281827bd09bSSatish Balay     }
282827bd09bSSatish Balay   col_indices[j]=-1;
283827bd09bSSatish Balay 
284827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
285827bd09bSSatish Balay   /* this looks like nsep[]=segments */
28652f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
28752f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
288827bd09bSSatish Balay   ivec_zero(stages,level+1);
289827bd09bSSatish Balay   ivec_copy(segs,nsep,level+1);
290827bd09bSSatish Balay   for (i=0; i<level; i++)
291827bd09bSSatish Balay     {segs[i+1] += segs[i];}
292827bd09bSSatish Balay   stages[0] = segs[0];
293827bd09bSSatish Balay 
294827bd09bSSatish Balay   /* temporary vectors  */
295a501084fSBarry Smith   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
296a501084fSBarry Smith   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
297a501084fSBarry Smith   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
298a501084fSBarry Smith   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
299a501084fSBarry Smith   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
300827bd09bSSatish Balay 
301827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
302827bd09bSSatish Balay   for (i=1, j=0; i<=level; i++)
303827bd09bSSatish Balay     {j+=nsep[i];}
304827bd09bSSatish Balay 
305827bd09bSSatish Balay   /* storage for sparse x values */
306827bd09bSSatish Balay   n_global = xxt_handle->info->n_global;
30752f87cdaSBarry Smith   xxt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
308a501084fSBarry Smith   x = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
309827bd09bSSatish Balay   xxt_nnz = 0;
310827bd09bSSatish Balay 
311827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
312827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
313827bd09bSSatish Balay   for (dim=i=j=0;i<m;i++)
314827bd09bSSatish Balay     {
315827bd09bSSatish Balay       /* time to move to the next level? */
316827bd09bSSatish Balay       while (i==segs[dim])
317827bd09bSSatish Balay 	{
318827bd09bSSatish Balay 	  if (dim==level)
319388eb383SBarry Smith 	    {SETERRQ(PETSC_ERR_PLIB,"dim about to exceed level\n"); break;}
320827bd09bSSatish Balay 
321827bd09bSSatish Balay 	  stages[dim++]=i;
322827bd09bSSatish Balay 	  end+=lnsep[dim];
323827bd09bSSatish Balay 	}
324827bd09bSSatish Balay       stages[dim]=i;
325827bd09bSSatish Balay 
326827bd09bSSatish Balay       /* which column are we firing? */
327827bd09bSSatish Balay       /* i.e. set v_l */
328827bd09bSSatish Balay       /* use new seps and do global min across hc to determine which one to fire */
329827bd09bSSatish Balay       (start<end) ? (col=fo[start]) : (col=INT_MAX);
330827bd09bSSatish Balay       giop_hc(&col,&work,1,op2,dim);
331827bd09bSSatish Balay 
332827bd09bSSatish Balay       /* shouldn't need this */
333827bd09bSSatish Balay       if (col==INT_MAX)
334827bd09bSSatish Balay 	{
335f1ed62a8SBarry Smith 	  ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
336827bd09bSSatish Balay 	  continue;
337827bd09bSSatish Balay 	}
338827bd09bSSatish Balay 
339827bd09bSSatish Balay       /* do I own it? I should */
340827bd09bSSatish Balay       rvec_zero(v ,a_m);
341827bd09bSSatish Balay       if (col==fo[start])
342827bd09bSSatish Balay 	{
343827bd09bSSatish Balay 	  start++;
344827bd09bSSatish Balay 	  idex=ivec_linear_search(col, a_local2global, a_n);
345827bd09bSSatish Balay 	  if (idex!=-1)
346827bd09bSSatish Balay 	    {v[idex] = 1.0; j++;}
347827bd09bSSatish Balay 	  else
348388eb383SBarry Smith 	    {SETERRQ(PETSC_ERR_PLIB,"NOT FOUND!\n");}
349827bd09bSSatish Balay 	}
350827bd09bSSatish Balay       else
351827bd09bSSatish Balay 	{
352827bd09bSSatish Balay 	  idex=ivec_linear_search(col, a_local2global, a_m);
353827bd09bSSatish Balay 	  if (idex!=-1)
354827bd09bSSatish Balay 	    {v[idex] = 1.0;}
355827bd09bSSatish Balay 	}
356827bd09bSSatish Balay 
357827bd09bSSatish Balay       /* perform u = A.v_l */
358827bd09bSSatish Balay       rvec_zero(u,n);
359827bd09bSSatish Balay       do_matvec(xxt_handle->mvi,v,u);
360827bd09bSSatish Balay 
361827bd09bSSatish Balay       /* uu =  X^T.u_l (local portion) */
362827bd09bSSatish Balay       /* technically only need to zero out first i entries */
363827bd09bSSatish Balay       /* later turn this into an XXT_solve call ? */
364827bd09bSSatish Balay       rvec_zero(uu,m);
365827bd09bSSatish Balay       x_ptr=x;
366827bd09bSSatish Balay       iptr = col_indices;
367827bd09bSSatish Balay       for (k=0; k<i; k++)
368827bd09bSSatish Balay 	{
369827bd09bSSatish Balay 	  off = *iptr++;
370827bd09bSSatish Balay 	  len = *iptr++;
3710805154bSBarry Smith           dlen = PetscBLASIntCast(len);
37271a0148aSBarry Smith 	  uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1);
373827bd09bSSatish Balay 	  x_ptr+=len;
374827bd09bSSatish Balay 	}
375827bd09bSSatish Balay 
376827bd09bSSatish Balay 
377827bd09bSSatish Balay       /* uu = X^T.u_l (comm portion) */
378827bd09bSSatish Balay       ssgl_radd  (uu, w, dim, stages);
379827bd09bSSatish Balay 
380827bd09bSSatish Balay       /* z = X.uu */
381827bd09bSSatish Balay       rvec_zero(z,n);
382827bd09bSSatish Balay       x_ptr=x;
383827bd09bSSatish Balay       iptr = col_indices;
384827bd09bSSatish Balay       for (k=0; k<i; k++)
385827bd09bSSatish Balay 	{
386827bd09bSSatish Balay 	  off = *iptr++;
387827bd09bSSatish Balay 	  len = *iptr++;
3880805154bSBarry Smith           dlen = PetscBLASIntCast(len);
38971a0148aSBarry Smith 	  BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
390827bd09bSSatish Balay 	  x_ptr+=len;
391827bd09bSSatish Balay 	}
392827bd09bSSatish Balay 
393827bd09bSSatish Balay       /* compute v_l = v_l - z */
394827bd09bSSatish Balay       rvec_zero(v+a_n,a_m-a_n);
3950805154bSBarry Smith       dlen = PetscBLASIntCast(n);
39671a0148aSBarry Smith       BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);
397827bd09bSSatish Balay 
398827bd09bSSatish Balay       /* compute u_l = A.v_l */
399827bd09bSSatish Balay       if (a_n!=a_m)
400827bd09bSSatish Balay 	{gs_gop_hc(gs_handle,v,"+\0",dim);}
401827bd09bSSatish Balay       rvec_zero(u,n);
402827bd09bSSatish Balay       do_matvec(xxt_handle->mvi,v,u);
403827bd09bSSatish Balay 
404827bd09bSSatish Balay       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
4050805154bSBarry Smith       dlen = PetscBLASIntCast(n);
40671a0148aSBarry Smith       alpha = BLASdot_(&dlen,u,&i1,v,&i1);
407827bd09bSSatish Balay       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
408827bd09bSSatish Balay       grop_hc(&alpha, &alpha_w, 1, op, dim);
409827bd09bSSatish Balay 
410a501084fSBarry Smith       alpha = (PetscScalar) sqrt((double)alpha);
411827bd09bSSatish Balay 
412827bd09bSSatish Balay       /* check for small alpha                             */
413827bd09bSSatish Balay       /* LATER use this to detect and determine null space */
414827bd09bSSatish Balay       if (fabs(alpha)<1.0e-14)
415388eb383SBarry Smith 	{SETERRQ1(PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);}
416827bd09bSSatish Balay 
417827bd09bSSatish Balay       /* compute v_l = v_l/sqrt(alpha) */
418827bd09bSSatish Balay       rvec_scale(v,1.0/alpha,n);
419827bd09bSSatish Balay 
420827bd09bSSatish Balay       /* add newly generated column, v_l, to X */
421827bd09bSSatish Balay       flag = 1;
422827bd09bSSatish Balay       off=len=0;
423827bd09bSSatish Balay       for (k=0; k<n; k++)
424827bd09bSSatish Balay 	{
425827bd09bSSatish Balay 	  if (v[k]!=0.0)
426827bd09bSSatish Balay 	    {
427827bd09bSSatish Balay 	      len=k;
428827bd09bSSatish Balay 	      if (flag)
429827bd09bSSatish Balay 		{off=k; flag=0;}
430827bd09bSSatish Balay 	    }
431827bd09bSSatish Balay 	}
432827bd09bSSatish Balay 
433827bd09bSSatish Balay       len -= (off-1);
434827bd09bSSatish Balay 
435827bd09bSSatish Balay       if (len>0)
436827bd09bSSatish Balay 	{
437827bd09bSSatish Balay 	  if ((xxt_nnz+len)>xxt_max_nnz)
438827bd09bSSatish Balay 	    {
439f1ed62a8SBarry Smith 	      ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
440827bd09bSSatish Balay 	      xxt_max_nnz *= 2;
441a501084fSBarry Smith 	      x_ptr = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
442827bd09bSSatish Balay 	      rvec_copy(x_ptr,x,xxt_nnz);
443a501084fSBarry Smith 	      free(x);
444827bd09bSSatish Balay 	      x = x_ptr;
445827bd09bSSatish Balay 	      x_ptr+=xxt_nnz;
446827bd09bSSatish Balay 	    }
447827bd09bSSatish Balay 	  xxt_nnz += len;
448827bd09bSSatish Balay 	  rvec_copy(x_ptr,v+off,len);
449827bd09bSSatish Balay 
450827bd09bSSatish Balay           /* keep track of number of zeros */
451827bd09bSSatish Balay 	  if (dim)
452827bd09bSSatish Balay 	    {
453827bd09bSSatish Balay 	      for (k=0; k<len; k++)
454827bd09bSSatish Balay 		{
455827bd09bSSatish Balay 		  if (x_ptr[k]==0.0)
456827bd09bSSatish Balay 		    {xxt_zero_nnz++;}
457827bd09bSSatish Balay 		}
458827bd09bSSatish Balay 	    }
459827bd09bSSatish Balay 	  else
460827bd09bSSatish Balay 	    {
461827bd09bSSatish Balay 	      for (k=0; k<len; k++)
462827bd09bSSatish Balay 		{
463827bd09bSSatish Balay 		  if (x_ptr[k]==0.0)
464827bd09bSSatish Balay 		    {xxt_zero_nnz_0++;}
465827bd09bSSatish Balay 		}
466827bd09bSSatish Balay 	    }
467827bd09bSSatish Balay 	  col_indices[2*i] = off;
468827bd09bSSatish Balay 	  col_sz[i] = col_indices[2*i+1] = len;
469827bd09bSSatish Balay 	  col_vals[i] = x_ptr;
470827bd09bSSatish Balay 	}
471827bd09bSSatish Balay       else
472827bd09bSSatish Balay 	{
473827bd09bSSatish Balay 	  col_indices[2*i] = 0;
474827bd09bSSatish Balay 	  col_sz[i] = col_indices[2*i+1] = 0;
475827bd09bSSatish Balay 	  col_vals[i] = x_ptr;
476827bd09bSSatish Balay 	}
477827bd09bSSatish Balay     }
478827bd09bSSatish Balay 
479827bd09bSSatish Balay   /* close off stages for execution phase */
480827bd09bSSatish Balay   while (dim!=level)
481827bd09bSSatish Balay     {
482827bd09bSSatish Balay       stages[dim++]=i;
48371a0148aSBarry Smith       ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
484827bd09bSSatish Balay     }
485827bd09bSSatish Balay   stages[dim]=i;
486827bd09bSSatish Balay 
487827bd09bSSatish Balay   xxt_handle->info->n=xxt_handle->mvi->n;
488827bd09bSSatish Balay   xxt_handle->info->m=m;
489827bd09bSSatish Balay   xxt_handle->info->nnz=xxt_nnz;
490827bd09bSSatish Balay   xxt_handle->info->max_nnz=xxt_max_nnz;
491827bd09bSSatish Balay   xxt_handle->info->msg_buf_sz=stages[level]-stages[0];
492a501084fSBarry Smith   xxt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
493a501084fSBarry Smith   xxt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
494827bd09bSSatish Balay   xxt_handle->info->x=x;
495827bd09bSSatish Balay   xxt_handle->info->col_vals=col_vals;
496827bd09bSSatish Balay   xxt_handle->info->col_sz=col_sz;
497827bd09bSSatish Balay   xxt_handle->info->col_indices=col_indices;
498827bd09bSSatish Balay   xxt_handle->info->stages=stages;
499827bd09bSSatish Balay   xxt_handle->info->nsolves=0;
500827bd09bSSatish Balay   xxt_handle->info->tot_solve_time=0.0;
501827bd09bSSatish Balay 
502a501084fSBarry Smith   free(segs);
503a501084fSBarry Smith   free(u);
504a501084fSBarry Smith   free(v);
505a501084fSBarry Smith   free(uu);
506a501084fSBarry Smith   free(z);
507a501084fSBarry Smith   free(w);
508827bd09bSSatish Balay 
509827bd09bSSatish Balay   return(0);
510827bd09bSSatish Balay }
511827bd09bSSatish Balay 
5127b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5130924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle,  PetscScalar *uc)
514827bd09bSSatish Balay {
51552f87cdaSBarry Smith    PetscInt off, len, *iptr;
51652f87cdaSBarry Smith   PetscInt level       =xxt_handle->level;
51752f87cdaSBarry Smith   PetscInt n           =xxt_handle->info->n;
51852f87cdaSBarry Smith   PetscInt m           =xxt_handle->info->m;
51952f87cdaSBarry Smith   PetscInt *stages     =xxt_handle->info->stages;
52052f87cdaSBarry Smith   PetscInt *col_indices=xxt_handle->info->col_indices;
521a501084fSBarry Smith   PetscScalar *x_ptr, *uu_ptr;
522a501084fSBarry Smith   PetscScalar *solve_uu=xxt_handle->info->solve_uu;
523a501084fSBarry Smith   PetscScalar *solve_w =xxt_handle->info->solve_w;
524a501084fSBarry Smith   PetscScalar *x       =xxt_handle->info->x;
52571a0148aSBarry Smith   PetscBLASInt i1 = 1,dlen;
526827bd09bSSatish Balay 
5270924e98cSBarry Smith   PetscFunctionBegin;
528827bd09bSSatish Balay   uu_ptr=solve_uu;
529827bd09bSSatish Balay   rvec_zero(uu_ptr,m);
530827bd09bSSatish Balay 
531827bd09bSSatish Balay   /* x  = X.Y^T.b */
532827bd09bSSatish Balay   /* uu = Y^T.b */
533827bd09bSSatish Balay   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
534827bd09bSSatish Balay     {
535827bd09bSSatish Balay       off=*iptr++; len=*iptr++;
5360805154bSBarry Smith       dlen = PetscBLASIntCast(len);
53771a0148aSBarry Smith       *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1);
538827bd09bSSatish Balay     }
539827bd09bSSatish Balay 
540827bd09bSSatish Balay   /* comunication of beta */
541827bd09bSSatish Balay   uu_ptr=solve_uu;
542827bd09bSSatish Balay   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}
543827bd09bSSatish Balay 
544827bd09bSSatish Balay   rvec_zero(uc,n);
545827bd09bSSatish Balay 
546827bd09bSSatish Balay   /* x = X.uu */
547827bd09bSSatish Balay   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
548827bd09bSSatish Balay     {
549827bd09bSSatish Balay       off=*iptr++; len=*iptr++;
5500805154bSBarry Smith       dlen = PetscBLASIntCast(len);
55171a0148aSBarry Smith       BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
552827bd09bSSatish Balay     }
5530924e98cSBarry Smith   PetscFunctionReturn(0);
554827bd09bSSatish Balay }
555827bd09bSSatish Balay 
5567b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5570924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle)
558827bd09bSSatish Balay {
55952f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
560827bd09bSSatish Balay 
5610924e98cSBarry Smith   PetscFunctionBegin;
562827bd09bSSatish Balay   if (xxt_handle==NULL)
56371a0148aSBarry Smith     {SETERRQ1(PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xxt_handle);}
564827bd09bSSatish Balay 
565827bd09bSSatish Balay   vals[0]=vals[1]=xxt_handle->id;
566827bd09bSSatish Balay   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
567827bd09bSSatish Balay   if ((vals[0]!=vals[1])||(xxt_handle->id<=0))
56871a0148aSBarry Smith     {SETERRQ3(PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %D/%D %D\n",vals[0],vals[1], xxt_handle->id);}
5690924e98cSBarry Smith   PetscFunctionReturn(0);
570827bd09bSSatish Balay }
571827bd09bSSatish Balay 
5727b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5730924e98cSBarry Smith static  PetscErrorCode det_separators(xxt_ADT xxt_handle)
574827bd09bSSatish Balay {
57552f87cdaSBarry Smith   PetscInt i, ct, id;
57652f87cdaSBarry Smith   PetscInt mask, edge, *iptr;
57752f87cdaSBarry Smith   PetscInt *dir, *used;
57852f87cdaSBarry Smith   PetscInt sum[4], w[4];
579a501084fSBarry Smith   PetscScalar rsum[4], rw[4];
58052f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
581a501084fSBarry Smith   PetscScalar *lhs, *rhs;
58252f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo=0;
583827bd09bSSatish Balay   gs_ADT gs_handle=xxt_handle->mvi->gs_handle;
58452f87cdaSBarry Smith   PetscInt *local2global=xxt_handle->mvi->local2global;
58552f87cdaSBarry Smith   PetscInt  n=xxt_handle->mvi->n;
58652f87cdaSBarry Smith   PetscInt  m=xxt_handle->mvi->m;
58752f87cdaSBarry Smith   PetscInt level=xxt_handle->level;
58852f87cdaSBarry Smith   PetscInt shared=FALSE;
589827bd09bSSatish Balay 
5900924e98cSBarry Smith   PetscFunctionBegin;
59152f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
59252f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
59352f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
59452f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
59552f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
596827bd09bSSatish Balay 
597827bd09bSSatish Balay   ivec_zero(dir  ,level+1);
598827bd09bSSatish Balay   ivec_zero(nsep ,level+1);
599827bd09bSSatish Balay   ivec_zero(lnsep,level+1);
600827bd09bSSatish Balay   ivec_set (fo   ,-1,n+1);
601827bd09bSSatish Balay   ivec_zero(used,n);
602827bd09bSSatish Balay 
603*8cda6cd7SBarry Smith   lhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
604*8cda6cd7SBarry Smith   rhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
605827bd09bSSatish Balay 
606827bd09bSSatish Balay   /* determine the # of unique dof */
607827bd09bSSatish Balay   rvec_zero(lhs,m);
608827bd09bSSatish Balay   rvec_set(lhs,1.0,n);
609827bd09bSSatish Balay   gs_gop_hc(gs_handle,lhs,"+\0",level);
610827bd09bSSatish Balay   rvec_zero(rsum,2);
611827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
612827bd09bSSatish Balay     {
613827bd09bSSatish Balay       if (lhs[i]!=0.0)
614827bd09bSSatish Balay 	{rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
615827bd09bSSatish Balay     }
616827bd09bSSatish Balay   grop_hc(rsum,rw,2,op,level);
617827bd09bSSatish Balay   rsum[0]+=0.1;
618827bd09bSSatish Balay   rsum[1]+=0.1;
619827bd09bSSatish Balay 
620827bd09bSSatish Balay   if (fabs(rsum[0]-rsum[1])>EPS)
621827bd09bSSatish Balay     {shared=TRUE;}
622827bd09bSSatish Balay 
62352f87cdaSBarry Smith   xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0];
62452f87cdaSBarry Smith   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0];
625827bd09bSSatish Balay 
626827bd09bSSatish Balay   /* determine separator sets top down */
627827bd09bSSatish Balay   if (shared)
628827bd09bSSatish Balay     {
629827bd09bSSatish Balay       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
630827bd09bSSatish Balay 	{
631827bd09bSSatish Balay 	  /* set rsh of hc, fire, and collect lhs responses */
632827bd09bSSatish Balay 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
633827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
634827bd09bSSatish Balay 
635827bd09bSSatish Balay 	  /* set lsh of hc, fire, and collect rhs responses */
636827bd09bSSatish Balay 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
637827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
638827bd09bSSatish Balay 
639827bd09bSSatish Balay 	  for (i=0;i<n;i++)
640827bd09bSSatish Balay 	    {
641827bd09bSSatish Balay 	      if (id< mask)
642827bd09bSSatish Balay 		{
643827bd09bSSatish Balay 		  if (lhs[i]!=0.0)
644827bd09bSSatish Balay 		    {lhs[i]=1.0;}
645827bd09bSSatish Balay 		}
646827bd09bSSatish Balay 	      if (id>=mask)
647827bd09bSSatish Balay 		{
648827bd09bSSatish Balay 		  if (rhs[i]!=0.0)
649827bd09bSSatish Balay 		    {rhs[i]=1.0;}
650827bd09bSSatish Balay 		}
651827bd09bSSatish Balay 	    }
652827bd09bSSatish Balay 
653827bd09bSSatish Balay 	  if (id< mask)
654827bd09bSSatish Balay 	    {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
655827bd09bSSatish Balay 	  else
656827bd09bSSatish Balay 	    {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}
657827bd09bSSatish Balay 
658827bd09bSSatish Balay 	  /* count number of dofs I own that have signal and not in sep set */
659827bd09bSSatish Balay 	  rvec_zero(rsum,4);
660827bd09bSSatish Balay 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
661827bd09bSSatish Balay 	    {
662827bd09bSSatish Balay 	      if (!used[i])
663827bd09bSSatish Balay 		{
664827bd09bSSatish Balay 		  /* number of unmarked dofs on node */
665827bd09bSSatish Balay 		  ct++;
666827bd09bSSatish Balay 		  /* number of dofs to be marked on lhs hc */
667827bd09bSSatish Balay 		  if (id< mask)
668827bd09bSSatish Balay 		    {
669827bd09bSSatish Balay 		      if (lhs[i]!=0.0)
670827bd09bSSatish Balay 			{sum[0]++; rsum[0]+=1.0/lhs[i];}
671827bd09bSSatish Balay 		    }
672827bd09bSSatish Balay 		  /* number of dofs to be marked on rhs hc */
673827bd09bSSatish Balay 		  if (id>=mask)
674827bd09bSSatish Balay 		    {
675827bd09bSSatish Balay 		      if (rhs[i]!=0.0)
676827bd09bSSatish Balay 			{sum[1]++; rsum[1]+=1.0/rhs[i];}
677827bd09bSSatish Balay 		    }
678827bd09bSSatish Balay 		}
679827bd09bSSatish Balay 	    }
680827bd09bSSatish Balay 
681827bd09bSSatish Balay 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
682827bd09bSSatish Balay 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
683827bd09bSSatish Balay 	  (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
684827bd09bSSatish Balay 	  giop_hc(sum,w,4,op,edge);
685827bd09bSSatish Balay 	  grop_hc(rsum,rw,4,op,edge);
686827bd09bSSatish Balay 	  rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
687827bd09bSSatish Balay 
688827bd09bSSatish Balay 	  if (id<mask)
689827bd09bSSatish Balay 	    {
690827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
691827bd09bSSatish Balay 	      for (ct=i=0;i<n;i++)
692827bd09bSSatish Balay 		{
693827bd09bSSatish Balay 		  if ((!used[i])&&(lhs[i]!=0.0))
694827bd09bSSatish Balay 		    {
695827bd09bSSatish Balay 		      ct++; nfo++;
696827bd09bSSatish Balay 
697827bd09bSSatish Balay 		      if (nfo>n)
698388eb383SBarry Smith 			{SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}
699827bd09bSSatish Balay 
700827bd09bSSatish Balay 		      *--iptr = local2global[i];
701827bd09bSSatish Balay 		      used[i]=edge;
702827bd09bSSatish Balay 		    }
703827bd09bSSatish Balay 		}
704827bd09bSSatish Balay 	      if (ct>1) {ivec_sort(iptr,ct);}
705827bd09bSSatish Balay 
706827bd09bSSatish Balay 	      lnsep[edge]=ct;
70752f87cdaSBarry Smith 	      nsep[edge]=(PetscInt) rsum[0];
708827bd09bSSatish Balay 	      dir [edge]=LEFT;
709827bd09bSSatish Balay 	    }
710827bd09bSSatish Balay 
711827bd09bSSatish Balay 	  if (id>=mask)
712827bd09bSSatish Balay 	    {
713827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
714827bd09bSSatish Balay 	      for (ct=i=0;i<n;i++)
715827bd09bSSatish Balay 		{
716827bd09bSSatish Balay 		  if ((!used[i])&&(rhs[i]!=0.0))
717827bd09bSSatish Balay 		    {
718827bd09bSSatish Balay 		      ct++; nfo++;
719827bd09bSSatish Balay 
720827bd09bSSatish Balay 		      if (nfo>n)
721388eb383SBarry Smith 			{SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}
722827bd09bSSatish Balay 
723827bd09bSSatish Balay 		      *--iptr = local2global[i];
724827bd09bSSatish Balay 		      used[i]=edge;
725827bd09bSSatish Balay 		    }
726827bd09bSSatish Balay 		}
727827bd09bSSatish Balay 	      if (ct>1) {ivec_sort(iptr,ct);}
728827bd09bSSatish Balay 
729827bd09bSSatish Balay 	      lnsep[edge]=ct;
73052f87cdaSBarry Smith 	      nsep[edge]= (PetscInt) rsum[1];
731827bd09bSSatish Balay 	      dir [edge]=RIGHT;
732827bd09bSSatish Balay 	    }
733827bd09bSSatish Balay 
734827bd09bSSatish Balay 	  /* LATER or we can recur on these to order seps at this level */
735827bd09bSSatish Balay 	  /* do we need full set of separators for this?                */
736827bd09bSSatish Balay 
737827bd09bSSatish Balay 	  /* fold rhs hc into lower */
738827bd09bSSatish Balay 	  if (id>=mask)
739827bd09bSSatish Balay 	    {id-=mask;}
740827bd09bSSatish Balay 	}
741827bd09bSSatish Balay     }
742827bd09bSSatish Balay   else
743827bd09bSSatish Balay     {
744827bd09bSSatish Balay       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
745827bd09bSSatish Balay 	{
746827bd09bSSatish Balay 	  /* set rsh of hc, fire, and collect lhs responses */
747827bd09bSSatish Balay 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
748827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
749827bd09bSSatish Balay 
750827bd09bSSatish Balay 	  /* set lsh of hc, fire, and collect rhs responses */
751827bd09bSSatish Balay 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
752827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
753827bd09bSSatish Balay 
754827bd09bSSatish Balay 	  /* count number of dofs I own that have signal and not in sep set */
755827bd09bSSatish Balay 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
756827bd09bSSatish Balay 	    {
757827bd09bSSatish Balay 	      if (!used[i])
758827bd09bSSatish Balay 		{
759827bd09bSSatish Balay 		  /* number of unmarked dofs on node */
760827bd09bSSatish Balay 		  ct++;
761827bd09bSSatish Balay 		  /* number of dofs to be marked on lhs hc */
762827bd09bSSatish Balay 		  if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
763827bd09bSSatish Balay 		  /* number of dofs to be marked on rhs hc */
764827bd09bSSatish Balay 		  if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
765827bd09bSSatish Balay 		}
766827bd09bSSatish Balay 	    }
767827bd09bSSatish Balay 
768827bd09bSSatish Balay 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
769827bd09bSSatish Balay 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
770827bd09bSSatish Balay 	  giop_hc(sum,w,4,op,edge);
771827bd09bSSatish Balay 
772827bd09bSSatish Balay 	  /* lhs hc wins */
773827bd09bSSatish Balay 	  if (sum[2]>=sum[3])
774827bd09bSSatish Balay 	    {
775827bd09bSSatish Balay 	      if (id<mask)
776827bd09bSSatish Balay 		{
777827bd09bSSatish Balay 		  /* mark dofs I own that have signal and not in sep set */
778827bd09bSSatish Balay 		  for (ct=i=0;i<n;i++)
779827bd09bSSatish Balay 		    {
780827bd09bSSatish Balay 		      if ((!used[i])&&(lhs[i]!=0.0))
781827bd09bSSatish Balay 			{
782827bd09bSSatish Balay 			  ct++; nfo++;
783827bd09bSSatish Balay 			  *--iptr = local2global[i];
784827bd09bSSatish Balay 			  used[i]=edge;
785827bd09bSSatish Balay 			}
786827bd09bSSatish Balay 		    }
787827bd09bSSatish Balay 		  if (ct>1) {ivec_sort(iptr,ct);}
788827bd09bSSatish Balay 		  lnsep[edge]=ct;
789827bd09bSSatish Balay 		}
790827bd09bSSatish Balay 	      nsep[edge]=sum[0];
791827bd09bSSatish Balay 	      dir [edge]=LEFT;
792827bd09bSSatish Balay 	    }
793827bd09bSSatish Balay 	  /* rhs hc wins */
794827bd09bSSatish Balay 	  else
795827bd09bSSatish Balay 	    {
796827bd09bSSatish Balay 	      if (id>=mask)
797827bd09bSSatish Balay 		{
798827bd09bSSatish Balay 		  /* mark dofs I own that have signal and not in sep set */
799827bd09bSSatish Balay 		  for (ct=i=0;i<n;i++)
800827bd09bSSatish Balay 		    {
801827bd09bSSatish Balay 		      if ((!used[i])&&(rhs[i]!=0.0))
802827bd09bSSatish Balay 			{
803827bd09bSSatish Balay 			  ct++; nfo++;
804827bd09bSSatish Balay 			  *--iptr = local2global[i];
805827bd09bSSatish Balay 			  used[i]=edge;
806827bd09bSSatish Balay 			}
807827bd09bSSatish Balay 		    }
808827bd09bSSatish Balay 		  if (ct>1) {ivec_sort(iptr,ct);}
809827bd09bSSatish Balay 		  lnsep[edge]=ct;
810827bd09bSSatish Balay 		}
811827bd09bSSatish Balay 	      nsep[edge]=sum[1];
812827bd09bSSatish Balay 	      dir [edge]=RIGHT;
813827bd09bSSatish Balay 	    }
814827bd09bSSatish Balay 	  /* LATER or we can recur on these to order seps at this level */
815827bd09bSSatish Balay 	  /* do we need full set of separators for this?                */
816827bd09bSSatish Balay 
817827bd09bSSatish Balay 	  /* fold rhs hc into lower */
818827bd09bSSatish Balay 	  if (id>=mask)
819827bd09bSSatish Balay 	    {id-=mask;}
820827bd09bSSatish Balay 	}
821827bd09bSSatish Balay     }
822827bd09bSSatish Balay 
823827bd09bSSatish Balay 
824827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
825827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
826827bd09bSSatish Balay     {
827827bd09bSSatish Balay       if (!used[i])
828827bd09bSSatish Balay 	{
829827bd09bSSatish Balay 	  ct++; nfo++;
830827bd09bSSatish Balay 	  *--iptr = local2global[i];
831827bd09bSSatish Balay 	  used[i]=edge;
832827bd09bSSatish Balay 	}
833827bd09bSSatish Balay     }
834827bd09bSSatish Balay   if (ct>1) {ivec_sort(iptr,ct);}
835827bd09bSSatish Balay   lnsep[edge]=ct;
836827bd09bSSatish Balay   nsep [edge]=ct;
837827bd09bSSatish Balay   dir  [edge]=LEFT;
838827bd09bSSatish Balay 
839827bd09bSSatish Balay   xxt_handle->info->nsep=nsep;
840827bd09bSSatish Balay   xxt_handle->info->lnsep=lnsep;
841827bd09bSSatish Balay   xxt_handle->info->fo=fo;
842827bd09bSSatish Balay   xxt_handle->info->nfo=nfo;
843827bd09bSSatish Balay 
844a501084fSBarry Smith   free(dir);
845a501084fSBarry Smith   free(lhs);
846a501084fSBarry Smith   free(rhs);
847a501084fSBarry Smith   free(used);
8480924e98cSBarry Smith   PetscFunctionReturn(0);
849827bd09bSSatish Balay }
850827bd09bSSatish Balay 
8517b1ae94cSBarry Smith /**************************************xxt.c***********************************/
85252f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
853827bd09bSSatish Balay {
854827bd09bSSatish Balay   mv_info *mvi;
855827bd09bSSatish Balay 
856827bd09bSSatish Balay 
857a501084fSBarry Smith   mvi = (mv_info*)malloc(sizeof(mv_info));
858827bd09bSSatish Balay   mvi->n=n;
859827bd09bSSatish Balay   mvi->m=m;
860827bd09bSSatish Balay   mvi->n_global=-1;
861827bd09bSSatish Balay   mvi->m_global=-1;
86252f87cdaSBarry Smith   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
863827bd09bSSatish Balay   ivec_copy(mvi->local2global,local2global,m);
864827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
865a501084fSBarry Smith   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
866827bd09bSSatish Balay   mvi->grid_data=grid_data;
867827bd09bSSatish Balay 
868827bd09bSSatish Balay   /* set xxt communication handle to perform restricted matvec */
869827bd09bSSatish Balay   mvi->gs_handle = gs_init(local2global, m, num_nodes);
870827bd09bSSatish Balay 
871827bd09bSSatish Balay   return(mvi);
872827bd09bSSatish Balay }
873827bd09bSSatish Balay 
8747b1ae94cSBarry Smith /**************************************xxt.c***********************************/
8750924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
876827bd09bSSatish Balay {
8770924e98cSBarry Smith   PetscFunctionBegin;
878827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
8790924e98cSBarry Smith   PetscFunctionReturn(0);
880827bd09bSSatish Balay }
881827bd09bSSatish Balay 
882827bd09bSSatish Balay 
883827bd09bSSatish Balay 
884