xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision 8b83055f287c8deb1a16325ae3faee24d0648b7a)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /*************************************xxt.c************************************
3827bd09bSSatish Balay Module Name: xxt
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 **************************************xxt.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 xxt_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    *col_sz, *col_indices;
30a501084fSBarry Smith   PetscScalar **col_vals, *x, *solve_uu, *solve_w;
3152f87cdaSBarry Smith   PetscInt    nsolves;
32a501084fSBarry Smith   PetscScalar tot_solve_time;
33827bd09bSSatish Balay } xxt_info;
34827bd09bSSatish Balay 
35827bd09bSSatish Balay typedef struct matvec_info {
3652f87cdaSBarry Smith   PetscInt     n, m, n_global, m_global;
3752f87cdaSBarry Smith   PetscInt     *local2global;
38ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle;
39a501084fSBarry Smith   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
40827bd09bSSatish Balay   void *grid_data;
41827bd09bSSatish Balay } mv_info;
42827bd09bSSatish Balay 
43827bd09bSSatish Balay struct xxt_CDT {
4452f87cdaSBarry Smith   PetscInt id;
4552f87cdaSBarry Smith   PetscInt ns;
4652f87cdaSBarry Smith   PetscInt level;
47827bd09bSSatish Balay   xxt_info *info;
48827bd09bSSatish Balay   mv_info  *mvi;
49827bd09bSSatish Balay };
50827bd09bSSatish Balay 
5152f87cdaSBarry Smith static PetscInt n_xxt        =0;
5252f87cdaSBarry Smith static PetscInt n_xxt_handles=0;
53827bd09bSSatish Balay 
54827bd09bSSatish Balay /* prototypes */
553fdc5746SBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *rhs);
563fdc5746SBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle);
573fdc5746SBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle);
583fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
5952f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle);
6052f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle);
615c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data);
62a501084fSBarry Smith 
637b1ae94cSBarry Smith /**************************************xxt.c***********************************/
647b1ae94cSBarry Smith xxt_ADT XXT_new(void)
65827bd09bSSatish Balay {
66827bd09bSSatish Balay   xxt_ADT xxt_handle;
67827bd09bSSatish Balay 
68827bd09bSSatish Balay   /* rolling count on n_xxt ... pot. problem here */
69827bd09bSSatish Balay   n_xxt_handles++;
70a501084fSBarry Smith   xxt_handle       = (xxt_ADT)malloc(sizeof(struct xxt_CDT));
71827bd09bSSatish Balay   xxt_handle->id   = ++n_xxt;
72827bd09bSSatish Balay   xxt_handle->info = NULL; xxt_handle->mvi  = NULL;
73827bd09bSSatish Balay 
74827bd09bSSatish Balay   return(xxt_handle);
75827bd09bSSatish Balay }
76827bd09bSSatish Balay 
777b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7852f87cdaSBarry Smith PetscInt XXT_factor(xxt_ADT xxt_handle,     /* prev. allocated xxt  handle */
7952f87cdaSBarry Smith                     PetscInt *local2global, /* global column mapping       */
8052f87cdaSBarry Smith                     PetscInt n,             /* local num rows              */
8152f87cdaSBarry Smith                     PetscInt m,             /* local num cols              */
825c8f6a95SKarl Rupp                     PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc         */
831147fc2aSKarl Rupp                     void *grid_data)        /* grid data for matvec        */
84827bd09bSSatish Balay {
85b1c944f5SJed Brown   PCTFS_comm_init();
86827bd09bSSatish Balay   check_handle(xxt_handle);
87827bd09bSSatish Balay 
88827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
89b1c944f5SJed Brown   if ((1<<(xxt_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);
90827bd09bSSatish Balay 
91827bd09bSSatish Balay   /* space for X info */
92a501084fSBarry Smith   xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info));
93827bd09bSSatish Balay 
94827bd09bSSatish Balay   /* set up matvec handles */
955c8f6a95SKarl Rupp   xxt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data);
96827bd09bSSatish Balay 
97827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
98827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
99827bd09bSSatish Balay   xxt_handle->ns=0;
100827bd09bSSatish Balay 
101827bd09bSSatish Balay   /* determine separators and generate firing order - NB xxt info set here */
102827bd09bSSatish Balay   det_separators(xxt_handle);
103827bd09bSSatish Balay 
104827bd09bSSatish Balay   return(do_xxt_factor(xxt_handle));
105827bd09bSSatish Balay }
106827bd09bSSatish Balay 
1077b1ae94cSBarry Smith /**************************************xxt.c***********************************/
1088cda6cd7SBarry Smith PetscInt XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b)
109827bd09bSSatish Balay {
110827bd09bSSatish Balay 
111b1c944f5SJed Brown   PCTFS_comm_init();
112827bd09bSSatish Balay   check_handle(xxt_handle);
113827bd09bSSatish Balay 
114827bd09bSSatish Balay   /* need to copy b into x? */
1152fa5cd67SKarl Rupp   if (b) PCTFS_rvec_copy(x,b,xxt_handle->mvi->n);
116827bd09bSSatish Balay   do_xxt_solve(xxt_handle,x);
117827bd09bSSatish Balay 
118827bd09bSSatish Balay   return(0);
119827bd09bSSatish Balay }
120827bd09bSSatish Balay 
1217b1ae94cSBarry Smith /**************************************xxt.c***********************************/
12252f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle)
123827bd09bSSatish Balay {
124827bd09bSSatish Balay 
125b1c944f5SJed Brown   PCTFS_comm_init();
126827bd09bSSatish Balay   check_handle(xxt_handle);
127827bd09bSSatish Balay   n_xxt_handles--;
128827bd09bSSatish Balay 
129a501084fSBarry Smith   free(xxt_handle->info->nsep);
130a501084fSBarry Smith   free(xxt_handle->info->lnsep);
131a501084fSBarry Smith   free(xxt_handle->info->fo);
132a501084fSBarry Smith   free(xxt_handle->info->stages);
133a501084fSBarry Smith   free(xxt_handle->info->solve_uu);
134a501084fSBarry Smith   free(xxt_handle->info->solve_w);
135a501084fSBarry Smith   free(xxt_handle->info->x);
136a501084fSBarry Smith   free(xxt_handle->info->col_vals);
137a501084fSBarry Smith   free(xxt_handle->info->col_sz);
138a501084fSBarry Smith   free(xxt_handle->info->col_indices);
139a501084fSBarry Smith   free(xxt_handle->info);
140a501084fSBarry Smith   free(xxt_handle->mvi->local2global);
141ca8e9878SJed Brown   PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle);
142a501084fSBarry Smith   free(xxt_handle->mvi);
143a501084fSBarry Smith   free(xxt_handle);
144827bd09bSSatish Balay 
145827bd09bSSatish Balay   /* if the check fails we nuke */
146a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
147827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
148827bd09bSSatish Balay   return(0);
149827bd09bSSatish Balay }
150827bd09bSSatish Balay 
1517b1ae94cSBarry Smith /**************************************xxt.c***********************************/
15252f87cdaSBarry Smith PetscInt XXT_stats(xxt_ADT xxt_handle)
153827bd09bSSatish Balay {
15452f87cdaSBarry Smith   PetscInt       op[]  = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
15552f87cdaSBarry Smith   PetscInt       fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
15652f87cdaSBarry Smith   PetscInt       vals[9],  work[9];
157a501084fSBarry Smith   PetscScalar    fvals[3], fwork[3];
15871a0148aSBarry Smith   PetscErrorCode ierr;
159827bd09bSSatish Balay 
160b1c944f5SJed Brown   PCTFS_comm_init();
161827bd09bSSatish Balay   check_handle(xxt_handle);
162827bd09bSSatish Balay 
163827bd09bSSatish Balay   /* if factorization not done there are no stats */
164db4deed7SKarl Rupp   if (!xxt_handle->info||!xxt_handle->mvi) {
165b1c944f5SJed Brown     if (!PCTFS_my_id) { ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");CHKERRQ(ierr); }
166827bd09bSSatish Balay     return 1;
167827bd09bSSatish Balay   }
168827bd09bSSatish Balay 
169827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
170827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
171827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
172b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
173827bd09bSSatish Balay 
1742fa5cd67SKarl Rupp   fvals[0]=fvals[1]=fvals[2] =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
175b1c944f5SJed Brown   PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
176827bd09bSSatish Balay 
177db4deed7SKarl Rupp   if (!PCTFS_my_id) {
178b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_nnz=%D\n",PCTFS_my_id,vals[0]);CHKERRQ(ierr);
179b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_nnz=%D\n",PCTFS_my_id,vals[1]);CHKERRQ(ierr);
180b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);CHKERRQ(ierr);
181b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_nnz=%D\n",PCTFS_my_id,vals[2]);CHKERRQ(ierr);
182b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(2d)  =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.5)));CHKERRQ(ierr);
183b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt   C(3d)  =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.6667)));CHKERRQ(ierr);
184b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_n  =%D\n",PCTFS_my_id,vals[3]);CHKERRQ(ierr);
185b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_n  =%D\n",PCTFS_my_id,vals[4]);CHKERRQ(ierr);
186b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_n  =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);CHKERRQ(ierr);
187b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_n  =%D\n",PCTFS_my_id,vals[5]);CHKERRQ(ierr);
188b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_buf=%D\n",PCTFS_my_id,vals[6]);CHKERRQ(ierr);
189b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_buf=%D\n",PCTFS_my_id,vals[7]);CHKERRQ(ierr);
190b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);CHKERRQ(ierr);
191b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_slv=%g\n",PCTFS_my_id,fvals[0]);CHKERRQ(ierr);
192b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_slv=%g\n",PCTFS_my_id,fvals[1]);CHKERRQ(ierr);
193b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);CHKERRQ(ierr);
194827bd09bSSatish Balay   }
195827bd09bSSatish Balay 
196827bd09bSSatish Balay   return(0);
197827bd09bSSatish Balay }
198827bd09bSSatish Balay 
199827bd09bSSatish Balay /*************************************xxt.c************************************
200827bd09bSSatish Balay 
201827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
202827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
203827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
204827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
205827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
206827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
207ca8e9878SJed Brown    PCTFS_gs_init/gop).
208827bd09bSSatish Balay 
209827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
210827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
211827bd09bSSatish Balay **************************************xxt.c***********************************/
21252f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle)
213827bd09bSSatish Balay {
2147b1ae94cSBarry Smith   return xxt_generate(xxt_handle);
215827bd09bSSatish Balay }
216827bd09bSSatish Balay 
2177b1ae94cSBarry Smith /**************************************xxt.c***********************************/
21852f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle)
219827bd09bSSatish Balay {
22052f87cdaSBarry Smith   PetscInt       i,j,k,idex;
22152f87cdaSBarry Smith   PetscInt       dim, col;
222a501084fSBarry Smith   PetscScalar    *u, *uu, *v, *z, *w, alpha, alpha_w;
22352f87cdaSBarry Smith   PetscInt       *segs;
22452f87cdaSBarry Smith   PetscInt       op[] = {GL_ADD,0};
22552f87cdaSBarry Smith   PetscInt       off, len;
226a501084fSBarry Smith   PetscScalar    *x_ptr;
22752f87cdaSBarry Smith   PetscInt       *iptr, flag;
22852f87cdaSBarry Smith   PetscInt       start =0, end, work;
22952f87cdaSBarry Smith   PetscInt       op2[] = {GL_MIN,0};
230ca8e9878SJed Brown   PCTFS_gs_ADT   PCTFS_gs_handle;
23152f87cdaSBarry Smith   PetscInt       *nsep, *lnsep, *fo;
23252f87cdaSBarry Smith   PetscInt       a_n            =xxt_handle->mvi->n;
23352f87cdaSBarry Smith   PetscInt       a_m            =xxt_handle->mvi->m;
23452f87cdaSBarry Smith   PetscInt       *a_local2global=xxt_handle->mvi->local2global;
23552f87cdaSBarry Smith   PetscInt       level;
23652f87cdaSBarry Smith   PetscInt       xxt_nnz=0, xxt_max_nnz=0;
23752f87cdaSBarry Smith   PetscInt       n, m;
23852f87cdaSBarry Smith   PetscInt       *col_sz, *col_indices, *stages;
239a501084fSBarry Smith   PetscScalar    **col_vals, *x;
24052f87cdaSBarry Smith   PetscInt       n_global;
24152f87cdaSBarry Smith   PetscInt       xxt_zero_nnz  =0;
24252f87cdaSBarry Smith   PetscInt       xxt_zero_nnz_0=0;
24371a0148aSBarry Smith   PetscBLASInt   i1            = 1,dlen;
244a501084fSBarry Smith   PetscScalar    dm1           = -1.0;
245d1528f56SBarry Smith   PetscErrorCode ierr;
246827bd09bSSatish Balay 
247827bd09bSSatish Balay   n               = xxt_handle->mvi->n;
248827bd09bSSatish Balay   nsep            = xxt_handle->info->nsep;
249827bd09bSSatish Balay   lnsep           = xxt_handle->info->lnsep;
250827bd09bSSatish Balay   fo              = xxt_handle->info->fo;
251827bd09bSSatish Balay   end             = lnsep[0];
252827bd09bSSatish Balay   level           = xxt_handle->level;
253ca8e9878SJed Brown   PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
254827bd09bSSatish Balay 
255827bd09bSSatish Balay   /* is there a null space? */
256827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
2572fa5cd67SKarl Rupp   for (i=0, j=0; i<=level; i++) j+=nsep[i];
258827bd09bSSatish Balay 
259827bd09bSSatish Balay   m = j-xxt_handle->ns;
26022d28d08SBarry Smith   if (m!=j) {
26122d28d08SBarry Smith     ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);CHKERRQ(ierr);
26222d28d08SBarry Smith   }
263827bd09bSSatish Balay 
264827bd09bSSatish Balay   /* get and initialize storage for x local         */
265827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
26652f87cdaSBarry Smith   col_sz      = (PetscInt*) malloc(m*sizeof(PetscInt));
26752f87cdaSBarry Smith   col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
268a501084fSBarry Smith   col_vals    = (PetscScalar**) malloc(m*sizeof(PetscScalar*));
269db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
270827bd09bSSatish Balay     col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
271827bd09bSSatish Balay     col_vals[i]   = NULL;
272827bd09bSSatish Balay   }
273827bd09bSSatish Balay   col_indices[j]=-1;
274827bd09bSSatish Balay 
275827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
276827bd09bSSatish Balay   /* this looks like nsep[]=segments */
27752f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
27852f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
279ca8e9878SJed Brown   PCTFS_ivec_zero(stages,level+1);
280ca8e9878SJed Brown   PCTFS_ivec_copy(segs,nsep,level+1);
2812fa5cd67SKarl Rupp   for (i=0; i<level; i++) segs[i+1] += segs[i];
282827bd09bSSatish Balay   stages[0] = segs[0];
283827bd09bSSatish Balay 
284827bd09bSSatish Balay   /* temporary vectors  */
285a501084fSBarry Smith   u  = (PetscScalar*) malloc(n*sizeof(PetscScalar));
286a501084fSBarry Smith   z  = (PetscScalar*) malloc(n*sizeof(PetscScalar));
287a501084fSBarry Smith   v  = (PetscScalar*) malloc(a_m*sizeof(PetscScalar));
288a501084fSBarry Smith   uu = (PetscScalar*) malloc(m*sizeof(PetscScalar));
289a501084fSBarry Smith   w  = (PetscScalar*) malloc(m*sizeof(PetscScalar));
290827bd09bSSatish Balay 
291827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
2922fa5cd67SKarl Rupp   for (i=1, j=0; i<=level; i++) j+=nsep[i];
293827bd09bSSatish Balay 
294827bd09bSSatish Balay   /* storage for sparse x values */
295827bd09bSSatish Balay   n_global    = xxt_handle->info->n_global;
29685ec1a3cSBarry Smith   xxt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes;
297a501084fSBarry Smith   x           = (PetscScalar*) malloc(xxt_max_nnz*sizeof(PetscScalar));
298827bd09bSSatish Balay   xxt_nnz     = 0;
299827bd09bSSatish Balay 
300827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
301827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
3022fa5cd67SKarl Rupp   for (dim=i=j=0; i<m; i++) {
303827bd09bSSatish Balay     /* time to move to the next level? */
304d4af0d30SBarry Smith     while (i==segs[dim]) {
305e32f2f54SBarry Smith       if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n");
306827bd09bSSatish Balay       stages[dim++]=i;
307827bd09bSSatish Balay       end         +=lnsep[dim];
308827bd09bSSatish Balay     }
309827bd09bSSatish Balay     stages[dim]=i;
310827bd09bSSatish Balay 
311827bd09bSSatish Balay     /* which column are we firing? */
312827bd09bSSatish Balay     /* i.e. set v_l */
313827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
314827bd09bSSatish Balay     (start<end) ? (col=fo[start]) : (col=INT_MAX);
315b1c944f5SJed Brown     PCTFS_giop_hc(&col,&work,1,op2,dim);
316827bd09bSSatish Balay 
317827bd09bSSatish Balay     /* shouldn't need this */
318db4deed7SKarl Rupp     if (col==INT_MAX) {
319f1ed62a8SBarry Smith       ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
320827bd09bSSatish Balay       continue;
321827bd09bSSatish Balay     }
322827bd09bSSatish Balay 
323827bd09bSSatish Balay     /* do I own it? I should */
324ca8e9878SJed Brown     PCTFS_rvec_zero(v,a_m);
325db4deed7SKarl Rupp     if (col==fo[start]) {
326827bd09bSSatish Balay       start++;
327ca8e9878SJed Brown       idex=PCTFS_ivec_linear_search(col, a_local2global, a_n);
328e7e72b3dSBarry Smith       if (idex!=-1) {
329e7e72b3dSBarry Smith         v[idex] = 1.0; j++;
330e7e72b3dSBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");
331db4deed7SKarl Rupp     } else {
332ca8e9878SJed Brown       idex=PCTFS_ivec_linear_search(col, a_local2global, a_m);
3332fa5cd67SKarl Rupp       if (idex!=-1) v[idex] = 1.0;
334827bd09bSSatish Balay     }
335827bd09bSSatish Balay 
336827bd09bSSatish Balay     /* perform u = A.v_l */
337ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
338827bd09bSSatish Balay     do_matvec(xxt_handle->mvi,v,u);
339827bd09bSSatish Balay 
340827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
341827bd09bSSatish Balay     /* technically only need to zero out first i entries */
342827bd09bSSatish Balay     /* later turn this into an XXT_solve call ? */
343ca8e9878SJed Brown     PCTFS_rvec_zero(uu,m);
344827bd09bSSatish Balay     x_ptr=x;
345827bd09bSSatish Balay     iptr = col_indices;
346db4deed7SKarl Rupp     for (k=0; k<i; k++) {
347827bd09bSSatish Balay       off   = *iptr++;
348827bd09bSSatish Balay       len   = *iptr++;
349c5df96a5SBarry Smith       ierr  = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
350*8b83055fSJed Brown       PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1));
351827bd09bSSatish Balay       x_ptr+=len;
352827bd09bSSatish Balay     }
353827bd09bSSatish Balay 
354827bd09bSSatish Balay 
355827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
356b1c944f5SJed Brown     PCTFS_ssgl_radd  (uu, w, dim, stages);
357827bd09bSSatish Balay 
358827bd09bSSatish Balay     /* z = X.uu */
359ca8e9878SJed Brown     PCTFS_rvec_zero(z,n);
360827bd09bSSatish Balay     x_ptr=x;
361827bd09bSSatish Balay     iptr = col_indices;
362db4deed7SKarl Rupp     for (k=0; k<i; k++) {
363827bd09bSSatish Balay       off  = *iptr++;
364827bd09bSSatish Balay       len  = *iptr++;
365c5df96a5SBarry Smith       ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
366*8b83055fSJed Brown       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1));
367827bd09bSSatish Balay       x_ptr+=len;
368827bd09bSSatish Balay     }
369827bd09bSSatish Balay 
370827bd09bSSatish Balay     /* compute v_l = v_l - z */
371ca8e9878SJed Brown     PCTFS_rvec_zero(v+a_n,a_m-a_n);
372c5df96a5SBarry Smith     ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr);
373*8b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1));
374827bd09bSSatish Balay 
375827bd09bSSatish Balay     /* compute u_l = A.v_l */
3762fa5cd67SKarl Rupp     if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim);
377ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
378827bd09bSSatish Balay     do_matvec(xxt_handle->mvi,v,u);
379827bd09bSSatish Balay 
380827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
381c5df96a5SBarry Smith     ierr  = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr);
382*8b83055fSJed Brown     PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,v,&i1));
383827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
384b1c944f5SJed Brown     PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
385827bd09bSSatish Balay 
3868f1a2a5eSBarry Smith     alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha);
387827bd09bSSatish Balay 
388827bd09bSSatish Balay     /* check for small alpha                             */
389827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
390c1235816SBarry Smith     if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);
391827bd09bSSatish Balay 
392827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
393ca8e9878SJed Brown     PCTFS_rvec_scale(v,1.0/alpha,n);
394827bd09bSSatish Balay 
395827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
396827bd09bSSatish Balay     flag = 1;
397827bd09bSSatish Balay     off=len=0;
398db4deed7SKarl Rupp     for (k=0; k<n; k++) {
399db4deed7SKarl Rupp       if (v[k]!=0.0) {
400827bd09bSSatish Balay         len=k;
401db4deed7SKarl Rupp         if (flag) { off=k; flag=0; }
402827bd09bSSatish Balay       }
403827bd09bSSatish Balay     }
404827bd09bSSatish Balay 
405827bd09bSSatish Balay     len -= (off-1);
406827bd09bSSatish Balay 
407db4deed7SKarl Rupp     if (len>0) {
408db4deed7SKarl Rupp       if ((xxt_nnz+len)>xxt_max_nnz) {
409f1ed62a8SBarry Smith         ierr         = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
410827bd09bSSatish Balay         xxt_max_nnz *= 2;
411a501084fSBarry Smith         x_ptr        = (PetscScalar*) malloc(xxt_max_nnz*sizeof(PetscScalar));
412ca8e9878SJed Brown         PCTFS_rvec_copy(x_ptr,x,xxt_nnz);
413a501084fSBarry Smith         free(x);
414827bd09bSSatish Balay         x     = x_ptr;
415827bd09bSSatish Balay         x_ptr+=xxt_nnz;
416827bd09bSSatish Balay       }
417827bd09bSSatish Balay       xxt_nnz += len;
418ca8e9878SJed Brown       PCTFS_rvec_copy(x_ptr,v+off,len);
419827bd09bSSatish Balay 
420827bd09bSSatish Balay       /* keep track of number of zeros */
421db4deed7SKarl Rupp       if (dim) {
422db4deed7SKarl Rupp         for (k=0; k<len; k++) {
4232fa5cd67SKarl Rupp           if (x_ptr[k]==0.0) xxt_zero_nnz++;
424827bd09bSSatish Balay         }
425db4deed7SKarl Rupp       } else {
426db4deed7SKarl Rupp         for (k=0; k<len; k++) {
4272fa5cd67SKarl Rupp           if (x_ptr[k]==0.0) xxt_zero_nnz_0++;
428827bd09bSSatish Balay         }
429827bd09bSSatish Balay       }
430827bd09bSSatish Balay       col_indices[2*i] = off;
431827bd09bSSatish Balay       col_sz[i] = col_indices[2*i+1] = len;
432827bd09bSSatish Balay       col_vals[i] = x_ptr;
433827bd09bSSatish Balay     }
434db4deed7SKarl Rupp     else {
435827bd09bSSatish Balay       col_indices[2*i] = 0;
436827bd09bSSatish Balay       col_sz[i]        = col_indices[2*i+1] = 0;
437827bd09bSSatish Balay       col_vals[i]      = x_ptr;
438827bd09bSSatish Balay     }
439827bd09bSSatish Balay   }
440827bd09bSSatish Balay 
441827bd09bSSatish Balay   /* close off stages for execution phase */
442db4deed7SKarl Rupp   while (dim!=level) {
443827bd09bSSatish Balay     stages[dim++] = i;
44471a0148aSBarry Smith     ierr          = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
445827bd09bSSatish Balay   }
446827bd09bSSatish Balay   stages[dim]=i;
447827bd09bSSatish Balay 
448827bd09bSSatish Balay   xxt_handle->info->n              = xxt_handle->mvi->n;
449827bd09bSSatish Balay   xxt_handle->info->m              = m;
450827bd09bSSatish Balay   xxt_handle->info->nnz            = xxt_nnz;
451827bd09bSSatish Balay   xxt_handle->info->max_nnz        = xxt_max_nnz;
452827bd09bSSatish Balay   xxt_handle->info->msg_buf_sz     = stages[level]-stages[0];
453a501084fSBarry Smith   xxt_handle->info->solve_uu       = (PetscScalar*) malloc(m*sizeof(PetscScalar));
454a501084fSBarry Smith   xxt_handle->info->solve_w        = (PetscScalar*) malloc(m*sizeof(PetscScalar));
455827bd09bSSatish Balay   xxt_handle->info->x              = x;
456827bd09bSSatish Balay   xxt_handle->info->col_vals       = col_vals;
457827bd09bSSatish Balay   xxt_handle->info->col_sz         = col_sz;
458827bd09bSSatish Balay   xxt_handle->info->col_indices    = col_indices;
459827bd09bSSatish Balay   xxt_handle->info->stages         = stages;
460827bd09bSSatish Balay   xxt_handle->info->nsolves        = 0;
461827bd09bSSatish Balay   xxt_handle->info->tot_solve_time = 0.0;
462827bd09bSSatish Balay 
463a501084fSBarry Smith   free(segs);
464a501084fSBarry Smith   free(u);
465a501084fSBarry Smith   free(v);
466a501084fSBarry Smith   free(uu);
467a501084fSBarry Smith   free(z);
468a501084fSBarry Smith   free(w);
469827bd09bSSatish Balay 
470827bd09bSSatish Balay   return(0);
471827bd09bSSatish Balay }
472827bd09bSSatish Balay 
4737b1ae94cSBarry Smith /**************************************xxt.c***********************************/
4740924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle,  PetscScalar *uc)
475827bd09bSSatish Balay {
47652f87cdaSBarry Smith   PetscInt       off, len, *iptr;
47752f87cdaSBarry Smith   PetscInt       level        = xxt_handle->level;
47852f87cdaSBarry Smith   PetscInt       n            = xxt_handle->info->n;
47952f87cdaSBarry Smith   PetscInt       m            = xxt_handle->info->m;
48052f87cdaSBarry Smith   PetscInt       *stages      = xxt_handle->info->stages;
48152f87cdaSBarry Smith   PetscInt       *col_indices = xxt_handle->info->col_indices;
482a501084fSBarry Smith   PetscScalar    *x_ptr, *uu_ptr;
483a501084fSBarry Smith   PetscScalar    *solve_uu = xxt_handle->info->solve_uu;
484a501084fSBarry Smith   PetscScalar    *solve_w  = xxt_handle->info->solve_w;
485a501084fSBarry Smith   PetscScalar    *x        = xxt_handle->info->x;
48671a0148aSBarry Smith   PetscBLASInt   i1        = 1,dlen;
487c5df96a5SBarry Smith   PetscErrorCode ierr;
488827bd09bSSatish Balay 
4890924e98cSBarry Smith   PetscFunctionBegin;
490827bd09bSSatish Balay   uu_ptr=solve_uu;
491ca8e9878SJed Brown   PCTFS_rvec_zero(uu_ptr,m);
492827bd09bSSatish Balay 
493827bd09bSSatish Balay   /* x  = X.Y^T.b */
494827bd09bSSatish Balay   /* uu = Y^T.b */
495db4deed7SKarl Rupp   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) {
496c5df96a5SBarry Smith     off       =*iptr++;
497c5df96a5SBarry Smith     len       =*iptr++;
498c5df96a5SBarry Smith     ierr      = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
499*8b83055fSJed Brown     PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1));
500827bd09bSSatish Balay   }
501827bd09bSSatish Balay 
502827bd09bSSatish Balay   /* comunication of beta */
503827bd09bSSatish Balay   uu_ptr=solve_uu;
5042fa5cd67SKarl Rupp   if (level) PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages);
505827bd09bSSatish Balay 
506ca8e9878SJed Brown   PCTFS_rvec_zero(uc,n);
507827bd09bSSatish Balay 
508827bd09bSSatish Balay   /* x = X.uu */
509db4deed7SKarl Rupp   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) {
510c5df96a5SBarry Smith     off  =*iptr++;
511c5df96a5SBarry Smith     len  =*iptr++;
512c5df96a5SBarry Smith     ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr);
513*8b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1));
514827bd09bSSatish Balay   }
5150924e98cSBarry Smith   PetscFunctionReturn(0);
516827bd09bSSatish Balay }
517827bd09bSSatish Balay 
5187b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5190924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle)
520827bd09bSSatish Balay {
52152f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
522827bd09bSSatish Balay 
5230924e98cSBarry Smith   PetscFunctionBegin;
524e32f2f54SBarry Smith   if (!xxt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xxt_handle);
525827bd09bSSatish Balay 
526827bd09bSSatish Balay   vals[0]=vals[1]=xxt_handle->id;
527b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
528e32f2f54SBarry Smith   if ((vals[0]!=vals[1])||(xxt_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], xxt_handle->id);
5290924e98cSBarry Smith   PetscFunctionReturn(0);
530827bd09bSSatish Balay }
531827bd09bSSatish Balay 
5327b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5330924e98cSBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle)
534827bd09bSSatish Balay {
53552f87cdaSBarry Smith   PetscInt     i, ct, id;
53652f87cdaSBarry Smith   PetscInt     mask, edge, *iptr;
53752f87cdaSBarry Smith   PetscInt     *dir, *used;
53852f87cdaSBarry Smith   PetscInt     sum[4], w[4];
539a501084fSBarry Smith   PetscScalar  rsum[4], rw[4];
54052f87cdaSBarry Smith   PetscInt     op[] = {GL_ADD,0};
541a501084fSBarry Smith   PetscScalar  *lhs, *rhs;
54252f87cdaSBarry Smith   PetscInt     *nsep, *lnsep, *fo, nfo=0;
543ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle;
54452f87cdaSBarry Smith   PetscInt     *local2global   = xxt_handle->mvi->local2global;
54552f87cdaSBarry Smith   PetscInt     n               = xxt_handle->mvi->n;
54652f87cdaSBarry Smith   PetscInt     m               = xxt_handle->mvi->m;
54752f87cdaSBarry Smith   PetscInt     level           = xxt_handle->level;
548ab824b78SBarry Smith   PetscInt     shared          = 0;
549827bd09bSSatish Balay 
5500924e98cSBarry Smith   PetscFunctionBegin;
55152f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
55252f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
55352f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
55452f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
55552f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
556827bd09bSSatish Balay 
557ca8e9878SJed Brown   PCTFS_ivec_zero(dir,level+1);
558ca8e9878SJed Brown   PCTFS_ivec_zero(nsep,level+1);
559ca8e9878SJed Brown   PCTFS_ivec_zero(lnsep,level+1);
560ca8e9878SJed Brown   PCTFS_ivec_set (fo,-1,n+1);
561ca8e9878SJed Brown   PCTFS_ivec_zero(used,n);
562827bd09bSSatish Balay 
5638cda6cd7SBarry Smith   lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
5648cda6cd7SBarry Smith   rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
565827bd09bSSatish Balay 
566827bd09bSSatish Balay   /* determine the # of unique dof */
567ca8e9878SJed Brown   PCTFS_rvec_zero(lhs,m);
568ca8e9878SJed Brown   PCTFS_rvec_set(lhs,1.0,n);
569ca8e9878SJed Brown   PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level);
570ca8e9878SJed Brown   PCTFS_rvec_zero(rsum,2);
571db4deed7SKarl Rupp   for (ct=i=0; i<n; i++) {
5722fa5cd67SKarl Rupp     if (lhs[i]!=0.0) {
5732fa5cd67SKarl Rupp       rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];
5742fa5cd67SKarl Rupp     }
575827bd09bSSatish Balay   }
576b1c944f5SJed Brown   PCTFS_grop_hc(rsum,rw,2,op,level);
577827bd09bSSatish Balay   rsum[0]+=0.1;
578827bd09bSSatish Balay   rsum[1]+=0.1;
579827bd09bSSatish Balay 
5802fa5cd67SKarl Rupp   if (fabs(rsum[0]-rsum[1])>EPS) shared=1;
581827bd09bSSatish Balay 
58252f87cdaSBarry Smith   xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0];
58352f87cdaSBarry Smith   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0];
584827bd09bSSatish Balay 
585827bd09bSSatish Balay   /* determine separator sets top down */
5862fa5cd67SKarl Rupp   if (shared) {
587db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
588db4deed7SKarl Rupp 
589827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
590ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
591ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
592827bd09bSSatish Balay 
593827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
594ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
595ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
596827bd09bSSatish Balay 
597db4deed7SKarl Rupp       for (i=0;i<n;i++) {
598db4deed7SKarl Rupp         if (id< mask) {
5992fa5cd67SKarl Rupp           if (lhs[i]!=0.0) lhs[i]=1.0;
600827bd09bSSatish Balay         }
601db4deed7SKarl Rupp 
602db4deed7SKarl Rupp         if (id>=mask) {
6032fa5cd67SKarl Rupp           if (rhs[i]!=0.0) rhs[i]=1.0;
604827bd09bSSatish Balay         }
605827bd09bSSatish Balay       }
606827bd09bSSatish Balay 
6072fa5cd67SKarl Rupp       if (id< mask) PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1);
6082fa5cd67SKarl Rupp       else          PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1);
609827bd09bSSatish Balay 
610827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
611ca8e9878SJed Brown       PCTFS_rvec_zero(rsum,4);
612db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
613db4deed7SKarl Rupp 
614db4deed7SKarl Rupp         if (!used[i]) {
615827bd09bSSatish Balay           /* number of unmarked dofs on node */
616827bd09bSSatish Balay           ct++;
617827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
618db4deed7SKarl Rupp           if (id< mask) {
619db4deed7SKarl Rupp             if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; }
620827bd09bSSatish Balay           }
621827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
622db4deed7SKarl Rupp           if (id>=mask) {
623db4deed7SKarl Rupp             if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; }
624827bd09bSSatish Balay           }
625827bd09bSSatish Balay         }
626827bd09bSSatish Balay       }
627827bd09bSSatish Balay 
628827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
629827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
630827bd09bSSatish Balay       (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
631b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
632b1c944f5SJed Brown       PCTFS_grop_hc(rsum,rw,4,op,edge);
633827bd09bSSatish Balay       rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
634827bd09bSSatish Balay 
635db4deed7SKarl Rupp       if (id<mask) {
636827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
637db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
638db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
639827bd09bSSatish Balay             ct++; nfo++;
640827bd09bSSatish Balay 
641c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
642827bd09bSSatish Balay 
643827bd09bSSatish Balay             *--iptr = local2global[i];
644827bd09bSSatish Balay             used[i] = edge;
645827bd09bSSatish Balay           }
646827bd09bSSatish Balay         }
6472fa5cd67SKarl Rupp         if (ct>1) PCTFS_ivec_sort(iptr,ct);
648827bd09bSSatish Balay 
649827bd09bSSatish Balay         lnsep[edge]=ct;
65052f87cdaSBarry Smith         nsep[edge]=(PetscInt) rsum[0];
651827bd09bSSatish Balay         dir [edge]=LEFT;
652827bd09bSSatish Balay       }
653827bd09bSSatish Balay 
654db4deed7SKarl Rupp       if (id>=mask) {
655827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
656db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
657db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
658827bd09bSSatish Balay             ct++; nfo++;
659827bd09bSSatish Balay 
660c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
661827bd09bSSatish Balay 
662827bd09bSSatish Balay             *--iptr = local2global[i];
663827bd09bSSatish Balay             used[i] = edge;
664827bd09bSSatish Balay           }
665827bd09bSSatish Balay         }
6662fa5cd67SKarl Rupp         if (ct>1) PCTFS_ivec_sort(iptr,ct);
667827bd09bSSatish Balay 
668827bd09bSSatish Balay         lnsep[edge] = ct;
66952f87cdaSBarry Smith         nsep[edge]  = (PetscInt) rsum[1];
670827bd09bSSatish Balay         dir [edge]  = RIGHT;
671827bd09bSSatish Balay       }
672827bd09bSSatish Balay 
673827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
674827bd09bSSatish Balay       /* do we need full set of separators for this?                */
675827bd09bSSatish Balay 
676827bd09bSSatish Balay       /* fold rhs hc into lower */
6772fa5cd67SKarl Rupp       if (id>=mask) id-=mask;
678827bd09bSSatish Balay     }
679db4deed7SKarl Rupp   } else {
680db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
681827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
682ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
683ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
684827bd09bSSatish Balay 
685827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
686ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
687ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
688827bd09bSSatish Balay 
689827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
690db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
691db4deed7SKarl Rupp         if (!used[i]) {
692827bd09bSSatish Balay           /* number of unmarked dofs on node */
693827bd09bSSatish Balay           ct++;
694827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
6952fa5cd67SKarl Rupp           if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++;
696827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
6972fa5cd67SKarl Rupp           if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++;
698827bd09bSSatish Balay         }
699827bd09bSSatish Balay       }
700827bd09bSSatish Balay 
701827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
702827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
703b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
704827bd09bSSatish Balay 
705827bd09bSSatish Balay       /* lhs hc wins */
706db4deed7SKarl Rupp       if (sum[2]>=sum[3]) {
707db4deed7SKarl Rupp         if (id<mask) {
708827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
709db4deed7SKarl Rupp           for (ct=i=0;i<n;i++) {
710db4deed7SKarl Rupp             if ((!used[i])&&(lhs[i]!=0.0)) {
711827bd09bSSatish Balay               ct++; nfo++;
712827bd09bSSatish Balay               *--iptr = local2global[i];
713827bd09bSSatish Balay               used[i]=edge;
714827bd09bSSatish Balay             }
715827bd09bSSatish Balay           }
7162fa5cd67SKarl Rupp           if (ct>1) PCTFS_ivec_sort(iptr,ct);
717827bd09bSSatish Balay           lnsep[edge]=ct;
718827bd09bSSatish Balay         }
719827bd09bSSatish Balay         nsep[edge]=sum[0];
720827bd09bSSatish Balay         dir [edge]=LEFT;
721db4deed7SKarl Rupp       } else { /* rhs hc wins */
722db4deed7SKarl Rupp         if (id>=mask) {
723827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
724db4deed7SKarl Rupp           for (ct=i=0;i<n;i++) {
725db4deed7SKarl Rupp             if ((!used[i])&&(rhs[i]!=0.0)) {
726827bd09bSSatish Balay               ct++; nfo++;
727827bd09bSSatish Balay               *--iptr = local2global[i];
728827bd09bSSatish Balay               used[i]=edge;
729827bd09bSSatish Balay             }
730827bd09bSSatish Balay           }
7312fa5cd67SKarl Rupp           if (ct>1) PCTFS_ivec_sort(iptr,ct);
732827bd09bSSatish Balay           lnsep[edge]=ct;
733827bd09bSSatish Balay         }
734827bd09bSSatish Balay         nsep[edge]=sum[1];
735827bd09bSSatish Balay         dir [edge]=RIGHT;
736827bd09bSSatish Balay       }
737827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
738827bd09bSSatish Balay       /* do we need full set of separators for this?                */
739827bd09bSSatish Balay 
740827bd09bSSatish Balay       /* fold rhs hc into lower */
7412fa5cd67SKarl Rupp       if (id>=mask) id-=mask;
742827bd09bSSatish Balay     }
743827bd09bSSatish Balay   }
744827bd09bSSatish Balay 
745827bd09bSSatish Balay 
746827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
747db4deed7SKarl Rupp   for (ct=i=0; i<n; i++) {
748db4deed7SKarl Rupp     if (!used[i]) {
749827bd09bSSatish Balay       ct++; nfo++;
750827bd09bSSatish Balay       *--iptr = local2global[i];
751827bd09bSSatish Balay       used[i] = edge;
752827bd09bSSatish Balay     }
753827bd09bSSatish Balay   }
7542fa5cd67SKarl Rupp   if (ct>1) PCTFS_ivec_sort(iptr,ct);
755827bd09bSSatish Balay   lnsep[edge]=ct;
756827bd09bSSatish Balay   nsep [edge]=ct;
757827bd09bSSatish Balay   dir  [edge]=LEFT;
758827bd09bSSatish Balay 
759827bd09bSSatish Balay   xxt_handle->info->nsep  = nsep;
760827bd09bSSatish Balay   xxt_handle->info->lnsep = lnsep;
761827bd09bSSatish Balay   xxt_handle->info->fo    = fo;
762827bd09bSSatish Balay   xxt_handle->info->nfo   = nfo;
763827bd09bSSatish Balay 
764a501084fSBarry Smith   free(dir);
765a501084fSBarry Smith   free(lhs);
766a501084fSBarry Smith   free(rhs);
767a501084fSBarry Smith   free(used);
7680924e98cSBarry Smith   PetscFunctionReturn(0);
769827bd09bSSatish Balay }
770827bd09bSSatish Balay 
7717b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7725c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global,PetscInt n,PetscInt m,PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*),void *grid_data)
773827bd09bSSatish Balay {
774827bd09bSSatish Balay   mv_info *mvi;
775827bd09bSSatish Balay 
776827bd09bSSatish Balay 
777a501084fSBarry Smith   mvi               = (mv_info*)malloc(sizeof(mv_info));
778827bd09bSSatish Balay   mvi->n            = n;
779827bd09bSSatish Balay   mvi->m            = m;
780827bd09bSSatish Balay   mvi->n_global     = -1;
781827bd09bSSatish Balay   mvi->m_global     = -1;
78252f87cdaSBarry Smith   mvi->local2global = (PetscInt*)malloc((m+1)*sizeof(PetscInt));
783ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global,local2global,m);
784827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
7855c8f6a95SKarl Rupp   mvi->matvec          = matvec;
786827bd09bSSatish Balay   mvi->grid_data       = grid_data;
787827bd09bSSatish Balay 
788827bd09bSSatish Balay   /* set xxt communication handle to perform restricted matvec */
789ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
790827bd09bSSatish Balay 
791827bd09bSSatish Balay   return(mvi);
792827bd09bSSatish Balay }
793827bd09bSSatish Balay 
7947b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7950924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
796827bd09bSSatish Balay {
7970924e98cSBarry Smith   PetscFunctionBegin;
798827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
7990924e98cSBarry Smith   PetscFunctionReturn(0);
800827bd09bSSatish Balay }
801827bd09bSSatish Balay 
802827bd09bSSatish Balay 
803827bd09bSSatish Balay 
804