xref: /petsc/src/ksp/pc/impls/tfs/xyt.c (revision d1528f56e7efe90e41719e8af1e4ce61b838432e)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2827bd09bSSatish Balay 
3827bd09bSSatish Balay /*************************************xyt.c************************************
4827bd09bSSatish Balay Module Name: xyt
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 **************************************xyt.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 xyt_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 *xcol_sz, *xcol_indices;
31a501084fSBarry Smith   PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
3252f87cdaSBarry Smith   PetscInt *ycol_sz, *ycol_indices;
33a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
3452f87cdaSBarry Smith   PetscInt nsolves;
35a501084fSBarry Smith   PetscScalar tot_solve_time;
36827bd09bSSatish Balay } xyt_info;
37827bd09bSSatish Balay 
38827bd09bSSatish Balay 
39827bd09bSSatish Balay typedef struct matvec_info {
4052f87cdaSBarry Smith   PetscInt n, m, n_global, m_global;
4152f87cdaSBarry Smith   PetscInt *local2global;
42827bd09bSSatish Balay   gs_ADT gs_handle;
43a501084fSBarry Smith   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
44827bd09bSSatish Balay   void *grid_data;
45827bd09bSSatish Balay } mv_info;
46827bd09bSSatish Balay 
47827bd09bSSatish Balay struct xyt_CDT{
4852f87cdaSBarry Smith   PetscInt id;
4952f87cdaSBarry Smith   PetscInt ns;
5052f87cdaSBarry Smith   PetscInt level;
51827bd09bSSatish Balay   xyt_info *info;
52827bd09bSSatish Balay   mv_info  *mvi;
53827bd09bSSatish Balay };
54827bd09bSSatish Balay 
5552f87cdaSBarry Smith static PetscInt n_xyt=0;
5652f87cdaSBarry Smith static PetscInt n_xyt_handles=0;
57827bd09bSSatish Balay 
58827bd09bSSatish Balay /* prototypes */
593fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
603fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle);
613fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle);
623fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
6352f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle);
6452f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle);
6552f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data);
66827bd09bSSatish Balay 
677b1ae94cSBarry Smith /**************************************xyt.c***********************************/
687b1ae94cSBarry Smith xyt_ADT XYT_new(void)
69827bd09bSSatish Balay {
70827bd09bSSatish Balay   xyt_ADT xyt_handle;
71827bd09bSSatish Balay 
72827bd09bSSatish Balay   /* rolling count on n_xyt ... pot. problem here */
73827bd09bSSatish Balay   n_xyt_handles++;
74a501084fSBarry Smith   xyt_handle       = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
75827bd09bSSatish Balay   xyt_handle->id   = ++n_xyt;
76827bd09bSSatish Balay   xyt_handle->info = NULL;
77827bd09bSSatish Balay   xyt_handle->mvi  = NULL;
78827bd09bSSatish Balay 
79827bd09bSSatish Balay   return(xyt_handle);
80827bd09bSSatish Balay }
81827bd09bSSatish Balay 
827b1ae94cSBarry Smith /**************************************xyt.c***********************************/
8352f87cdaSBarry Smith PetscInt XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt  handle */
8452f87cdaSBarry Smith 	   PetscInt *local2global,  /* global column mapping       */
8552f87cdaSBarry Smith 	   PetscInt n,              /* local num rows              */
8652f87cdaSBarry Smith 	   PetscInt m,              /* local num cols              */
87827bd09bSSatish Balay 	   void *matvec,       /* b_loc=A_local.x_loc         */
88827bd09bSSatish Balay 	   void *grid_data     /* grid data for matvec        */
89827bd09bSSatish Balay 	   )
90827bd09bSSatish Balay {
91827bd09bSSatish Balay 
92a40341efSBarry Smith   comm_init();
93827bd09bSSatish Balay   check_handle(xyt_handle);
94827bd09bSSatish Balay 
95827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
96827bd09bSSatish Balay   if ((1<<(xyt_handle->level=i_log2_num_nodes))!=num_nodes)
97388eb383SBarry Smith     {SETERRQ2(PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %d != %d\n",1<<i_log2_num_nodes,num_nodes);}
98827bd09bSSatish Balay 
99827bd09bSSatish Balay   /* space for X info */
100a501084fSBarry Smith   xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));
101827bd09bSSatish Balay 
102827bd09bSSatish Balay   /* set up matvec handles */
103827bd09bSSatish Balay   xyt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);
104827bd09bSSatish Balay 
105827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
106827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
107827bd09bSSatish Balay   xyt_handle->ns=0;
108827bd09bSSatish Balay 
109827bd09bSSatish Balay   /* determine separators and generate firing order - NB xyt info set here */
110827bd09bSSatish Balay   det_separators(xyt_handle);
111827bd09bSSatish Balay 
112827bd09bSSatish Balay   return(do_xyt_factor(xyt_handle));
113827bd09bSSatish Balay }
114827bd09bSSatish Balay 
1157b1ae94cSBarry Smith /**************************************xyt.c***********************************/
11652f87cdaSBarry Smith PetscInt XYT_solve(xyt_ADT xyt_handle, double *x, double *b)
117827bd09bSSatish Balay {
118a40341efSBarry Smith   comm_init();
119827bd09bSSatish Balay   check_handle(xyt_handle);
120827bd09bSSatish Balay 
121827bd09bSSatish Balay   /* need to copy b into x? */
122827bd09bSSatish Balay   if (b)
123827bd09bSSatish Balay     {rvec_copy(x,b,xyt_handle->mvi->n);}
124827bd09bSSatish Balay   do_xyt_solve(xyt_handle,x);
125827bd09bSSatish Balay 
126827bd09bSSatish Balay   return(0);
127827bd09bSSatish Balay }
128827bd09bSSatish Balay 
1297b1ae94cSBarry Smith /**************************************xyt.c***********************************/
13052f87cdaSBarry Smith PetscInt XYT_free(xyt_ADT xyt_handle)
131827bd09bSSatish Balay {
132a40341efSBarry Smith   comm_init();
133827bd09bSSatish Balay   check_handle(xyt_handle);
134827bd09bSSatish Balay   n_xyt_handles--;
135827bd09bSSatish Balay 
136a501084fSBarry Smith   free(xyt_handle->info->nsep);
137a501084fSBarry Smith   free(xyt_handle->info->lnsep);
138a501084fSBarry Smith   free(xyt_handle->info->fo);
139a501084fSBarry Smith   free(xyt_handle->info->stages);
140a501084fSBarry Smith   free(xyt_handle->info->solve_uu);
141a501084fSBarry Smith   free(xyt_handle->info->solve_w);
142a501084fSBarry Smith   free(xyt_handle->info->x);
143a501084fSBarry Smith   free(xyt_handle->info->xcol_vals);
144a501084fSBarry Smith   free(xyt_handle->info->xcol_sz);
145a501084fSBarry Smith   free(xyt_handle->info->xcol_indices);
146a501084fSBarry Smith   free(xyt_handle->info->y);
147a501084fSBarry Smith   free(xyt_handle->info->ycol_vals);
148a501084fSBarry Smith   free(xyt_handle->info->ycol_sz);
149a501084fSBarry Smith   free(xyt_handle->info->ycol_indices);
150a501084fSBarry Smith   free(xyt_handle->info);
151a501084fSBarry Smith   free(xyt_handle->mvi->local2global);
152827bd09bSSatish Balay   gs_free(xyt_handle->mvi->gs_handle);
153a501084fSBarry Smith   free(xyt_handle->mvi);
154a501084fSBarry Smith   free(xyt_handle);
155827bd09bSSatish Balay 
156827bd09bSSatish Balay 
157827bd09bSSatish Balay   /* if the check fails we nuke */
158a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
159827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
160827bd09bSSatish Balay   return(0);
161827bd09bSSatish Balay }
162827bd09bSSatish Balay 
1637b1ae94cSBarry Smith /**************************************xyt.c***********************************/
16452f87cdaSBarry Smith PetscInt XYT_stats(xyt_ADT xyt_handle)
165827bd09bSSatish Balay {
16652f87cdaSBarry Smith   PetscInt  op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
16752f87cdaSBarry Smith   PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
16852f87cdaSBarry Smith   PetscInt   vals[9],  work[9];
169a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
170827bd09bSSatish Balay 
171827bd09bSSatish Balay 
172a40341efSBarry Smith   comm_init();
173827bd09bSSatish Balay   check_handle(xyt_handle);
174827bd09bSSatish Balay 
175827bd09bSSatish Balay   /* if factorization not done there are no stats */
176827bd09bSSatish Balay   if (!xyt_handle->info||!xyt_handle->mvi)
177827bd09bSSatish Balay     {
178827bd09bSSatish Balay       if (!my_id)
179827bd09bSSatish Balay 	{printf("XYT_stats() :: no stats available!\n");}
180827bd09bSSatish Balay       return 1;
181827bd09bSSatish Balay     }
182827bd09bSSatish Balay 
183827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
184827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
185827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
186827bd09bSSatish Balay   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
187827bd09bSSatish Balay 
188827bd09bSSatish Balay   fvals[0]=fvals[1]=fvals[2]
189827bd09bSSatish Balay     =xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
190827bd09bSSatish Balay   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
191827bd09bSSatish Balay 
192827bd09bSSatish Balay   if (!my_id)
193827bd09bSSatish Balay     {
194827bd09bSSatish Balay       printf("%d :: min   xyt_nnz=%d\n",my_id,vals[0]);
195827bd09bSSatish Balay       printf("%d :: max   xyt_nnz=%d\n",my_id,vals[1]);
196827bd09bSSatish Balay       printf("%d :: avg   xyt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
197827bd09bSSatish Balay       printf("%d :: tot   xyt_nnz=%d\n",my_id,vals[2]);
198827bd09bSSatish Balay       printf("%d :: xyt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
199827bd09bSSatish Balay       printf("%d :: xyt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
200827bd09bSSatish Balay       printf("%d :: min   xyt_n  =%d\n",my_id,vals[3]);
201827bd09bSSatish Balay       printf("%d :: max   xyt_n  =%d\n",my_id,vals[4]);
202827bd09bSSatish Balay       printf("%d :: avg   xyt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
203827bd09bSSatish Balay       printf("%d :: tot   xyt_n  =%d\n",my_id,vals[5]);
204827bd09bSSatish Balay       printf("%d :: min   xyt_buf=%d\n",my_id,vals[6]);
205827bd09bSSatish Balay       printf("%d :: max   xyt_buf=%d\n",my_id,vals[7]);
206827bd09bSSatish Balay       printf("%d :: avg   xyt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
207827bd09bSSatish Balay       printf("%d :: min   xyt_slv=%g\n",my_id,fvals[0]);
208827bd09bSSatish Balay       printf("%d :: max   xyt_slv=%g\n",my_id,fvals[1]);
209827bd09bSSatish Balay       printf("%d :: avg   xyt_slv=%g\n",my_id,fvals[2]/num_nodes);
210827bd09bSSatish Balay     }
211827bd09bSSatish Balay 
212827bd09bSSatish Balay   return(0);
213827bd09bSSatish Balay }
214827bd09bSSatish Balay 
215827bd09bSSatish Balay 
216827bd09bSSatish Balay /*************************************xyt.c************************************
217827bd09bSSatish Balay 
218827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
219827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
220827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
221827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
222827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
223827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
224827bd09bSSatish Balay    gs_init/gop).
225827bd09bSSatish Balay 
226827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
227827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
228827bd09bSSatish Balay **************************************xyt.c***********************************/
22952f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle)
230827bd09bSSatish Balay {
2317b1ae94cSBarry Smith   return xyt_generate(xyt_handle);
232827bd09bSSatish Balay }
233827bd09bSSatish Balay 
2347b1ae94cSBarry Smith /**************************************xyt.c***********************************/
23552f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle)
236827bd09bSSatish Balay {
23752f87cdaSBarry Smith   PetscInt i,j,k,idx;
23852f87cdaSBarry Smith   PetscInt dim, col;
239a501084fSBarry Smith   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
24052f87cdaSBarry Smith   PetscInt *segs;
24152f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
24252f87cdaSBarry Smith   PetscInt off, len;
243a501084fSBarry Smith   PetscScalar *x_ptr, *y_ptr;
24452f87cdaSBarry Smith   PetscInt *iptr, flag;
24552f87cdaSBarry Smith   PetscInt start=0, end, work;
24652f87cdaSBarry Smith   PetscInt op2[] = {GL_MIN,0};
247827bd09bSSatish Balay   gs_ADT gs_handle;
24852f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo;
24952f87cdaSBarry Smith   PetscInt a_n=xyt_handle->mvi->n;
25052f87cdaSBarry Smith   PetscInt a_m=xyt_handle->mvi->m;
25152f87cdaSBarry Smith   PetscInt *a_local2global=xyt_handle->mvi->local2global;
25252f87cdaSBarry Smith   PetscInt level;
25352f87cdaSBarry Smith   PetscInt n, m;
25452f87cdaSBarry Smith   PetscInt *xcol_sz, *xcol_indices, *stages;
255a501084fSBarry Smith   PetscScalar **xcol_vals, *x;
25652f87cdaSBarry Smith   PetscInt *ycol_sz, *ycol_indices;
257a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
25852f87cdaSBarry Smith   PetscInt n_global;
25952f87cdaSBarry Smith   PetscInt xt_nnz=0, xt_max_nnz=0;
26052f87cdaSBarry Smith   PetscInt yt_nnz=0, yt_max_nnz=0;
26152f87cdaSBarry Smith   PetscInt xt_zero_nnz  =0;
26252f87cdaSBarry Smith   PetscInt xt_zero_nnz_0=0;
26352f87cdaSBarry Smith   PetscInt yt_zero_nnz  =0;
26452f87cdaSBarry Smith   PetscInt yt_zero_nnz_0=0;
265a501084fSBarry Smith   PetscBLASInt i1 = 1;
266a501084fSBarry Smith   PetscScalar dm1 = -1.0;
267*d1528f56SBarry Smith   PetscErrorCode ierr;
268827bd09bSSatish Balay 
269827bd09bSSatish Balay   n=xyt_handle->mvi->n;
270827bd09bSSatish Balay   nsep=xyt_handle->info->nsep;
271827bd09bSSatish Balay   lnsep=xyt_handle->info->lnsep;
272827bd09bSSatish Balay   fo=xyt_handle->info->fo;
273827bd09bSSatish Balay   end=lnsep[0];
274827bd09bSSatish Balay   level=xyt_handle->level;
275827bd09bSSatish Balay   gs_handle=xyt_handle->mvi->gs_handle;
276827bd09bSSatish Balay 
277827bd09bSSatish Balay   /* is there a null space? */
278827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
279827bd09bSSatish Balay   for (i=0, j=0; i<=level; i++)
280827bd09bSSatish Balay     {j+=nsep[i];}
281827bd09bSSatish Balay 
282827bd09bSSatish Balay   m = j-xyt_handle->ns;
283827bd09bSSatish Balay   if (m!=j)
284827bd09bSSatish Balay     {printf("xyt_generate() :: null space exists %d %d %d\n",m,j,xyt_handle->ns);}
285827bd09bSSatish Balay 
286*d1528f56SBarry Smith   ierr = PetscInfo2(0,"xyt_generate() :: X(%d,%d)\n",n,m);CHKERRQ(ierr);
287827bd09bSSatish Balay 
288827bd09bSSatish Balay   /* get and initialize storage for x local         */
289827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
29052f87cdaSBarry Smith   xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
29152f87cdaSBarry Smith   xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
292a501084fSBarry Smith   xcol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
293827bd09bSSatish Balay   for (i=j=0; i<m; i++, j+=2)
294827bd09bSSatish Balay     {
295827bd09bSSatish Balay       xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
296827bd09bSSatish Balay       xcol_vals[i] = NULL;
297827bd09bSSatish Balay     }
298827bd09bSSatish Balay   xcol_indices[j]=-1;
299827bd09bSSatish Balay 
300827bd09bSSatish Balay   /* get and initialize storage for y local         */
301827bd09bSSatish Balay   /* note that y local is nxm and stored by columns */
30252f87cdaSBarry Smith   ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
30352f87cdaSBarry Smith   ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
304a501084fSBarry Smith   ycol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
305827bd09bSSatish Balay   for (i=j=0; i<m; i++, j+=2)
306827bd09bSSatish Balay     {
307827bd09bSSatish Balay       ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
308827bd09bSSatish Balay       ycol_vals[i] = NULL;
309827bd09bSSatish Balay     }
310827bd09bSSatish Balay   ycol_indices[j]=-1;
311827bd09bSSatish Balay 
312827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
313827bd09bSSatish Balay   /* this looks like nsep[]=segments */
31452f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
31552f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
316827bd09bSSatish Balay   ivec_zero(stages,level+1);
317827bd09bSSatish Balay   ivec_copy(segs,nsep,level+1);
318827bd09bSSatish Balay   for (i=0; i<level; i++)
319827bd09bSSatish Balay     {segs[i+1] += segs[i];}
320827bd09bSSatish Balay   stages[0] = segs[0];
321827bd09bSSatish Balay 
322827bd09bSSatish Balay   /* temporary vectors  */
323a501084fSBarry Smith   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
324a501084fSBarry Smith   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
325a501084fSBarry Smith   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
326a501084fSBarry Smith   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
327a501084fSBarry Smith   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
328827bd09bSSatish Balay 
329827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
330827bd09bSSatish Balay   for (i=1, j=0; i<=level; i++)
331827bd09bSSatish Balay     {j+=nsep[i];}
332827bd09bSSatish Balay 
333827bd09bSSatish Balay   /* storage for sparse x values */
334827bd09bSSatish Balay   n_global = xyt_handle->info->n_global;
33552f87cdaSBarry Smith   xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
336a501084fSBarry Smith   x = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
337a501084fSBarry Smith   y = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
338827bd09bSSatish Balay 
339827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
340827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
341827bd09bSSatish Balay   for (dim=i=j=0;i<m;i++)
342827bd09bSSatish Balay     {
343827bd09bSSatish Balay       /* time to move to the next level? */
344827bd09bSSatish Balay       while (i==segs[dim])
345827bd09bSSatish Balay 	{
346827bd09bSSatish Balay 	  if (dim==level)
347388eb383SBarry Smith 	    {SETERRQ(PETSC_ERR_PLIB,"dim about to exceed level\n"); break;}
348827bd09bSSatish Balay 
349827bd09bSSatish Balay 	  stages[dim++]=i;
350827bd09bSSatish Balay 	  end+=lnsep[dim];
351827bd09bSSatish Balay 	}
352827bd09bSSatish Balay       stages[dim]=i;
353827bd09bSSatish Balay 
354827bd09bSSatish Balay       /* which column are we firing? */
355827bd09bSSatish Balay       /* i.e. set v_l */
356827bd09bSSatish Balay       /* use new seps and do global min across hc to determine which one to fire */
357827bd09bSSatish Balay       (start<end) ? (col=fo[start]) : (col=INT_MAX);
358827bd09bSSatish Balay       giop_hc(&col,&work,1,op2,dim);
359827bd09bSSatish Balay 
360827bd09bSSatish Balay       /* shouldn't need this */
361827bd09bSSatish Balay       if (col==INT_MAX)
362827bd09bSSatish Balay 	{
363*d1528f56SBarry Smith 	  ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");
364827bd09bSSatish Balay 	  continue;
365827bd09bSSatish Balay 	}
366827bd09bSSatish Balay 
367827bd09bSSatish Balay       /* do I own it? I should */
368827bd09bSSatish Balay       rvec_zero(v ,a_m);
369827bd09bSSatish Balay       if (col==fo[start])
370827bd09bSSatish Balay 	{
371827bd09bSSatish Balay 	  start++;
372827bd09bSSatish Balay 	  idx=ivec_linear_search(col, a_local2global, a_n);
373827bd09bSSatish Balay 	  if (idx!=-1)
374827bd09bSSatish Balay 	    {v[idx] = 1.0; j++;}
375827bd09bSSatish Balay 	  else
376388eb383SBarry Smith 	    {SETERRQ(PETSC_ERR_PLIB,"NOT FOUND!\n");}
377827bd09bSSatish Balay 	}
378827bd09bSSatish Balay       else
379827bd09bSSatish Balay 	{
380827bd09bSSatish Balay 	  idx=ivec_linear_search(col, a_local2global, a_m);
381827bd09bSSatish Balay 	  if (idx!=-1)
382827bd09bSSatish Balay 	    {v[idx] = 1.0;}
383827bd09bSSatish Balay 	}
384827bd09bSSatish Balay 
385827bd09bSSatish Balay       /* perform u = A.v_l */
386827bd09bSSatish Balay       rvec_zero(u,n);
387827bd09bSSatish Balay       do_matvec(xyt_handle->mvi,v,u);
388827bd09bSSatish Balay 
389827bd09bSSatish Balay       /* uu =  X^T.u_l (local portion) */
390827bd09bSSatish Balay       /* technically only need to zero out first i entries */
391827bd09bSSatish Balay       /* later turn this into an XYT_solve call ? */
392827bd09bSSatish Balay       rvec_zero(uu,m);
393827bd09bSSatish Balay       y_ptr=y;
394827bd09bSSatish Balay       iptr = ycol_indices;
395827bd09bSSatish Balay       for (k=0; k<i; k++)
396827bd09bSSatish Balay 	{
397827bd09bSSatish Balay 	  off = *iptr++;
398827bd09bSSatish Balay 	  len = *iptr++;
399827bd09bSSatish Balay 
40071044d3cSBarry Smith 	  uu[k] = BLASdot_(&len,u+off,&i1,y_ptr,&i1);
401827bd09bSSatish Balay 	  y_ptr+=len;
402827bd09bSSatish Balay 	}
403827bd09bSSatish Balay 
404827bd09bSSatish Balay       /* uu = X^T.u_l (comm portion) */
405827bd09bSSatish Balay       ssgl_radd  (uu, w, dim, stages);
406827bd09bSSatish Balay 
407827bd09bSSatish Balay       /* z = X.uu */
408827bd09bSSatish Balay       rvec_zero(z,n);
409827bd09bSSatish Balay       x_ptr=x;
410827bd09bSSatish Balay       iptr = xcol_indices;
411827bd09bSSatish Balay       for (k=0; k<i; k++)
412827bd09bSSatish Balay 	{
413827bd09bSSatish Balay 	  off = *iptr++;
414827bd09bSSatish Balay 	  len = *iptr++;
415827bd09bSSatish Balay 
41671044d3cSBarry Smith 	  BLASaxpy_(&len,&uu[k],x_ptr,&i1,z+off,&i1);
417827bd09bSSatish Balay 	  x_ptr+=len;
418827bd09bSSatish Balay 	}
419827bd09bSSatish Balay 
420827bd09bSSatish Balay       /* compute v_l = v_l - z */
421827bd09bSSatish Balay       rvec_zero(v+a_n,a_m-a_n);
42271044d3cSBarry Smith       BLASaxpy_(&n,&dm1,z,&i1,v,&i1);
423827bd09bSSatish Balay 
424827bd09bSSatish Balay       /* compute u_l = A.v_l */
425827bd09bSSatish Balay       if (a_n!=a_m)
426827bd09bSSatish Balay 	{gs_gop_hc(gs_handle,v,"+\0",dim);}
427827bd09bSSatish Balay       rvec_zero(u,n);
428827bd09bSSatish Balay      do_matvec(xyt_handle->mvi,v,u);
429827bd09bSSatish Balay 
430827bd09bSSatish Balay       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
43171044d3cSBarry Smith       alpha = BLASdot_(&n,u,&i1,u,&i1);
432827bd09bSSatish Balay       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
433827bd09bSSatish Balay       grop_hc(&alpha, &alpha_w, 1, op, dim);
434827bd09bSSatish Balay 
435a501084fSBarry Smith       alpha = (PetscScalar) sqrt((double)alpha);
436827bd09bSSatish Balay 
437827bd09bSSatish Balay       /* check for small alpha                             */
438827bd09bSSatish Balay       /* LATER use this to detect and determine null space */
439827bd09bSSatish Balay       if (fabs(alpha)<1.0e-14)
440388eb383SBarry Smith 	{SETERRQ1(PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);}
441827bd09bSSatish Balay 
442827bd09bSSatish Balay       /* compute v_l = v_l/sqrt(alpha) */
443827bd09bSSatish Balay       rvec_scale(v,1.0/alpha,n);
444827bd09bSSatish Balay       rvec_scale(u,1.0/alpha,n);
445827bd09bSSatish Balay 
446827bd09bSSatish Balay       /* add newly generated column, v_l, to X */
447827bd09bSSatish Balay       flag = 1;
448827bd09bSSatish Balay       off=len=0;
449827bd09bSSatish Balay       for (k=0; k<n; k++)
450827bd09bSSatish Balay 	{
451827bd09bSSatish Balay 	  if (v[k]!=0.0)
452827bd09bSSatish Balay 	    {
453827bd09bSSatish Balay 	      len=k;
454827bd09bSSatish Balay 	      if (flag)
455827bd09bSSatish Balay 		{off=k; flag=0;}
456827bd09bSSatish Balay 	    }
457827bd09bSSatish Balay 	}
458827bd09bSSatish Balay 
459827bd09bSSatish Balay       len -= (off-1);
460827bd09bSSatish Balay 
461827bd09bSSatish Balay       if (len>0)
462827bd09bSSatish Balay 	{
463827bd09bSSatish Balay 	  if ((xt_nnz+len)>xt_max_nnz)
464827bd09bSSatish Balay 	    {
465*d1528f56SBarry Smith 	      ierr = PetscInfo(0,"increasing space for X by 2x!\n");
466827bd09bSSatish Balay 	      xt_max_nnz *= 2;
467a501084fSBarry Smith 	      x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
468827bd09bSSatish Balay 	      rvec_copy(x_ptr,x,xt_nnz);
469a501084fSBarry Smith 	      free(x);
470827bd09bSSatish Balay 	      x = x_ptr;
471827bd09bSSatish Balay 	      x_ptr+=xt_nnz;
472827bd09bSSatish Balay 	    }
473827bd09bSSatish Balay 	  xt_nnz += len;
474827bd09bSSatish Balay 	  rvec_copy(x_ptr,v+off,len);
475827bd09bSSatish Balay 
476827bd09bSSatish Balay           /* keep track of number of zeros */
477827bd09bSSatish Balay 	  if (dim)
478827bd09bSSatish Balay 	    {
479827bd09bSSatish Balay 	      for (k=0; k<len; k++)
480827bd09bSSatish Balay 		{
481827bd09bSSatish Balay 		  if (x_ptr[k]==0.0)
482827bd09bSSatish Balay 		    {xt_zero_nnz++;}
483827bd09bSSatish Balay 		}
484827bd09bSSatish Balay 	    }
485827bd09bSSatish Balay 	  else
486827bd09bSSatish Balay 	    {
487827bd09bSSatish Balay 	      for (k=0; k<len; k++)
488827bd09bSSatish Balay 		{
489827bd09bSSatish Balay 		  if (x_ptr[k]==0.0)
490827bd09bSSatish Balay 		    {xt_zero_nnz_0++;}
491827bd09bSSatish Balay 		}
492827bd09bSSatish Balay 	    }
493827bd09bSSatish Balay 	  xcol_indices[2*i] = off;
494827bd09bSSatish Balay 	  xcol_sz[i] = xcol_indices[2*i+1] = len;
495827bd09bSSatish Balay 	  xcol_vals[i] = x_ptr;
496827bd09bSSatish Balay 	}
497827bd09bSSatish Balay       else
498827bd09bSSatish Balay 	{
499827bd09bSSatish Balay 	  xcol_indices[2*i] = 0;
500827bd09bSSatish Balay 	  xcol_sz[i] = xcol_indices[2*i+1] = 0;
501827bd09bSSatish Balay 	  xcol_vals[i] = x_ptr;
502827bd09bSSatish Balay 	}
503827bd09bSSatish Balay 
504827bd09bSSatish Balay 
505827bd09bSSatish Balay       /* add newly generated column, u_l, to Y */
506827bd09bSSatish Balay       flag = 1;
507827bd09bSSatish Balay       off=len=0;
508827bd09bSSatish Balay       for (k=0; k<n; k++)
509827bd09bSSatish Balay 	{
510827bd09bSSatish Balay 	  if (u[k]!=0.0)
511827bd09bSSatish Balay 	    {
512827bd09bSSatish Balay 	      len=k;
513827bd09bSSatish Balay 	      if (flag)
514827bd09bSSatish Balay 		{off=k; flag=0;}
515827bd09bSSatish Balay 	    }
516827bd09bSSatish Balay 	}
517827bd09bSSatish Balay 
518827bd09bSSatish Balay       len -= (off-1);
519827bd09bSSatish Balay 
520827bd09bSSatish Balay       if (len>0)
521827bd09bSSatish Balay 	{
522827bd09bSSatish Balay 	  if ((yt_nnz+len)>yt_max_nnz)
523827bd09bSSatish Balay 	    {
524*d1528f56SBarry Smith 	      ierr = PetscInfo(0,"increasing space for Y by 2x!\n");
525827bd09bSSatish Balay 	      yt_max_nnz *= 2;
526a501084fSBarry Smith 	      y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
527827bd09bSSatish Balay 	      rvec_copy(y_ptr,y,yt_nnz);
528a501084fSBarry Smith 	      free(y);
529827bd09bSSatish Balay 	      y = y_ptr;
530827bd09bSSatish Balay 	      y_ptr+=yt_nnz;
531827bd09bSSatish Balay 	    }
532827bd09bSSatish Balay 	  yt_nnz += len;
533827bd09bSSatish Balay 	  rvec_copy(y_ptr,u+off,len);
534827bd09bSSatish Balay 
535827bd09bSSatish Balay           /* keep track of number of zeros */
536827bd09bSSatish Balay 	  if (dim)
537827bd09bSSatish Balay 	    {
538827bd09bSSatish Balay 	      for (k=0; k<len; k++)
539827bd09bSSatish Balay 		{
540827bd09bSSatish Balay 		  if (y_ptr[k]==0.0)
541827bd09bSSatish Balay 		    {yt_zero_nnz++;}
542827bd09bSSatish Balay 		}
543827bd09bSSatish Balay 	    }
544827bd09bSSatish Balay 	  else
545827bd09bSSatish Balay 	    {
546827bd09bSSatish Balay 	      for (k=0; k<len; k++)
547827bd09bSSatish Balay 		{
548827bd09bSSatish Balay 		  if (y_ptr[k]==0.0)
549827bd09bSSatish Balay 		    {yt_zero_nnz_0++;}
550827bd09bSSatish Balay 		}
551827bd09bSSatish Balay 	    }
552827bd09bSSatish Balay 	  ycol_indices[2*i] = off;
553827bd09bSSatish Balay 	  ycol_sz[i] = ycol_indices[2*i+1] = len;
554827bd09bSSatish Balay 	  ycol_vals[i] = y_ptr;
555827bd09bSSatish Balay 	}
556827bd09bSSatish Balay       else
557827bd09bSSatish Balay 	{
558827bd09bSSatish Balay 	  ycol_indices[2*i] = 0;
559827bd09bSSatish Balay 	  ycol_sz[i] = ycol_indices[2*i+1] = 0;
560827bd09bSSatish Balay 	  ycol_vals[i] = y_ptr;
561827bd09bSSatish Balay 	}
562827bd09bSSatish Balay     }
563827bd09bSSatish Balay 
564827bd09bSSatish Balay   /* close off stages for execution phase */
565827bd09bSSatish Balay   while (dim!=level)
566827bd09bSSatish Balay     {
567827bd09bSSatish Balay       stages[dim++]=i;
568*d1528f56SBarry Smith       ierr = PetscInfo2(0,"disconnected!!! dim(%d)!=level(%d)\n",dim,level);CHKERRQ(ierr);
569827bd09bSSatish Balay     }
570827bd09bSSatish Balay   stages[dim]=i;
571827bd09bSSatish Balay 
572827bd09bSSatish Balay   xyt_handle->info->n=xyt_handle->mvi->n;
573827bd09bSSatish Balay   xyt_handle->info->m=m;
574827bd09bSSatish Balay   xyt_handle->info->nnz=xt_nnz + yt_nnz;
575827bd09bSSatish Balay   xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz;
576827bd09bSSatish Balay   xyt_handle->info->msg_buf_sz=stages[level]-stages[0];
577a501084fSBarry Smith   xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
578a501084fSBarry Smith   xyt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
579827bd09bSSatish Balay   xyt_handle->info->x=x;
580827bd09bSSatish Balay   xyt_handle->info->xcol_vals=xcol_vals;
581827bd09bSSatish Balay   xyt_handle->info->xcol_sz=xcol_sz;
582827bd09bSSatish Balay   xyt_handle->info->xcol_indices=xcol_indices;
583827bd09bSSatish Balay   xyt_handle->info->stages=stages;
584827bd09bSSatish Balay   xyt_handle->info->y=y;
585827bd09bSSatish Balay   xyt_handle->info->ycol_vals=ycol_vals;
586827bd09bSSatish Balay   xyt_handle->info->ycol_sz=ycol_sz;
587827bd09bSSatish Balay   xyt_handle->info->ycol_indices=ycol_indices;
588827bd09bSSatish Balay 
589a501084fSBarry Smith   free(segs);
590a501084fSBarry Smith   free(u);
591a501084fSBarry Smith   free(v);
592a501084fSBarry Smith   free(uu);
593a501084fSBarry Smith   free(z);
594a501084fSBarry Smith   free(w);
595827bd09bSSatish Balay 
596827bd09bSSatish Balay   return(0);
597827bd09bSSatish Balay }
598827bd09bSSatish Balay 
5997b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6000924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle,  PetscScalar *uc)
601827bd09bSSatish Balay {
60252f87cdaSBarry Smith   PetscInt off, len, *iptr;
60352f87cdaSBarry Smith   PetscInt level       =xyt_handle->level;
60452f87cdaSBarry Smith   PetscInt n           =xyt_handle->info->n;
60552f87cdaSBarry Smith   PetscInt m           =xyt_handle->info->m;
60652f87cdaSBarry Smith   PetscInt *stages     =xyt_handle->info->stages;
60752f87cdaSBarry Smith   PetscInt *xcol_indices=xyt_handle->info->xcol_indices;
60852f87cdaSBarry Smith   PetscInt *ycol_indices=xyt_handle->info->ycol_indices;
609a501084fSBarry Smith    PetscScalar *x_ptr, *y_ptr, *uu_ptr;
610a501084fSBarry Smith   PetscScalar *solve_uu=xyt_handle->info->solve_uu;
611a501084fSBarry Smith   PetscScalar *solve_w =xyt_handle->info->solve_w;
612a501084fSBarry Smith   PetscScalar *x       =xyt_handle->info->x;
613a501084fSBarry Smith   PetscScalar *y       =xyt_handle->info->y;
614a501084fSBarry Smith   PetscBLASInt i1 = 1;
615827bd09bSSatish Balay 
6160924e98cSBarry Smith   PetscFunctionBegin;
617827bd09bSSatish Balay   uu_ptr=solve_uu;
618827bd09bSSatish Balay   rvec_zero(uu_ptr,m);
619827bd09bSSatish Balay 
620827bd09bSSatish Balay   /* x  = X.Y^T.b */
621827bd09bSSatish Balay   /* uu = Y^T.b */
622827bd09bSSatish Balay   for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len)
623827bd09bSSatish Balay     {
624827bd09bSSatish Balay       off=*iptr++; len=*iptr++;
62571044d3cSBarry Smith       *uu_ptr++ = BLASdot_(&len,uc+off,&i1,y_ptr,&i1);
626827bd09bSSatish Balay     }
627827bd09bSSatish Balay 
628827bd09bSSatish Balay   /* comunication of beta */
629827bd09bSSatish Balay   uu_ptr=solve_uu;
630827bd09bSSatish Balay   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}
631827bd09bSSatish Balay 
632827bd09bSSatish Balay   rvec_zero(uc,n);
633827bd09bSSatish Balay 
634827bd09bSSatish Balay   /* x = X.uu */
635827bd09bSSatish Balay   for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len)
636827bd09bSSatish Balay     {
637827bd09bSSatish Balay       off=*iptr++; len=*iptr++;
63871044d3cSBarry Smith       BLASaxpy_(&len,uu_ptr++,x_ptr,&i1,uc+off,&i1);
639827bd09bSSatish Balay     }
6400924e98cSBarry Smith   PetscFunctionReturn(0);
641827bd09bSSatish Balay }
642827bd09bSSatish Balay 
6437b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6440924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle)
645827bd09bSSatish Balay {
64652f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
647827bd09bSSatish Balay 
6480924e98cSBarry Smith   PetscFunctionBegin;
649827bd09bSSatish Balay   if (xyt_handle==NULL)
650388eb383SBarry Smith     {SETERRQ1(PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %d\n",xyt_handle);}
651827bd09bSSatish Balay 
652827bd09bSSatish Balay   vals[0]=vals[1]=xyt_handle->id;
653827bd09bSSatish Balay   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
654827bd09bSSatish Balay   if ((vals[0]!=vals[1])||(xyt_handle->id<=0))
655388eb383SBarry Smith     {SETERRQ3(PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %d/%d %d\n", vals[0],vals[1], xyt_handle->id);}
6560924e98cSBarry Smith   PetscFunctionReturn(0);
657827bd09bSSatish Balay }
658827bd09bSSatish Balay 
6597b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6600924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle)
661827bd09bSSatish Balay {
66252f87cdaSBarry Smith   PetscInt i, ct, id;
66352f87cdaSBarry Smith   PetscInt mask, edge, *iptr;
66452f87cdaSBarry Smith   PetscInt *dir, *used;
66552f87cdaSBarry Smith   PetscInt sum[4], w[4];
666a501084fSBarry Smith   PetscScalar rsum[4], rw[4];
66752f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
668a501084fSBarry Smith   PetscScalar *lhs, *rhs;
66952f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo=0;
670827bd09bSSatish Balay   gs_ADT gs_handle=xyt_handle->mvi->gs_handle;
67152f87cdaSBarry Smith   PetscInt *local2global=xyt_handle->mvi->local2global;
67252f87cdaSBarry Smith   PetscInt  n=xyt_handle->mvi->n;
67352f87cdaSBarry Smith   PetscInt  m=xyt_handle->mvi->m;
67452f87cdaSBarry Smith   PetscInt level=xyt_handle->level;
67552f87cdaSBarry Smith   PetscInt shared=FALSE;
676*d1528f56SBarry Smith   PetscErrorCode ierr;
677827bd09bSSatish Balay 
6780924e98cSBarry Smith   PetscFunctionBegin;
67952f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
68052f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
68152f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
68252f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
68352f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
684827bd09bSSatish Balay 
685827bd09bSSatish Balay   ivec_zero(dir  ,level+1);
686827bd09bSSatish Balay   ivec_zero(nsep ,level+1);
687827bd09bSSatish Balay   ivec_zero(lnsep,level+1);
688827bd09bSSatish Balay   ivec_set (fo   ,-1,n+1);
689827bd09bSSatish Balay   ivec_zero(used,n);
690827bd09bSSatish Balay 
691a501084fSBarry Smith   lhs  = (double*)malloc(sizeof(PetscScalar)*m);
692a501084fSBarry Smith   rhs  = (double*)malloc(sizeof(PetscScalar)*m);
693827bd09bSSatish Balay 
694827bd09bSSatish Balay   /* determine the # of unique dof */
695827bd09bSSatish Balay   rvec_zero(lhs,m);
696827bd09bSSatish Balay   rvec_set(lhs,1.0,n);
697827bd09bSSatish Balay   gs_gop_hc(gs_handle,lhs,"+\0",level);
698*d1528f56SBarry Smith   ierr = PetscInfo(0,"done first gs_gop_hc\n");CHKERRQ(ierr);
699827bd09bSSatish Balay   rvec_zero(rsum,2);
700827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
701827bd09bSSatish Balay     {
702827bd09bSSatish Balay       if (lhs[i]!=0.0)
703827bd09bSSatish Balay 	{rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
704827bd09bSSatish Balay 
705827bd09bSSatish Balay       if (lhs[i]!=1.0)
706827bd09bSSatish Balay 	{
707827bd09bSSatish Balay           shared=TRUE;
708827bd09bSSatish Balay         }
709827bd09bSSatish Balay     }
710827bd09bSSatish Balay 
711827bd09bSSatish Balay   grop_hc(rsum,rw,2,op,level);
712827bd09bSSatish Balay   rsum[0]+=0.1;
713827bd09bSSatish Balay   rsum[1]+=0.1;
714827bd09bSSatish Balay 
715827bd09bSSatish Balay   /*
716827bd09bSSatish Balay       if (!my_id)
717827bd09bSSatish Balay       {
718827bd09bSSatish Balay       printf("xyt n unique = %d (%g)\n",(int) rsum[0], rsum[0]);
719827bd09bSSatish Balay       printf("xyt n shared = %d (%g)\n",(int) rsum[1], rsum[1]);
720827bd09bSSatish Balay       }
721827bd09bSSatish Balay   */
722827bd09bSSatish Balay 
72352f87cdaSBarry Smith   xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0];
72452f87cdaSBarry Smith   xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0];
725827bd09bSSatish Balay 
726827bd09bSSatish Balay   /* determine separator sets top down */
727827bd09bSSatish Balay   if (shared)
728827bd09bSSatish Balay     {
729827bd09bSSatish Balay       /* solution is to do as in the symmetric shared case but then */
730827bd09bSSatish Balay       /* pick the sub-hc with the most free dofs and do a mat-vec   */
731827bd09bSSatish Balay       /* and pick up the responses on the other sub-hc from the     */
732827bd09bSSatish Balay       /* initial separator set obtained from the symm. shared case  */
733388eb383SBarry Smith       SETERRQ(PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n");
734827bd09bSSatish Balay       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
735827bd09bSSatish Balay 	{
736827bd09bSSatish Balay 	  /* set rsh of hc, fire, and collect lhs responses */
737827bd09bSSatish Balay 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
738827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
739827bd09bSSatish Balay 
740827bd09bSSatish Balay 	  /* set lsh of hc, fire, and collect rhs responses */
741827bd09bSSatish Balay 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
742827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
743827bd09bSSatish Balay 
744827bd09bSSatish Balay 	  for (i=0;i<n;i++)
745827bd09bSSatish Balay 	    {
746827bd09bSSatish Balay 	      if (id< mask)
747827bd09bSSatish Balay 		{
748827bd09bSSatish Balay 		  if (lhs[i]!=0.0)
749827bd09bSSatish Balay 		    {lhs[i]=1.0;}
750827bd09bSSatish Balay 		}
751827bd09bSSatish Balay 	      if (id>=mask)
752827bd09bSSatish Balay 		{
753827bd09bSSatish Balay 		  if (rhs[i]!=0.0)
754827bd09bSSatish Balay 		    {rhs[i]=1.0;}
755827bd09bSSatish Balay 		}
756827bd09bSSatish Balay 	    }
757827bd09bSSatish Balay 
758827bd09bSSatish Balay 	  if (id< mask)
759827bd09bSSatish Balay 	    {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
760827bd09bSSatish Balay 	  else
761827bd09bSSatish Balay 	    {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}
762827bd09bSSatish Balay 
763827bd09bSSatish Balay 	  /* count number of dofs I own that have signal and not in sep set */
764827bd09bSSatish Balay 	  rvec_zero(rsum,4);
765827bd09bSSatish Balay 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
766827bd09bSSatish Balay 	    {
767827bd09bSSatish Balay 	      if (!used[i])
768827bd09bSSatish Balay 		{
769827bd09bSSatish Balay 		  /* number of unmarked dofs on node */
770827bd09bSSatish Balay 		  ct++;
771827bd09bSSatish Balay 		  /* number of dofs to be marked on lhs hc */
772827bd09bSSatish Balay 		  if (id< mask)
773827bd09bSSatish Balay 		    {
774827bd09bSSatish Balay 		      if (lhs[i]!=0.0)
775827bd09bSSatish Balay 			{sum[0]++; rsum[0]+=1.0/lhs[i];}
776827bd09bSSatish Balay 		    }
777827bd09bSSatish Balay 		  /* number of dofs to be marked on rhs hc */
778827bd09bSSatish Balay 		  if (id>=mask)
779827bd09bSSatish Balay 		    {
780827bd09bSSatish Balay 		      if (rhs[i]!=0.0)
781827bd09bSSatish Balay 			{sum[1]++; rsum[1]+=1.0/rhs[i];}
782827bd09bSSatish Balay 		    }
783827bd09bSSatish Balay 		}
784827bd09bSSatish Balay 	    }
785827bd09bSSatish Balay 
786827bd09bSSatish Balay 	  /* go for load balance - choose half with most unmarked dofs, bias LHS */
787827bd09bSSatish Balay 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
788827bd09bSSatish Balay 	  (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
789827bd09bSSatish Balay 	  giop_hc(sum,w,4,op,edge);
790827bd09bSSatish Balay 	  grop_hc(rsum,rw,4,op,edge);
791827bd09bSSatish Balay 	  rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
792827bd09bSSatish Balay 
793827bd09bSSatish Balay 	  if (id<mask)
794827bd09bSSatish Balay 	    {
795827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
796827bd09bSSatish Balay 	      for (ct=i=0;i<n;i++)
797827bd09bSSatish Balay 		{
798827bd09bSSatish Balay 		  if ((!used[i])&&(lhs[i]!=0.0))
799827bd09bSSatish Balay 		    {
800827bd09bSSatish Balay 		      ct++; nfo++;
801827bd09bSSatish Balay 
802827bd09bSSatish Balay 		      if (nfo>n)
803388eb383SBarry Smith 			{SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}
804827bd09bSSatish Balay 
805827bd09bSSatish Balay 		      *--iptr = local2global[i];
806827bd09bSSatish Balay 		      used[i]=edge;
807827bd09bSSatish Balay 		    }
808827bd09bSSatish Balay 		}
809827bd09bSSatish Balay 	      if (ct>1) {ivec_sort(iptr,ct);}
810827bd09bSSatish Balay 
811827bd09bSSatish Balay 	      lnsep[edge]=ct;
81252f87cdaSBarry Smith 	      nsep[edge]=(PetscInt) rsum[0];
813827bd09bSSatish Balay 	      dir [edge]=LEFT;
814827bd09bSSatish Balay 	    }
815827bd09bSSatish Balay 
816827bd09bSSatish Balay 	  if (id>=mask)
817827bd09bSSatish Balay 	    {
818827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
819827bd09bSSatish Balay 	      for (ct=i=0;i<n;i++)
820827bd09bSSatish Balay 		{
821827bd09bSSatish Balay 		  if ((!used[i])&&(rhs[i]!=0.0))
822827bd09bSSatish Balay 		    {
823827bd09bSSatish Balay 		      ct++; nfo++;
824827bd09bSSatish Balay 
825827bd09bSSatish Balay 		      if (nfo>n)
826388eb383SBarry Smith 			{SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}
827827bd09bSSatish Balay 
828827bd09bSSatish Balay 		      *--iptr = local2global[i];
829827bd09bSSatish Balay 		      used[i]=edge;
830827bd09bSSatish Balay 		    }
831827bd09bSSatish Balay 		}
832827bd09bSSatish Balay 	      if (ct>1) {ivec_sort(iptr,ct);}
833827bd09bSSatish Balay 
834827bd09bSSatish Balay 	      lnsep[edge]=ct;
83552f87cdaSBarry Smith 	      nsep[edge]= (PetscInt) rsum[1];
836827bd09bSSatish Balay 	      dir [edge]=RIGHT;
837827bd09bSSatish Balay 	    }
838827bd09bSSatish Balay 
839827bd09bSSatish Balay 	  /* LATER or we can recur on these to order seps at this level */
840827bd09bSSatish Balay 	  /* do we need full set of separators for this?                */
841827bd09bSSatish Balay 
842827bd09bSSatish Balay 	  /* fold rhs hc into lower */
843827bd09bSSatish Balay 	  if (id>=mask)
844827bd09bSSatish Balay 	    {id-=mask;}
845827bd09bSSatish Balay 	}
846827bd09bSSatish Balay     }
847827bd09bSSatish Balay   else
848827bd09bSSatish Balay     {
849827bd09bSSatish Balay       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
850827bd09bSSatish Balay 	{
851827bd09bSSatish Balay 	  /* set rsh of hc, fire, and collect lhs responses */
852827bd09bSSatish Balay 	  (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
853827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,lhs,"+\0",edge);
854827bd09bSSatish Balay 
855827bd09bSSatish Balay 	  /* set lsh of hc, fire, and collect rhs responses */
856827bd09bSSatish Balay 	  (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
857827bd09bSSatish Balay 	  gs_gop_hc(gs_handle,rhs,"+\0",edge);
858827bd09bSSatish Balay 
859827bd09bSSatish Balay 	  /* count number of dofs I own that have signal and not in sep set */
860827bd09bSSatish Balay 	  for (ivec_zero(sum,4),ct=i=0;i<n;i++)
861827bd09bSSatish Balay 	    {
862827bd09bSSatish Balay 	      if (!used[i])
863827bd09bSSatish Balay 		{
864827bd09bSSatish Balay 		  /* number of unmarked dofs on node */
865827bd09bSSatish Balay 		  ct++;
866827bd09bSSatish Balay 		  /* number of dofs to be marked on lhs hc */
867827bd09bSSatish Balay 		  if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
868827bd09bSSatish Balay 		  /* number of dofs to be marked on rhs hc */
869827bd09bSSatish Balay 		  if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
870827bd09bSSatish Balay 		}
871827bd09bSSatish Balay 	    }
872827bd09bSSatish Balay 
873827bd09bSSatish Balay 	  /* for the non-symmetric case we need separators of width 2 */
874827bd09bSSatish Balay 	  /* so take both sides */
875827bd09bSSatish Balay 	  (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
876827bd09bSSatish Balay 	  giop_hc(sum,w,4,op,edge);
877827bd09bSSatish Balay 
878827bd09bSSatish Balay 	  ct=0;
879827bd09bSSatish Balay 	  if (id<mask)
880827bd09bSSatish Balay 	    {
881827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
882827bd09bSSatish Balay 	      for (i=0;i<n;i++)
883827bd09bSSatish Balay 		{
884827bd09bSSatish Balay 		  if ((!used[i])&&(lhs[i]!=0.0))
885827bd09bSSatish Balay 		    {
886827bd09bSSatish Balay 		      ct++; nfo++;
887827bd09bSSatish Balay 		      *--iptr = local2global[i];
888827bd09bSSatish Balay 		      used[i]=edge;
889827bd09bSSatish Balay 		    }
890827bd09bSSatish Balay 		}
891827bd09bSSatish Balay 	      /* LSH hc summation of ct should be sum[0] */
892827bd09bSSatish Balay 	    }
893827bd09bSSatish Balay 	  else
894827bd09bSSatish Balay 	    {
895827bd09bSSatish Balay 	      /* mark dofs I own that have signal and not in sep set */
896827bd09bSSatish Balay 	      for (i=0;i<n;i++)
897827bd09bSSatish Balay 		{
898827bd09bSSatish Balay 		  if ((!used[i])&&(rhs[i]!=0.0))
899827bd09bSSatish Balay 		    {
900827bd09bSSatish Balay 		      ct++; nfo++;
901827bd09bSSatish Balay 		      *--iptr = local2global[i];
902827bd09bSSatish Balay 		      used[i]=edge;
903827bd09bSSatish Balay 		    }
904827bd09bSSatish Balay 		}
905827bd09bSSatish Balay 	      /* RSH hc summation of ct should be sum[1] */
906827bd09bSSatish Balay 	    }
907827bd09bSSatish Balay 
908827bd09bSSatish Balay 	  if (ct>1) {ivec_sort(iptr,ct);}
909827bd09bSSatish Balay 	  lnsep[edge]=ct;
910827bd09bSSatish Balay 	  nsep[edge]=sum[0]+sum[1];
911827bd09bSSatish Balay 	  dir [edge]=BOTH;
912827bd09bSSatish Balay 
913827bd09bSSatish Balay 	  /* LATER or we can recur on these to order seps at this level */
914827bd09bSSatish Balay 	  /* do we need full set of separators for this?                */
915827bd09bSSatish Balay 
916827bd09bSSatish Balay 	  /* fold rhs hc into lower */
917827bd09bSSatish Balay 	  if (id>=mask)
918827bd09bSSatish Balay 	    {id-=mask;}
919827bd09bSSatish Balay 	}
920827bd09bSSatish Balay     }
921827bd09bSSatish Balay 
922827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
923827bd09bSSatish Balay   for (ct=i=0;i<n;i++)
924827bd09bSSatish Balay     {
925827bd09bSSatish Balay       if (!used[i])
926827bd09bSSatish Balay 	{
927827bd09bSSatish Balay 	  ct++; nfo++;
928827bd09bSSatish Balay 	  *--iptr = local2global[i];
929827bd09bSSatish Balay 	  used[i]=edge;
930827bd09bSSatish Balay 	}
931827bd09bSSatish Balay     }
932827bd09bSSatish Balay   if (ct>1) {ivec_sort(iptr,ct);}
933827bd09bSSatish Balay   lnsep[edge]=ct;
934827bd09bSSatish Balay   nsep [edge]=ct;
935827bd09bSSatish Balay   dir  [edge]=BOTH;
936827bd09bSSatish Balay 
937827bd09bSSatish Balay   xyt_handle->info->nsep=nsep;
938827bd09bSSatish Balay   xyt_handle->info->lnsep=lnsep;
939827bd09bSSatish Balay   xyt_handle->info->fo=fo;
940827bd09bSSatish Balay   xyt_handle->info->nfo=nfo;
941827bd09bSSatish Balay 
942a501084fSBarry Smith   free(dir);
943a501084fSBarry Smith   free(lhs);
944a501084fSBarry Smith   free(rhs);
945a501084fSBarry Smith   free(used);
9460924e98cSBarry Smith   PetscFunctionReturn(0);
947827bd09bSSatish Balay }
948827bd09bSSatish Balay 
9497b1ae94cSBarry Smith /**************************************xyt.c***********************************/
95052f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
951827bd09bSSatish Balay {
952827bd09bSSatish Balay   mv_info *mvi;
953827bd09bSSatish Balay 
954827bd09bSSatish Balay 
955a501084fSBarry Smith   mvi = (mv_info*)malloc(sizeof(mv_info));
956827bd09bSSatish Balay   mvi->n=n;
957827bd09bSSatish Balay   mvi->m=m;
958827bd09bSSatish Balay   mvi->n_global=-1;
959827bd09bSSatish Balay   mvi->m_global=-1;
96052f87cdaSBarry Smith   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
961827bd09bSSatish Balay   ivec_copy(mvi->local2global,local2global,m);
962827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
963a501084fSBarry Smith   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
964827bd09bSSatish Balay   mvi->grid_data=grid_data;
965827bd09bSSatish Balay 
966827bd09bSSatish Balay   /* set xyt communication handle to perform restricted matvec */
967827bd09bSSatish Balay   mvi->gs_handle = gs_init(local2global, m, num_nodes);
968827bd09bSSatish Balay 
969827bd09bSSatish Balay   return(mvi);
970827bd09bSSatish Balay }
971827bd09bSSatish Balay 
9727b1ae94cSBarry Smith /**************************************xyt.c***********************************/
9733fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
974827bd09bSSatish Balay {
9750924e98cSBarry Smith   PetscFunctionBegin;
976827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
9770924e98cSBarry Smith   PetscFunctionReturn(0);
978827bd09bSSatish Balay }
979827bd09bSSatish Balay 
980827bd09bSSatish Balay 
981827bd09bSSatish Balay 
982