xref: /petsc/src/ksp/pc/impls/tfs/xyt.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /*************************************xyt.c************************************
3827bd09bSSatish Balay Module Name: xyt
4827bd09bSSatish Balay Module Info:
5827bd09bSSatish Balay 
6827bd09bSSatish Balay author:  Henry M. Tufo III
7827bd09bSSatish Balay e-mail:  hmt@asci.uchicago.edu
8827bd09bSSatish Balay contact:
9827bd09bSSatish Balay +--------------------------------+--------------------------------+
10827bd09bSSatish Balay |MCS Division - Building 221     |Department of Computer Science  |
11827bd09bSSatish Balay |Argonne National Laboratory     |Ryerson 152                     |
12827bd09bSSatish Balay |9700 S. Cass Avenue             |The University of Chicago       |
13827bd09bSSatish Balay |Argonne, IL  60439              |Chicago, IL  60637              |
14827bd09bSSatish Balay |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
15827bd09bSSatish Balay +--------------------------------+--------------------------------+
16827bd09bSSatish Balay 
17827bd09bSSatish Balay Last Modification: 3.20.01
18827bd09bSSatish Balay **************************************xyt.c***********************************/
19c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
20827bd09bSSatish Balay 
21827bd09bSSatish Balay #define LEFT  -1
22827bd09bSSatish Balay #define RIGHT  1
23827bd09bSSatish Balay #define BOTH   0
24827bd09bSSatish Balay 
25827bd09bSSatish Balay typedef struct xyt_solver_info {
2652f87cdaSBarry Smith   PetscInt    n, m, n_global, m_global;
2752f87cdaSBarry Smith   PetscInt    nnz, max_nnz, msg_buf_sz;
2852f87cdaSBarry Smith   PetscInt    *nsep, *lnsep, *fo, nfo, *stages;
2952f87cdaSBarry Smith   PetscInt    *xcol_sz, *xcol_indices;
30a501084fSBarry Smith   PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
3152f87cdaSBarry Smith   PetscInt    *ycol_sz, *ycol_indices;
32a501084fSBarry Smith   PetscScalar **ycol_vals, *y;
3352f87cdaSBarry Smith   PetscInt    nsolves;
34a501084fSBarry Smith   PetscScalar tot_solve_time;
35827bd09bSSatish Balay } xyt_info;
36827bd09bSSatish Balay 
37827bd09bSSatish Balay 
38827bd09bSSatish Balay typedef struct matvec_info {
3952f87cdaSBarry Smith   PetscInt     n, m, n_global, m_global;
4052f87cdaSBarry Smith   PetscInt     *local2global;
41ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle;
42a501084fSBarry Smith   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
43827bd09bSSatish Balay   void *grid_data;
44827bd09bSSatish Balay } mv_info;
45827bd09bSSatish Balay 
46827bd09bSSatish Balay struct xyt_CDT {
4752f87cdaSBarry Smith   PetscInt id;
4852f87cdaSBarry Smith   PetscInt ns;
4952f87cdaSBarry Smith   PetscInt level;
50827bd09bSSatish Balay   xyt_info *info;
51827bd09bSSatish Balay   mv_info  *mvi;
52827bd09bSSatish Balay };
53827bd09bSSatish Balay 
5452f87cdaSBarry Smith static PetscInt n_xyt        =0;
5552f87cdaSBarry Smith static PetscInt n_xyt_handles=0;
56827bd09bSSatish Balay 
57827bd09bSSatish Balay /* prototypes */
583fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
593fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle);
603fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle);
613fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
6252f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle);
6352f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle);
645c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data);
65827bd09bSSatish Balay 
667b1ae94cSBarry Smith /**************************************xyt.c***********************************/
677b1ae94cSBarry Smith xyt_ADT XYT_new(void)
68827bd09bSSatish Balay {
69827bd09bSSatish Balay   xyt_ADT xyt_handle;
70827bd09bSSatish Balay 
71827bd09bSSatish Balay   /* rolling count on n_xyt ... pot. problem here */
72827bd09bSSatish Balay   n_xyt_handles++;
73a501084fSBarry Smith   xyt_handle       = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
74827bd09bSSatish Balay   xyt_handle->id   = ++n_xyt;
75827bd09bSSatish Balay   xyt_handle->info = NULL;
76827bd09bSSatish Balay   xyt_handle->mvi  = NULL;
77827bd09bSSatish Balay 
78827bd09bSSatish Balay   return(xyt_handle);
79827bd09bSSatish Balay }
80827bd09bSSatish Balay 
817b1ae94cSBarry Smith /**************************************xyt.c***********************************/
8252f87cdaSBarry Smith PetscInt XYT_factor(xyt_ADT xyt_handle,     /* prev. allocated xyt  handle */
8352f87cdaSBarry Smith                     PetscInt *local2global, /* global column mapping       */
8452f87cdaSBarry Smith                     PetscInt n,             /* local num rows              */
8552f87cdaSBarry Smith                     PetscInt m,             /* local num cols              */
865c8f6a95SKarl Rupp                     PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc         */
871147fc2aSKarl Rupp                     void *grid_data)        /* grid data for matvec        */
88827bd09bSSatish Balay {
89827bd09bSSatish Balay 
90b1c944f5SJed Brown   PCTFS_comm_init();
91827bd09bSSatish Balay   check_handle(xyt_handle);
92827bd09bSSatish Balay 
93827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
94b1c944f5SJed Brown   if ((1<<(xyt_handle->level=PCTFS_i_log2_num_nodes))!=PCTFS_num_nodes) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %D != %D\n",1<<PCTFS_i_log2_num_nodes,PCTFS_num_nodes);
95827bd09bSSatish Balay 
96827bd09bSSatish Balay   /* space for X info */
97a501084fSBarry Smith   xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));
98827bd09bSSatish Balay 
99827bd09bSSatish Balay   /* set up matvec handles */
1005c8f6a95SKarl Rupp   xyt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data);
101827bd09bSSatish Balay 
102827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
103827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
104827bd09bSSatish Balay   xyt_handle->ns=0;
105827bd09bSSatish Balay 
106827bd09bSSatish Balay   /* determine separators and generate firing order - NB xyt info set here */
107827bd09bSSatish Balay   det_separators(xyt_handle);
108827bd09bSSatish Balay 
109827bd09bSSatish Balay   return(do_xyt_factor(xyt_handle));
110827bd09bSSatish Balay }
111827bd09bSSatish Balay 
1127b1ae94cSBarry Smith /**************************************xyt.c***********************************/
1138cda6cd7SBarry Smith PetscInt XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b)
114827bd09bSSatish Balay {
115b1c944f5SJed Brown   PCTFS_comm_init();
116827bd09bSSatish Balay   check_handle(xyt_handle);
117827bd09bSSatish Balay 
118827bd09bSSatish Balay   /* need to copy b into x? */
119*2fa5cd67SKarl Rupp   if (b) PCTFS_rvec_copy(x,b,xyt_handle->mvi->n);do_xyt_solve(xyt_handle,x);
120827bd09bSSatish Balay 
121827bd09bSSatish Balay   return(0);
122827bd09bSSatish Balay }
123827bd09bSSatish Balay 
1247b1ae94cSBarry Smith /**************************************xyt.c***********************************/
12552f87cdaSBarry Smith PetscInt XYT_free(xyt_ADT xyt_handle)
126827bd09bSSatish Balay {
127b1c944f5SJed Brown   PCTFS_comm_init();
128827bd09bSSatish Balay   check_handle(xyt_handle);
129827bd09bSSatish Balay   n_xyt_handles--;
130827bd09bSSatish Balay 
131a501084fSBarry Smith   free(xyt_handle->info->nsep);
132a501084fSBarry Smith   free(xyt_handle->info->lnsep);
133a501084fSBarry Smith   free(xyt_handle->info->fo);
134a501084fSBarry Smith   free(xyt_handle->info->stages);
135a501084fSBarry Smith   free(xyt_handle->info->solve_uu);
136a501084fSBarry Smith   free(xyt_handle->info->solve_w);
137a501084fSBarry Smith   free(xyt_handle->info->x);
138a501084fSBarry Smith   free(xyt_handle->info->xcol_vals);
139a501084fSBarry Smith   free(xyt_handle->info->xcol_sz);
140a501084fSBarry Smith   free(xyt_handle->info->xcol_indices);
141a501084fSBarry Smith   free(xyt_handle->info->y);
142a501084fSBarry Smith   free(xyt_handle->info->ycol_vals);
143a501084fSBarry Smith   free(xyt_handle->info->ycol_sz);
144a501084fSBarry Smith   free(xyt_handle->info->ycol_indices);
145a501084fSBarry Smith   free(xyt_handle->info);
146a501084fSBarry Smith   free(xyt_handle->mvi->local2global);
147ca8e9878SJed Brown   PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle);
148a501084fSBarry Smith   free(xyt_handle->mvi);
149a501084fSBarry Smith   free(xyt_handle);
150827bd09bSSatish Balay 
151827bd09bSSatish Balay 
152827bd09bSSatish Balay   /* if the check fails we nuke */
153a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
154827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
155827bd09bSSatish Balay   return(0);
156827bd09bSSatish Balay }
157827bd09bSSatish Balay 
1587b1ae94cSBarry Smith /**************************************xyt.c***********************************/
15952f87cdaSBarry Smith PetscInt XYT_stats(xyt_ADT xyt_handle)
160827bd09bSSatish Balay {
16152f87cdaSBarry Smith   PetscInt    op[]  = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
16252f87cdaSBarry Smith   PetscInt    fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
16352f87cdaSBarry Smith   PetscInt    vals[9],  work[9];
164a501084fSBarry Smith   PetscScalar fvals[3], fwork[3];
165827bd09bSSatish Balay 
166b1c944f5SJed Brown   PCTFS_comm_init();
167827bd09bSSatish Balay   check_handle(xyt_handle);
168827bd09bSSatish Balay 
169827bd09bSSatish Balay   /* if factorization not done there are no stats */
1707d0a6c19SBarry Smith   if (!xyt_handle->info||!xyt_handle->mvi) {
171b1c944f5SJed Brown     if (!PCTFS_my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n");
172827bd09bSSatish Balay     return 1;
173827bd09bSSatish Balay   }
174827bd09bSSatish Balay 
175827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
176827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
177827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
178b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
179827bd09bSSatish Balay 
180*2fa5cd67SKarl Rupp   fvals[0]=fvals[1]=fvals[2]=xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
181b1c944f5SJed Brown   PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
182827bd09bSSatish Balay 
183b1c944f5SJed Brown   if (!PCTFS_my_id) {
184b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_nnz=%D\n",PCTFS_my_id,vals[0]);
185b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_nnz=%D\n",PCTFS_my_id,vals[1]);
186b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);
187b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_nnz=%D\n",PCTFS_my_id,vals[2]);
188b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(2d)  =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.5)));
189b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(3d)  =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
190b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_n  =%D\n",PCTFS_my_id,vals[3]);
191b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_n  =%D\n",PCTFS_my_id,vals[4]);
192b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_n  =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);
193b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_n  =%D\n",PCTFS_my_id,vals[5]);
194b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_buf=%D\n",PCTFS_my_id,vals[6]);
195b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_buf=%D\n",PCTFS_my_id,vals[7]);
196b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);
197b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_slv=%g\n",PCTFS_my_id,fvals[0]);
198b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_slv=%g\n",PCTFS_my_id,fvals[1]);
199b1c944f5SJed Brown     PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);
200827bd09bSSatish Balay   }
201827bd09bSSatish Balay 
202827bd09bSSatish Balay   return(0);
203827bd09bSSatish Balay }
204827bd09bSSatish Balay 
205827bd09bSSatish Balay 
206827bd09bSSatish Balay /*************************************xyt.c************************************
207827bd09bSSatish Balay 
208827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
209827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
210827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
211827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
212827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
213827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
214ca8e9878SJed Brown    PCTFS_gs_init/gop).
215827bd09bSSatish Balay 
216827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
217827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
218827bd09bSSatish Balay **************************************xyt.c***********************************/
21952f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle)
220827bd09bSSatish Balay {
2217b1ae94cSBarry Smith   return xyt_generate(xyt_handle);
222827bd09bSSatish Balay }
223827bd09bSSatish Balay 
2247b1ae94cSBarry Smith /**************************************xyt.c***********************************/
22552f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle)
226827bd09bSSatish Balay {
22752f87cdaSBarry Smith   PetscInt       i,j,k,idx;
22852f87cdaSBarry Smith   PetscInt       dim, col;
229a501084fSBarry Smith   PetscScalar    *u, *uu, *v, *z, *w, alpha, alpha_w;
23052f87cdaSBarry Smith   PetscInt       *segs;
23152f87cdaSBarry Smith   PetscInt       op[] = {GL_ADD,0};
23252f87cdaSBarry Smith   PetscInt       off, len;
233a501084fSBarry Smith   PetscScalar    *x_ptr, *y_ptr;
23452f87cdaSBarry Smith   PetscInt       *iptr, flag;
23552f87cdaSBarry Smith   PetscInt       start =0, end, work;
23652f87cdaSBarry Smith   PetscInt       op2[] = {GL_MIN,0};
237ca8e9878SJed Brown   PCTFS_gs_ADT   PCTFS_gs_handle;
23852f87cdaSBarry Smith   PetscInt       *nsep, *lnsep, *fo;
23952f87cdaSBarry Smith   PetscInt       a_n            =xyt_handle->mvi->n;
24052f87cdaSBarry Smith   PetscInt       a_m            =xyt_handle->mvi->m;
24152f87cdaSBarry Smith   PetscInt       *a_local2global=xyt_handle->mvi->local2global;
24252f87cdaSBarry Smith   PetscInt       level;
24352f87cdaSBarry Smith   PetscInt       n, m;
24452f87cdaSBarry Smith   PetscInt       *xcol_sz, *xcol_indices, *stages;
245a501084fSBarry Smith   PetscScalar    **xcol_vals, *x;
24652f87cdaSBarry Smith   PetscInt       *ycol_sz, *ycol_indices;
247a501084fSBarry Smith   PetscScalar    **ycol_vals, *y;
24852f87cdaSBarry Smith   PetscInt       n_global;
24952f87cdaSBarry Smith   PetscInt       xt_nnz       =0, xt_max_nnz=0;
25052f87cdaSBarry Smith   PetscInt       yt_nnz       =0, yt_max_nnz=0;
25152f87cdaSBarry Smith   PetscInt       xt_zero_nnz  =0;
25252f87cdaSBarry Smith   PetscInt       xt_zero_nnz_0=0;
25352f87cdaSBarry Smith   PetscInt       yt_zero_nnz  =0;
25452f87cdaSBarry Smith   PetscInt       yt_zero_nnz_0=0;
2554a0de3f6SBarry Smith   PetscBLASInt   i1           = 1,dlen;
256a501084fSBarry Smith   PetscScalar    dm1          = -1.0;
257d1528f56SBarry Smith   PetscErrorCode ierr;
258827bd09bSSatish Balay 
259827bd09bSSatish Balay   n              =xyt_handle->mvi->n;
260827bd09bSSatish Balay   nsep           =xyt_handle->info->nsep;
261827bd09bSSatish Balay   lnsep          =xyt_handle->info->lnsep;
262827bd09bSSatish Balay   fo             =xyt_handle->info->fo;
263827bd09bSSatish Balay   end            =lnsep[0];
264827bd09bSSatish Balay   level          =xyt_handle->level;
265ca8e9878SJed Brown   PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle;
266827bd09bSSatish Balay 
267827bd09bSSatish Balay   /* is there a null space? */
268827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
269*2fa5cd67SKarl Rupp   for (i=0, j=0; i<=level; i++) j+=nsep[i];
270827bd09bSSatish Balay 
271827bd09bSSatish Balay   m = j-xyt_handle->ns;
27222d28d08SBarry Smith   if (m!=j) {
27322d28d08SBarry Smith     ierr = PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);CHKERRQ(ierr);
27422d28d08SBarry Smith   }
275827bd09bSSatish Balay 
2764a0de3f6SBarry Smith   ierr = PetscInfo2(0,"xyt_generate() :: X(%D,%D)\n",n,m);CHKERRQ(ierr);
277827bd09bSSatish Balay 
278827bd09bSSatish Balay   /* get and initialize storage for x local         */
279827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
28052f87cdaSBarry Smith   xcol_sz      = (PetscInt*) malloc(m*sizeof(PetscInt));
28152f87cdaSBarry Smith   xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
282a501084fSBarry Smith   xcol_vals    = (PetscScalar**) malloc(m*sizeof(PetscScalar*));
283db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
284827bd09bSSatish Balay     xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
285827bd09bSSatish Balay     xcol_vals[i]   = NULL;
286827bd09bSSatish Balay   }
287827bd09bSSatish Balay   xcol_indices[j]=-1;
288827bd09bSSatish Balay 
289827bd09bSSatish Balay   /* get and initialize storage for y local         */
290827bd09bSSatish Balay   /* note that y local is nxm and stored by columns */
29152f87cdaSBarry Smith   ycol_sz      = (PetscInt*) malloc(m*sizeof(PetscInt));
29252f87cdaSBarry Smith   ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
293a501084fSBarry Smith   ycol_vals    = (PetscScalar**) malloc(m*sizeof(PetscScalar*));
294db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
295827bd09bSSatish Balay     ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
296827bd09bSSatish Balay     ycol_vals[i]   = NULL;
297827bd09bSSatish Balay   }
298827bd09bSSatish Balay   ycol_indices[j]=-1;
299827bd09bSSatish Balay 
300827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
301827bd09bSSatish Balay   /* this looks like nsep[]=segments */
30252f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
30352f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
304ca8e9878SJed Brown   PCTFS_ivec_zero(stages,level+1);
305ca8e9878SJed Brown   PCTFS_ivec_copy(segs,nsep,level+1);
306*2fa5cd67SKarl Rupp   for (i=0; i<level; i++) segs[i+1] += segs[i];
307827bd09bSSatish Balay   stages[0] = segs[0];
308827bd09bSSatish Balay 
309827bd09bSSatish Balay   /* temporary vectors  */
310a501084fSBarry Smith   u  = (PetscScalar*) malloc(n*sizeof(PetscScalar));
311a501084fSBarry Smith   z  = (PetscScalar*) malloc(n*sizeof(PetscScalar));
312a501084fSBarry Smith   v  = (PetscScalar*) malloc(a_m*sizeof(PetscScalar));
313a501084fSBarry Smith   uu = (PetscScalar*) malloc(m*sizeof(PetscScalar));
314a501084fSBarry Smith   w  = (PetscScalar*) malloc(m*sizeof(PetscScalar));
315827bd09bSSatish Balay 
316827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
317*2fa5cd67SKarl Rupp   for (i=1, j=0; i<=level; i++) j+=nsep[i];
318827bd09bSSatish Balay 
319827bd09bSSatish Balay   /* storage for sparse x values */
320827bd09bSSatish Balay   n_global   = xyt_handle->info->n_global;
321b1c944f5SJed Brown   xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes;
322a501084fSBarry Smith   x          = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar));
323a501084fSBarry Smith   y          = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar));
324827bd09bSSatish Balay 
325827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
326827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
327db4deed7SKarl Rupp   for (dim=i=j=0; i<m; i++) {
328827bd09bSSatish Balay     /* time to move to the next level? */
329d4af0d30SBarry Smith     while (i==segs[dim]) {
330e32f2f54SBarry Smith       if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n");
331827bd09bSSatish Balay       stages[dim++]=i;
332827bd09bSSatish Balay       end         +=lnsep[dim];
333827bd09bSSatish Balay     }
334827bd09bSSatish Balay     stages[dim]=i;
335827bd09bSSatish Balay 
336827bd09bSSatish Balay     /* which column are we firing? */
337827bd09bSSatish Balay     /* i.e. set v_l */
338827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
339827bd09bSSatish Balay     (start<end) ? (col=fo[start]) : (col=INT_MAX);
340b1c944f5SJed Brown     PCTFS_giop_hc(&col,&work,1,op2,dim);
341827bd09bSSatish Balay 
342827bd09bSSatish Balay     /* shouldn't need this */
343db4deed7SKarl Rupp     if (col==INT_MAX) {
344f1ed62a8SBarry Smith       ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
345827bd09bSSatish Balay       continue;
346827bd09bSSatish Balay     }
347827bd09bSSatish Balay 
348827bd09bSSatish Balay     /* do I own it? I should */
349ca8e9878SJed Brown     PCTFS_rvec_zero(v,a_m);
350*2fa5cd67SKarl Rupp     if (col==fo[start]) {
351827bd09bSSatish Balay       start++;
352ca8e9878SJed Brown       idx=PCTFS_ivec_linear_search(col, a_local2global, a_n);
353*2fa5cd67SKarl Rupp       if (idx!=-1) {
354*2fa5cd67SKarl Rupp         v[idx] = 1.0;
355*2fa5cd67SKarl Rupp         j++;
356*2fa5cd67SKarl Rupp       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");
357db4deed7SKarl Rupp     } else {
358ca8e9878SJed Brown       idx=PCTFS_ivec_linear_search(col, a_local2global, a_m);
359*2fa5cd67SKarl Rupp       if (idx!=-1) v[idx] = 1.0;
360827bd09bSSatish Balay     }
361827bd09bSSatish Balay 
362827bd09bSSatish Balay     /* perform u = A.v_l */
363ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
364827bd09bSSatish Balay     do_matvec(xyt_handle->mvi,v,u);
365827bd09bSSatish Balay 
366827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
367827bd09bSSatish Balay     /* technically only need to zero out first i entries */
368827bd09bSSatish Balay     /* later turn this into an XYT_solve call ? */
369ca8e9878SJed Brown     PCTFS_rvec_zero(uu,m);
370827bd09bSSatish Balay     y_ptr=y;
371827bd09bSSatish Balay     iptr = ycol_indices;
372db4deed7SKarl Rupp     for (k=0; k<i; k++) {
373827bd09bSSatish Balay       off   = *iptr++;
374827bd09bSSatish Balay       len   = *iptr++;
375c5df96a5SBarry Smith       ierr  = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
3764a0de3f6SBarry Smith       uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1);
377827bd09bSSatish Balay       y_ptr+=len;
378827bd09bSSatish Balay     }
379827bd09bSSatish Balay 
380827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
381b1c944f5SJed Brown     PCTFS_ssgl_radd  (uu, w, dim, stages);
382827bd09bSSatish Balay 
383827bd09bSSatish Balay     /* z = X.uu */
384ca8e9878SJed Brown     PCTFS_rvec_zero(z,n);
385827bd09bSSatish Balay     x_ptr=x;
386827bd09bSSatish Balay     iptr = xcol_indices;
387db4deed7SKarl Rupp     for (k=0; k<i; k++) {
388827bd09bSSatish Balay       off  = *iptr++;
389827bd09bSSatish Balay       len  = *iptr++;
390c5df96a5SBarry Smith       ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
3914a0de3f6SBarry Smith       BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
392827bd09bSSatish Balay       x_ptr+=len;
393827bd09bSSatish Balay     }
394827bd09bSSatish Balay 
395827bd09bSSatish Balay     /* compute v_l = v_l - z */
396ca8e9878SJed Brown     PCTFS_rvec_zero(v+a_n,a_m-a_n);
397c5df96a5SBarry Smith     ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr);
3984a0de3f6SBarry Smith     BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);
399827bd09bSSatish Balay 
400827bd09bSSatish Balay     /* compute u_l = A.v_l */
401*2fa5cd67SKarl Rupp     if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim);
402ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
403827bd09bSSatish Balay     do_matvec(xyt_handle->mvi,v,u);
404827bd09bSSatish Balay 
405827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
406c5df96a5SBarry Smith     ierr  = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr);
4074a0de3f6SBarry Smith     alpha = BLASdot_(&dlen,u,&i1,u,&i1);
408827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
409b1c944f5SJed Brown     PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
410827bd09bSSatish Balay 
4118f1a2a5eSBarry Smith     alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha);
412827bd09bSSatish Balay 
413827bd09bSSatish Balay     /* check for small alpha                             */
414827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
415e32f2f54SBarry Smith     if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);
416827bd09bSSatish Balay 
417827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
418ca8e9878SJed Brown     PCTFS_rvec_scale(v,1.0/alpha,n);
419ca8e9878SJed Brown     PCTFS_rvec_scale(u,1.0/alpha,n);
420827bd09bSSatish Balay 
421827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
422827bd09bSSatish Balay     flag = 1;
423827bd09bSSatish Balay     off  =len=0;
424db4deed7SKarl Rupp     for (k=0; k<n; k++) {
425db4deed7SKarl Rupp       if (v[k]!=0.0) {
426827bd09bSSatish Balay         len=k;
427db4deed7SKarl Rupp         if (flag) {off=k; flag=0;}
428827bd09bSSatish Balay       }
429827bd09bSSatish Balay     }
430827bd09bSSatish Balay 
431827bd09bSSatish Balay     len -= (off-1);
432827bd09bSSatish Balay 
433db4deed7SKarl Rupp     if (len>0) {
434db4deed7SKarl Rupp       if ((xt_nnz+len)>xt_max_nnz) {
435f1ed62a8SBarry Smith         ierr        = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
436827bd09bSSatish Balay         xt_max_nnz *= 2;
437a501084fSBarry Smith         x_ptr       = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar));
438ca8e9878SJed Brown         PCTFS_rvec_copy(x_ptr,x,xt_nnz);
439a501084fSBarry Smith         free(x);
440827bd09bSSatish Balay         x     = x_ptr;
441827bd09bSSatish Balay         x_ptr+=xt_nnz;
442827bd09bSSatish Balay       }
443827bd09bSSatish Balay       xt_nnz += len;
444ca8e9878SJed Brown       PCTFS_rvec_copy(x_ptr,v+off,len);
445827bd09bSSatish Balay 
446827bd09bSSatish Balay       /* keep track of number of zeros */
447db4deed7SKarl Rupp       if (dim) {
448db4deed7SKarl Rupp         for (k=0; k<len; k++) {
449*2fa5cd67SKarl Rupp           if (x_ptr[k]==0.0) xt_zero_nnz++;
450827bd09bSSatish Balay         }
451db4deed7SKarl Rupp       } else {
452db4deed7SKarl Rupp         for (k=0; k<len; k++) {
453*2fa5cd67SKarl Rupp           if (x_ptr[k]==0.0) xt_zero_nnz_0++;
454827bd09bSSatish Balay         }
455827bd09bSSatish Balay       }
456827bd09bSSatish Balay       xcol_indices[2*i] = off;
457827bd09bSSatish Balay       xcol_sz[i]        = xcol_indices[2*i+1] = len;
458827bd09bSSatish Balay       xcol_vals[i]      = x_ptr;
459db4deed7SKarl Rupp     } else {
460827bd09bSSatish Balay       xcol_indices[2*i] = 0;
461827bd09bSSatish Balay       xcol_sz[i]        = xcol_indices[2*i+1] = 0;
462827bd09bSSatish Balay       xcol_vals[i]      = x_ptr;
463827bd09bSSatish Balay     }
464827bd09bSSatish Balay 
465827bd09bSSatish Balay 
466827bd09bSSatish Balay     /* add newly generated column, u_l, to Y */
467827bd09bSSatish Balay     flag = 1;
468827bd09bSSatish Balay     off  =len=0;
469db4deed7SKarl Rupp     for (k=0; k<n; k++) {
470db4deed7SKarl Rupp       if (u[k]!=0.0) {
471827bd09bSSatish Balay         len=k;
472db4deed7SKarl Rupp         if (flag) { off=k; flag=0; }
473827bd09bSSatish Balay       }
474827bd09bSSatish Balay     }
475827bd09bSSatish Balay 
476827bd09bSSatish Balay     len -= (off-1);
477827bd09bSSatish Balay 
478db4deed7SKarl Rupp     if (len>0) {
479db4deed7SKarl Rupp       if ((yt_nnz+len)>yt_max_nnz) {
480f1ed62a8SBarry Smith         ierr        = PetscInfo(0,"increasing space for Y by 2x!\n");CHKERRQ(ierr);
481827bd09bSSatish Balay         yt_max_nnz *= 2;
482a501084fSBarry Smith         y_ptr       = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar));
483ca8e9878SJed Brown         PCTFS_rvec_copy(y_ptr,y,yt_nnz);
484a501084fSBarry Smith         free(y);
485827bd09bSSatish Balay         y     = y_ptr;
486827bd09bSSatish Balay         y_ptr+=yt_nnz;
487827bd09bSSatish Balay       }
488827bd09bSSatish Balay       yt_nnz += len;
489ca8e9878SJed Brown       PCTFS_rvec_copy(y_ptr,u+off,len);
490827bd09bSSatish Balay 
491827bd09bSSatish Balay       /* keep track of number of zeros */
492db4deed7SKarl Rupp       if (dim) {
493db4deed7SKarl Rupp         for (k=0; k<len; k++) {
494*2fa5cd67SKarl Rupp           if (y_ptr[k]==0.0) yt_zero_nnz++;
495827bd09bSSatish Balay         }
496db4deed7SKarl Rupp       } else {
497db4deed7SKarl Rupp         for (k=0; k<len; k++) {
498*2fa5cd67SKarl Rupp           if (y_ptr[k]==0.0) yt_zero_nnz_0++;
499827bd09bSSatish Balay         }
500827bd09bSSatish Balay       }
501827bd09bSSatish Balay       ycol_indices[2*i] = off;
502827bd09bSSatish Balay       ycol_sz[i]        = ycol_indices[2*i+1] = len;
503827bd09bSSatish Balay       ycol_vals[i]      = y_ptr;
504db4deed7SKarl Rupp     } else {
505827bd09bSSatish Balay       ycol_indices[2*i] = 0;
506827bd09bSSatish Balay       ycol_sz[i]        = ycol_indices[2*i+1] = 0;
507827bd09bSSatish Balay       ycol_vals[i]      = y_ptr;
508827bd09bSSatish Balay     }
509827bd09bSSatish Balay   }
510827bd09bSSatish Balay 
511827bd09bSSatish Balay   /* close off stages for execution phase */
512db4deed7SKarl Rupp   while (dim!=level) {
513827bd09bSSatish Balay     stages[dim++]=i;
5144a0de3f6SBarry Smith     ierr         = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
515827bd09bSSatish Balay   }
516827bd09bSSatish Balay   stages[dim]=i;
517827bd09bSSatish Balay 
518827bd09bSSatish Balay   xyt_handle->info->n           =xyt_handle->mvi->n;
519827bd09bSSatish Balay   xyt_handle->info->m           =m;
520827bd09bSSatish Balay   xyt_handle->info->nnz         =xt_nnz + yt_nnz;
521827bd09bSSatish Balay   xyt_handle->info->max_nnz     =xt_max_nnz + yt_max_nnz;
522827bd09bSSatish Balay   xyt_handle->info->msg_buf_sz  =stages[level]-stages[0];
523a501084fSBarry Smith   xyt_handle->info->solve_uu    = (PetscScalar*) malloc(m*sizeof(PetscScalar));
524a501084fSBarry Smith   xyt_handle->info->solve_w     = (PetscScalar*) malloc(m*sizeof(PetscScalar));
525827bd09bSSatish Balay   xyt_handle->info->x           =x;
526827bd09bSSatish Balay   xyt_handle->info->xcol_vals   =xcol_vals;
527827bd09bSSatish Balay   xyt_handle->info->xcol_sz     =xcol_sz;
528827bd09bSSatish Balay   xyt_handle->info->xcol_indices=xcol_indices;
529827bd09bSSatish Balay   xyt_handle->info->stages      =stages;
530827bd09bSSatish Balay   xyt_handle->info->y           =y;
531827bd09bSSatish Balay   xyt_handle->info->ycol_vals   =ycol_vals;
532827bd09bSSatish Balay   xyt_handle->info->ycol_sz     =ycol_sz;
533827bd09bSSatish Balay   xyt_handle->info->ycol_indices=ycol_indices;
534827bd09bSSatish Balay 
535a501084fSBarry Smith   free(segs);
536a501084fSBarry Smith   free(u);
537a501084fSBarry Smith   free(v);
538a501084fSBarry Smith   free(uu);
539a501084fSBarry Smith   free(z);
540a501084fSBarry Smith   free(w);
541827bd09bSSatish Balay 
542827bd09bSSatish Balay   return(0);
543827bd09bSSatish Balay }
544827bd09bSSatish Balay 
5457b1ae94cSBarry Smith /**************************************xyt.c***********************************/
5460924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle,  PetscScalar *uc)
547827bd09bSSatish Balay {
54852f87cdaSBarry Smith   PetscInt       off, len, *iptr;
54952f87cdaSBarry Smith   PetscInt       level        =xyt_handle->level;
55052f87cdaSBarry Smith   PetscInt       n            =xyt_handle->info->n;
55152f87cdaSBarry Smith   PetscInt       m            =xyt_handle->info->m;
55252f87cdaSBarry Smith   PetscInt       *stages      =xyt_handle->info->stages;
55352f87cdaSBarry Smith   PetscInt       *xcol_indices=xyt_handle->info->xcol_indices;
55452f87cdaSBarry Smith   PetscInt       *ycol_indices=xyt_handle->info->ycol_indices;
555a501084fSBarry Smith   PetscScalar    *x_ptr, *y_ptr, *uu_ptr;
556a501084fSBarry Smith   PetscScalar    *solve_uu=xyt_handle->info->solve_uu;
557a501084fSBarry Smith   PetscScalar    *solve_w =xyt_handle->info->solve_w;
558a501084fSBarry Smith   PetscScalar    *x       =xyt_handle->info->x;
559a501084fSBarry Smith   PetscScalar    *y       =xyt_handle->info->y;
5604a0de3f6SBarry Smith   PetscBLASInt   i1       = 1,dlen;
561c5df96a5SBarry Smith   PetscErrorCode ierr;
562827bd09bSSatish Balay 
5630924e98cSBarry Smith   PetscFunctionBegin;
564827bd09bSSatish Balay   uu_ptr=solve_uu;
565ca8e9878SJed Brown   PCTFS_rvec_zero(uu_ptr,m);
566827bd09bSSatish Balay 
567827bd09bSSatish Balay   /* x  = X.Y^T.b */
568827bd09bSSatish Balay   /* uu = Y^T.b */
569827bd09bSSatish Balay   for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len)
570827bd09bSSatish Balay   {
571c5df96a5SBarry Smith     off       =*iptr++;
572c5df96a5SBarry Smith     len       =*iptr++;
573c5df96a5SBarry Smith     ierr      = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
5744a0de3f6SBarry Smith     *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1);
575827bd09bSSatish Balay   }
576827bd09bSSatish Balay 
577827bd09bSSatish Balay   /* comunication of beta */
578827bd09bSSatish Balay   uu_ptr=solve_uu;
579*2fa5cd67SKarl Rupp   if (level) PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages);
580ca8e9878SJed Brown   PCTFS_rvec_zero(uc,n);
581827bd09bSSatish Balay 
582827bd09bSSatish Balay   /* x = X.uu */
583db4deed7SKarl Rupp   for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) {
584c5df96a5SBarry Smith     off  =*iptr++;
585c5df96a5SBarry Smith     len  =*iptr++;
586c5df96a5SBarry Smith     ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
5874a0de3f6SBarry Smith     BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
588827bd09bSSatish Balay   }
5890924e98cSBarry Smith   PetscFunctionReturn(0);
590827bd09bSSatish Balay }
591827bd09bSSatish Balay 
5927b1ae94cSBarry Smith /**************************************xyt.c***********************************/
5930924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle)
594827bd09bSSatish Balay {
59552f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
596827bd09bSSatish Balay 
5970924e98cSBarry Smith   PetscFunctionBegin;
598e32f2f54SBarry Smith   if (!xyt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle);
599827bd09bSSatish Balay 
600827bd09bSSatish Balay   vals[0]=vals[1]=xyt_handle->id;
601b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
602e32f2f54SBarry Smith   if ((vals[0]!=vals[1])||(xyt_handle->id<=0)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %D/%D %D\n", vals[0],vals[1], xyt_handle->id);
6030924e98cSBarry Smith   PetscFunctionReturn(0);
604827bd09bSSatish Balay }
605827bd09bSSatish Balay 
6067b1ae94cSBarry Smith /**************************************xyt.c***********************************/
6070924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle)
608827bd09bSSatish Balay {
60952f87cdaSBarry Smith   PetscInt       i, ct, id;
61052f87cdaSBarry Smith   PetscInt       mask, edge, *iptr;
61152f87cdaSBarry Smith   PetscInt       *dir, *used;
61252f87cdaSBarry Smith   PetscInt       sum[4], w[4];
613a501084fSBarry Smith   PetscScalar    rsum[4], rw[4];
61452f87cdaSBarry Smith   PetscInt       op[] = {GL_ADD,0};
615a501084fSBarry Smith   PetscScalar    *lhs, *rhs;
61652f87cdaSBarry Smith   PetscInt       *nsep, *lnsep, *fo, nfo=0;
617ca8e9878SJed Brown   PCTFS_gs_ADT   PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle;
61852f87cdaSBarry Smith   PetscInt       *local2global  =xyt_handle->mvi->local2global;
61952f87cdaSBarry Smith   PetscInt       n              =xyt_handle->mvi->n;
62052f87cdaSBarry Smith   PetscInt       m              =xyt_handle->mvi->m;
62152f87cdaSBarry Smith   PetscInt       level          =xyt_handle->level;
622ab824b78SBarry Smith   PetscInt       shared         =0;
623d1528f56SBarry Smith   PetscErrorCode ierr;
624827bd09bSSatish Balay 
6250924e98cSBarry Smith   PetscFunctionBegin;
62652f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
62752f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
62852f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
62952f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
63052f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
631827bd09bSSatish Balay 
632ca8e9878SJed Brown   PCTFS_ivec_zero(dir,level+1);
633ca8e9878SJed Brown   PCTFS_ivec_zero(nsep,level+1);
634ca8e9878SJed Brown   PCTFS_ivec_zero(lnsep,level+1);
635ca8e9878SJed Brown   PCTFS_ivec_set (fo,-1,n+1);
636ca8e9878SJed Brown   PCTFS_ivec_zero(used,n);
637827bd09bSSatish Balay 
6388cda6cd7SBarry Smith   lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
6398cda6cd7SBarry Smith   rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
640827bd09bSSatish Balay 
641827bd09bSSatish Balay   /* determine the # of unique dof */
642ca8e9878SJed Brown   PCTFS_rvec_zero(lhs,m);
643ca8e9878SJed Brown   PCTFS_rvec_set(lhs,1.0,n);
644ca8e9878SJed Brown   PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level);
645ca8e9878SJed Brown   ierr = PetscInfo(0,"done first PCTFS_gs_gop_hc\n");CHKERRQ(ierr);
646ca8e9878SJed Brown   PCTFS_rvec_zero(rsum,2);
647827bd09bSSatish Balay   for (ct=i=0; i<n; i++)
648827bd09bSSatish Balay   {
649db4deed7SKarl Rupp     if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; }
650827bd09bSSatish Balay 
651*2fa5cd67SKarl Rupp     if (lhs[i]!=1.0) shared=1;
652827bd09bSSatish Balay   }
653827bd09bSSatish Balay 
654b1c944f5SJed Brown   PCTFS_grop_hc(rsum,rw,2,op,level);
655827bd09bSSatish Balay   rsum[0]+=0.1;
656827bd09bSSatish Balay   rsum[1]+=0.1;
657827bd09bSSatish Balay 
65852f87cdaSBarry Smith   xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0];
65952f87cdaSBarry Smith   xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0];
660827bd09bSSatish Balay 
661827bd09bSSatish Balay   /* determine separator sets top down */
662*2fa5cd67SKarl Rupp   if (shared) {
663827bd09bSSatish Balay     /* solution is to do as in the symmetric shared case but then */
664827bd09bSSatish Balay     /* pick the sub-hc with the most free dofs and do a mat-vec   */
665827bd09bSSatish Balay     /* and pick up the responses on the other sub-hc from the     */
666827bd09bSSatish Balay     /* initial separator set obtained from the symm. shared case  */
667e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n");
668db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level; edge>0; edge--,mask>>=1) {
669db4deed7SKarl Rupp 
670827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
671ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
672ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
673827bd09bSSatish Balay 
674827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
675ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
676ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
677827bd09bSSatish Balay 
678db4deed7SKarl Rupp       for (i=0; i<n; i++) {
679db4deed7SKarl Rupp         if (id< mask) {
680*2fa5cd67SKarl Rupp           if (lhs[i]!=0.0) lhs[i]=1.0;
681827bd09bSSatish Balay         }
682db4deed7SKarl Rupp         if (id>=mask) {
683*2fa5cd67SKarl Rupp           if (rhs[i]!=0.0) rhs[i]=1.0;
684827bd09bSSatish Balay         }
685827bd09bSSatish Balay       }
686827bd09bSSatish Balay 
687*2fa5cd67SKarl Rupp       if (id< mask) PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1);
688*2fa5cd67SKarl Rupp       else PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1);
689827bd09bSSatish Balay 
690827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
691ca8e9878SJed Brown       PCTFS_rvec_zero(rsum,4);
692db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
693db4deed7SKarl Rupp         if (!used[i]) {
694827bd09bSSatish Balay           /* number of unmarked dofs on node */
695827bd09bSSatish Balay           ct++;
696827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
697db4deed7SKarl Rupp           if (id< mask) {
698db4deed7SKarl Rupp             if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; }
699827bd09bSSatish Balay           }
700827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
701db4deed7SKarl Rupp           if (id>=mask) {
702db4deed7SKarl Rupp             if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; }
703827bd09bSSatish Balay           }
704827bd09bSSatish Balay         }
705827bd09bSSatish Balay       }
706827bd09bSSatish Balay 
707827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
708827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
709827bd09bSSatish Balay       (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
710b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
711b1c944f5SJed Brown       PCTFS_grop_hc(rsum,rw,4,op,edge);
712827bd09bSSatish Balay       rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
713827bd09bSSatish Balay 
714db4deed7SKarl Rupp       if (id<mask) {
715827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
716db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
717db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
718827bd09bSSatish Balay             ct++; nfo++;
719827bd09bSSatish Balay 
720c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
721827bd09bSSatish Balay 
722827bd09bSSatish Balay             *--iptr = local2global[i];
723827bd09bSSatish Balay             used[i]=edge;
724827bd09bSSatish Balay           }
725827bd09bSSatish Balay         }
726*2fa5cd67SKarl Rupp         if (ct>1) PCTFS_ivec_sort(iptr,ct);
727827bd09bSSatish Balay 
728827bd09bSSatish Balay         lnsep[edge]=ct;
72952f87cdaSBarry Smith         nsep[edge] =(PetscInt) rsum[0];
730827bd09bSSatish Balay         dir [edge] =LEFT;
731827bd09bSSatish Balay       }
732827bd09bSSatish Balay 
733db4deed7SKarl Rupp       if (id>=mask) {
734827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
735db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
736db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
737827bd09bSSatish Balay             ct++; nfo++;
738827bd09bSSatish Balay 
739c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
740827bd09bSSatish Balay 
741827bd09bSSatish Balay             *--iptr = local2global[i];
742827bd09bSSatish Balay             used[i] =edge;
743827bd09bSSatish Balay           }
744827bd09bSSatish Balay         }
745*2fa5cd67SKarl Rupp         if (ct>1) PCTFS_ivec_sort(iptr,ct);
746827bd09bSSatish Balay 
747827bd09bSSatish Balay         lnsep[edge] = ct;
74852f87cdaSBarry Smith         nsep[edge]  = (PetscInt) rsum[1];
749827bd09bSSatish Balay         dir [edge]  = RIGHT;
750827bd09bSSatish Balay       }
751827bd09bSSatish Balay 
752827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
753827bd09bSSatish Balay       /* do we need full set of separators for this?                */
754827bd09bSSatish Balay 
755827bd09bSSatish Balay       /* fold rhs hc into lower */
756*2fa5cd67SKarl Rupp       if (id>=mask) id-=mask;
757827bd09bSSatish Balay     }
758db4deed7SKarl Rupp   } else {
759db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
760827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
761ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
762ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
763827bd09bSSatish Balay 
764827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
765ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
766ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
767827bd09bSSatish Balay 
768827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
769db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
770db4deed7SKarl Rupp         if (!used[i]) {
771827bd09bSSatish Balay           /* number of unmarked dofs on node */
772827bd09bSSatish Balay           ct++;
773827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
774*2fa5cd67SKarl Rupp           if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++;
775827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
776*2fa5cd67SKarl Rupp           if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++;
777827bd09bSSatish Balay         }
778827bd09bSSatish Balay       }
779827bd09bSSatish Balay 
780827bd09bSSatish Balay       /* for the non-symmetric case we need separators of width 2 */
781827bd09bSSatish Balay       /* so take both sides */
782827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
783b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
784827bd09bSSatish Balay 
785827bd09bSSatish Balay       ct=0;
786db4deed7SKarl Rupp       if (id<mask) {
787827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
788db4deed7SKarl Rupp         for (i=0;i<n;i++) {
789db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
790827bd09bSSatish Balay             ct++; nfo++;
791827bd09bSSatish Balay             *--iptr = local2global[i];
792827bd09bSSatish Balay             used[i] =edge;
793827bd09bSSatish Balay           }
794827bd09bSSatish Balay         }
795827bd09bSSatish Balay         /* LSH hc summation of ct should be sum[0] */
796db4deed7SKarl Rupp       } else {
797827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
798db4deed7SKarl Rupp         for (i=0; i<n; i++) {
799db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
800827bd09bSSatish Balay             ct++; nfo++;
801827bd09bSSatish Balay             *--iptr = local2global[i];
802827bd09bSSatish Balay             used[i] = edge;
803827bd09bSSatish Balay           }
804827bd09bSSatish Balay         }
805827bd09bSSatish Balay         /* RSH hc summation of ct should be sum[1] */
806827bd09bSSatish Balay       }
807827bd09bSSatish Balay 
808*2fa5cd67SKarl Rupp       if (ct>1) PCTFS_ivec_sort(iptr,ct);
809827bd09bSSatish Balay       lnsep[edge]=ct;
810827bd09bSSatish Balay       nsep[edge] =sum[0]+sum[1];
811827bd09bSSatish Balay       dir [edge] =BOTH;
812827bd09bSSatish Balay 
813827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
814827bd09bSSatish Balay       /* do we need full set of separators for this?                */
815827bd09bSSatish Balay 
816827bd09bSSatish Balay       /* fold rhs hc into lower */
817*2fa5cd67SKarl Rupp       if (id>=mask) id-=mask;
818827bd09bSSatish Balay     }
819827bd09bSSatish Balay   }
820827bd09bSSatish Balay 
821827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
822db4deed7SKarl Rupp   for (ct=i=0;i<n;i++) {
823db4deed7SKarl Rupp     if (!used[i]) {
824827bd09bSSatish Balay       ct++; nfo++;
825827bd09bSSatish Balay       *--iptr = local2global[i];
826827bd09bSSatish Balay       used[i] = edge;
827827bd09bSSatish Balay     }
828827bd09bSSatish Balay   }
829*2fa5cd67SKarl Rupp   if (ct>1) PCTFS_ivec_sort(iptr,ct);
830827bd09bSSatish Balay   lnsep[edge]=ct;
831827bd09bSSatish Balay   nsep [edge]=ct;
832827bd09bSSatish Balay   dir  [edge]=BOTH;
833827bd09bSSatish Balay 
834827bd09bSSatish Balay   xyt_handle->info->nsep  = nsep;
835827bd09bSSatish Balay   xyt_handle->info->lnsep = lnsep;
836827bd09bSSatish Balay   xyt_handle->info->fo    = fo;
837827bd09bSSatish Balay   xyt_handle->info->nfo   = nfo;
838827bd09bSSatish Balay 
839a501084fSBarry Smith   free(dir);
840a501084fSBarry Smith   free(lhs);
841a501084fSBarry Smith   free(rhs);
842a501084fSBarry Smith   free(used);
8430924e98cSBarry Smith   PetscFunctionReturn(0);
844827bd09bSSatish Balay }
845827bd09bSSatish Balay 
8467b1ae94cSBarry Smith /**************************************xyt.c***********************************/
8475c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data)
848827bd09bSSatish Balay {
849827bd09bSSatish Balay   mv_info *mvi;
850827bd09bSSatish Balay 
851827bd09bSSatish Balay 
852a501084fSBarry Smith   mvi              = (mv_info*)malloc(sizeof(mv_info));
853827bd09bSSatish Balay   mvi->n           = n;
854827bd09bSSatish Balay   mvi->m           = m;
855827bd09bSSatish Balay   mvi->n_global    = -1;
856827bd09bSSatish Balay   mvi->m_global    = -1;
85752f87cdaSBarry Smith   mvi->local2global= (PetscInt*)malloc((m+1)*sizeof(PetscInt));
858*2fa5cd67SKarl Rupp 
859ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global,local2global,m);
860827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
8615c8f6a95SKarl Rupp   mvi->matvec          = matvec;
862827bd09bSSatish Balay   mvi->grid_data       = grid_data;
863827bd09bSSatish Balay 
864827bd09bSSatish Balay   /* set xyt communication handle to perform restricted matvec */
865ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
866827bd09bSSatish Balay 
867827bd09bSSatish Balay   return(mvi);
868827bd09bSSatish Balay }
869827bd09bSSatish Balay 
8707b1ae94cSBarry Smith /**************************************xyt.c***********************************/
8713fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
872827bd09bSSatish Balay {
8730924e98cSBarry Smith   PetscFunctionBegin;
874827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
8750924e98cSBarry Smith   PetscFunctionReturn(0);
876827bd09bSSatish Balay }
877827bd09bSSatish Balay 
878827bd09bSSatish Balay 
879827bd09bSSatish Balay 
880