xref: /petsc/src/ksp/pc/impls/tfs/xxt.c (revision db4deed73f06ae556a0f2e02a0fd6e016e156849)
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);
6152f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, 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              */
82827bd09bSSatish Balay            void *matvec,       /* b_loc=A_local.x_loc         */
83827bd09bSSatish Balay            void *grid_data     /* grid data for matvec        */
84827bd09bSSatish Balay            )
85827bd09bSSatish Balay {
86b1c944f5SJed Brown   PCTFS_comm_init();
87827bd09bSSatish Balay   check_handle(xxt_handle);
88827bd09bSSatish Balay 
89827bd09bSSatish Balay   /* only 2^k for now and all nodes participating */
90b1c944f5SJed 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);
91827bd09bSSatish Balay 
92827bd09bSSatish Balay   /* space for X info */
93a501084fSBarry Smith   xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info));
94827bd09bSSatish Balay 
95827bd09bSSatish Balay   /* set up matvec handles */
96827bd09bSSatish Balay   xxt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);
97827bd09bSSatish Balay 
98827bd09bSSatish Balay   /* matrix is assumed to be of full rank */
99827bd09bSSatish Balay   /* LATER we can reset to indicate rank def. */
100827bd09bSSatish Balay   xxt_handle->ns=0;
101827bd09bSSatish Balay 
102827bd09bSSatish Balay   /* determine separators and generate firing order - NB xxt info set here */
103827bd09bSSatish Balay   det_separators(xxt_handle);
104827bd09bSSatish Balay 
105827bd09bSSatish Balay   return(do_xxt_factor(xxt_handle));
106827bd09bSSatish Balay }
107827bd09bSSatish Balay 
1087b1ae94cSBarry Smith /**************************************xxt.c***********************************/
1098cda6cd7SBarry Smith PetscInt XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b)
110827bd09bSSatish Balay {
111827bd09bSSatish Balay 
112b1c944f5SJed Brown   PCTFS_comm_init();
113827bd09bSSatish Balay   check_handle(xxt_handle);
114827bd09bSSatish Balay 
115827bd09bSSatish Balay   /* need to copy b into x? */
116*db4deed7SKarl Rupp   if (b) {PCTFS_rvec_copy(x,b,xxt_handle->mvi->n);}
117827bd09bSSatish Balay   do_xxt_solve(xxt_handle,x);
118827bd09bSSatish Balay 
119827bd09bSSatish Balay   return(0);
120827bd09bSSatish Balay }
121827bd09bSSatish Balay 
1227b1ae94cSBarry Smith /**************************************xxt.c***********************************/
12352f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle)
124827bd09bSSatish Balay {
125827bd09bSSatish Balay 
126b1c944f5SJed Brown   PCTFS_comm_init();
127827bd09bSSatish Balay   check_handle(xxt_handle);
128827bd09bSSatish Balay   n_xxt_handles--;
129827bd09bSSatish Balay 
130a501084fSBarry Smith   free(xxt_handle->info->nsep);
131a501084fSBarry Smith   free(xxt_handle->info->lnsep);
132a501084fSBarry Smith   free(xxt_handle->info->fo);
133a501084fSBarry Smith   free(xxt_handle->info->stages);
134a501084fSBarry Smith   free(xxt_handle->info->solve_uu);
135a501084fSBarry Smith   free(xxt_handle->info->solve_w);
136a501084fSBarry Smith   free(xxt_handle->info->x);
137a501084fSBarry Smith   free(xxt_handle->info->col_vals);
138a501084fSBarry Smith   free(xxt_handle->info->col_sz);
139a501084fSBarry Smith   free(xxt_handle->info->col_indices);
140a501084fSBarry Smith   free(xxt_handle->info);
141a501084fSBarry Smith   free(xxt_handle->mvi->local2global);
142ca8e9878SJed Brown    PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle);
143a501084fSBarry Smith   free(xxt_handle->mvi);
144a501084fSBarry Smith   free(xxt_handle);
145827bd09bSSatish Balay 
146827bd09bSSatish Balay   /* if the check fails we nuke */
147a501084fSBarry Smith   /* if NULL pointer passed to free we nuke */
148827bd09bSSatish Balay   /* if the calls to free fail that's not my problem */
149827bd09bSSatish Balay   return(0);
150827bd09bSSatish Balay }
151827bd09bSSatish Balay 
1527b1ae94cSBarry Smith /**************************************xxt.c***********************************/
15352f87cdaSBarry Smith PetscInt XXT_stats(xxt_ADT xxt_handle)
154827bd09bSSatish Balay {
15552f87cdaSBarry Smith   PetscInt       op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
15652f87cdaSBarry Smith   PetscInt       fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
15752f87cdaSBarry Smith   PetscInt       vals[9],  work[9];
158a501084fSBarry Smith   PetscScalar    fvals[3], fwork[3];
15971a0148aSBarry Smith   PetscErrorCode ierr;
160827bd09bSSatish Balay 
161b1c944f5SJed Brown   PCTFS_comm_init();
162827bd09bSSatish Balay   check_handle(xxt_handle);
163827bd09bSSatish Balay 
164827bd09bSSatish Balay   /* if factorization not done there are no stats */
165*db4deed7SKarl Rupp   if (!xxt_handle->info||!xxt_handle->mvi) {
166b1c944f5SJed Brown     if (!PCTFS_my_id) { ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");CHKERRQ(ierr); }
167827bd09bSSatish Balay     return 1;
168827bd09bSSatish Balay   }
169827bd09bSSatish Balay 
170827bd09bSSatish Balay   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
171827bd09bSSatish Balay   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
172827bd09bSSatish Balay   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
173b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
174827bd09bSSatish Balay 
175827bd09bSSatish Balay   fvals[0]=fvals[1]=fvals[2]
176827bd09bSSatish Balay     =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
177b1c944f5SJed Brown   PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
178827bd09bSSatish Balay 
179*db4deed7SKarl Rupp   if (!PCTFS_my_id) {
180b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_nnz=%D\n",PCTFS_my_id,vals[0]);CHKERRQ(ierr);
181b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_nnz=%D\n",PCTFS_my_id,vals[1]);CHKERRQ(ierr);
182b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);CHKERRQ(ierr);
183b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_nnz=%D\n",PCTFS_my_id,vals[2]);CHKERRQ(ierr);
184b1c944f5SJed 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);
185b1c944f5SJed 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);
186b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_n  =%D\n",PCTFS_my_id,vals[3]);CHKERRQ(ierr);
187b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_n  =%D\n",PCTFS_my_id,vals[4]);CHKERRQ(ierr);
188b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_n  =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);CHKERRQ(ierr);
189b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xxt_n  =%D\n",PCTFS_my_id,vals[5]);CHKERRQ(ierr);
190b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_buf=%D\n",PCTFS_my_id,vals[6]);CHKERRQ(ierr);
191b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_buf=%D\n",PCTFS_my_id,vals[7]);CHKERRQ(ierr);
192b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);CHKERRQ(ierr);
193b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xxt_slv=%g\n",PCTFS_my_id,fvals[0]);CHKERRQ(ierr);
194b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xxt_slv=%g\n",PCTFS_my_id,fvals[1]);CHKERRQ(ierr);
195b1c944f5SJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xxt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);CHKERRQ(ierr);
196827bd09bSSatish Balay   }
197827bd09bSSatish Balay 
198827bd09bSSatish Balay   return(0);
199827bd09bSSatish Balay }
200827bd09bSSatish Balay 
201827bd09bSSatish Balay /*************************************xxt.c************************************
202827bd09bSSatish Balay 
203827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which
204827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m.
205827bd09bSSatish Balay    o my_ml holds address of ML struct associated w/A_local and coarse grid
206827bd09bSSatish Balay    o local2global holds global number of column i (i=0,...,m-1)
207827bd09bSSatish Balay    o local2global holds global number of row    i (i=0,...,n-1)
208827bd09bSSatish Balay    o mylocmatvec performs A_local . vec_local (note that gs is performed using
209ca8e9878SJed Brown    PCTFS_gs_init/gop).
210827bd09bSSatish Balay 
211827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
212827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out)
213827bd09bSSatish Balay **************************************xxt.c***********************************/
21452f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle)
215827bd09bSSatish Balay {
2167b1ae94cSBarry Smith   return xxt_generate(xxt_handle);
217827bd09bSSatish Balay }
218827bd09bSSatish Balay 
2197b1ae94cSBarry Smith /**************************************xxt.c***********************************/
22052f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle)
221827bd09bSSatish Balay {
22252f87cdaSBarry Smith   PetscInt i,j,k,idex;
22352f87cdaSBarry Smith   PetscInt dim, col;
224a501084fSBarry Smith   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
22552f87cdaSBarry Smith   PetscInt *segs;
22652f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
22752f87cdaSBarry Smith   PetscInt off, len;
228a501084fSBarry Smith   PetscScalar *x_ptr;
22952f87cdaSBarry Smith   PetscInt *iptr, flag;
23052f87cdaSBarry Smith   PetscInt start=0, end, work;
23152f87cdaSBarry Smith   PetscInt op2[] = {GL_MIN,0};
232ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle;
23352f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo;
23452f87cdaSBarry Smith   PetscInt a_n=xxt_handle->mvi->n;
23552f87cdaSBarry Smith   PetscInt a_m=xxt_handle->mvi->m;
23652f87cdaSBarry Smith   PetscInt *a_local2global=xxt_handle->mvi->local2global;
23752f87cdaSBarry Smith   PetscInt level;
23852f87cdaSBarry Smith   PetscInt xxt_nnz=0, xxt_max_nnz=0;
23952f87cdaSBarry Smith   PetscInt n, m;
24052f87cdaSBarry Smith   PetscInt *col_sz, *col_indices, *stages;
241a501084fSBarry Smith   PetscScalar **col_vals, *x;
24252f87cdaSBarry Smith   PetscInt n_global;
24352f87cdaSBarry Smith   PetscInt xxt_zero_nnz=0;
24452f87cdaSBarry Smith   PetscInt xxt_zero_nnz_0=0;
24571a0148aSBarry Smith   PetscBLASInt i1 = 1,dlen;
246a501084fSBarry Smith   PetscScalar dm1 = -1.0;
247d1528f56SBarry Smith   PetscErrorCode ierr;
248827bd09bSSatish Balay 
249827bd09bSSatish Balay   n=xxt_handle->mvi->n;
250827bd09bSSatish Balay   nsep=xxt_handle->info->nsep;
251827bd09bSSatish Balay   lnsep=xxt_handle->info->lnsep;
252827bd09bSSatish Balay   fo=xxt_handle->info->fo;
253827bd09bSSatish Balay   end=lnsep[0];
254827bd09bSSatish Balay   level=xxt_handle->level;
255ca8e9878SJed Brown   PCTFS_gs_handle=xxt_handle->mvi->PCTFS_gs_handle;
256827bd09bSSatish Balay 
257827bd09bSSatish Balay   /* is there a null space? */
258827bd09bSSatish Balay   /* LATER add in ability to detect null space by checking alpha */
259827bd09bSSatish Balay   for (i=0, j=0; i<=level; i++)
260827bd09bSSatish Balay     {j+=nsep[i];}
261827bd09bSSatish Balay 
262827bd09bSSatish Balay   m = j-xxt_handle->ns;
26322d28d08SBarry Smith   if (m!=j) {
26422d28d08SBarry Smith     ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);CHKERRQ(ierr);
26522d28d08SBarry Smith   }
266827bd09bSSatish Balay 
267827bd09bSSatish Balay   /* get and initialize storage for x local         */
268827bd09bSSatish Balay   /* note that x local is nxm and stored by columns */
26952f87cdaSBarry Smith   col_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
27052f87cdaSBarry Smith   col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
271a501084fSBarry Smith   col_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
272*db4deed7SKarl Rupp   for (i=j=0; i<m; i++, j+=2) {
273827bd09bSSatish Balay     col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
274827bd09bSSatish Balay     col_vals[i] = NULL;
275827bd09bSSatish Balay   }
276827bd09bSSatish Balay   col_indices[j]=-1;
277827bd09bSSatish Balay 
278827bd09bSSatish Balay   /* size of separators for each sub-hc working from bottom of tree to top */
279827bd09bSSatish Balay   /* this looks like nsep[]=segments */
28052f87cdaSBarry Smith   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
28152f87cdaSBarry Smith   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
282ca8e9878SJed Brown   PCTFS_ivec_zero(stages,level+1);
283ca8e9878SJed Brown   PCTFS_ivec_copy(segs,nsep,level+1);
284*db4deed7SKarl Rupp   for (i=0; i<level; i++) { segs[i+1] += segs[i]; }
285827bd09bSSatish Balay   stages[0] = segs[0];
286827bd09bSSatish Balay 
287827bd09bSSatish Balay   /* temporary vectors  */
288a501084fSBarry Smith   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
289a501084fSBarry Smith   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
290a501084fSBarry Smith   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
291a501084fSBarry Smith   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
292a501084fSBarry Smith   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
293827bd09bSSatish Balay 
294827bd09bSSatish Balay   /* extra nnz due to replication of vertices across separators */
295*db4deed7SKarl Rupp   for (i=1, j=0; i<=level; i++) {j+=nsep[i];}
296827bd09bSSatish Balay 
297827bd09bSSatish Balay   /* storage for sparse x values */
298827bd09bSSatish Balay   n_global = xxt_handle->info->n_global;
299b1c944f5SJed Brown   xxt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes;
300a501084fSBarry Smith   x = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
301827bd09bSSatish Balay   xxt_nnz = 0;
302827bd09bSSatish Balay 
303827bd09bSSatish Balay   /* LATER - can embed next sep to fire in gs */
304827bd09bSSatish Balay   /* time to make the donuts - generate X factor */
305827bd09bSSatish Balay   for (dim=i=j=0;i<m;i++)
306827bd09bSSatish Balay   {
307827bd09bSSatish Balay     /* time to move to the next level? */
308d4af0d30SBarry Smith     while (i==segs[dim]) {
309e32f2f54SBarry Smith       if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n");
310827bd09bSSatish Balay       stages[dim++]=i;
311827bd09bSSatish Balay       end+=lnsep[dim];
312827bd09bSSatish Balay     }
313827bd09bSSatish Balay     stages[dim]=i;
314827bd09bSSatish Balay 
315827bd09bSSatish Balay     /* which column are we firing? */
316827bd09bSSatish Balay     /* i.e. set v_l */
317827bd09bSSatish Balay     /* use new seps and do global min across hc to determine which one to fire */
318827bd09bSSatish Balay     (start<end) ? (col=fo[start]) : (col=INT_MAX);
319b1c944f5SJed Brown     PCTFS_giop_hc(&col,&work,1,op2,dim);
320827bd09bSSatish Balay 
321827bd09bSSatish Balay     /* shouldn't need this */
322*db4deed7SKarl Rupp     if (col==INT_MAX) {
323f1ed62a8SBarry Smith       ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr);
324827bd09bSSatish Balay       continue;
325827bd09bSSatish Balay     }
326827bd09bSSatish Balay 
327827bd09bSSatish Balay     /* do I own it? I should */
328ca8e9878SJed Brown     PCTFS_rvec_zero(v ,a_m);
329*db4deed7SKarl Rupp     if (col==fo[start]) {
330827bd09bSSatish Balay       start++;
331ca8e9878SJed Brown       idex=PCTFS_ivec_linear_search(col, a_local2global, a_n);
332e7e72b3dSBarry Smith       if (idex!=-1) {
333e7e72b3dSBarry Smith         v[idex] = 1.0; j++;
334e7e72b3dSBarry Smith       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n");
335*db4deed7SKarl Rupp     } else {
336ca8e9878SJed Brown       idex=PCTFS_ivec_linear_search(col, a_local2global, a_m);
337827bd09bSSatish Balay       if (idex!=-1)
338827bd09bSSatish Balay         {v[idex] = 1.0;}
339827bd09bSSatish Balay     }
340827bd09bSSatish Balay 
341827bd09bSSatish Balay     /* perform u = A.v_l */
342ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
343827bd09bSSatish Balay     do_matvec(xxt_handle->mvi,v,u);
344827bd09bSSatish Balay 
345827bd09bSSatish Balay     /* uu =  X^T.u_l (local portion) */
346827bd09bSSatish Balay     /* technically only need to zero out first i entries */
347827bd09bSSatish Balay     /* later turn this into an XXT_solve call ? */
348ca8e9878SJed Brown     PCTFS_rvec_zero(uu,m);
349827bd09bSSatish Balay     x_ptr=x;
350827bd09bSSatish Balay     iptr = col_indices;
351*db4deed7SKarl Rupp     for (k=0; k<i; k++) {
352827bd09bSSatish Balay       off = *iptr++;
353827bd09bSSatish Balay       len = *iptr++;
3540805154bSBarry Smith       dlen = PetscBLASIntCast(len);
35571a0148aSBarry Smith       uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1);
356827bd09bSSatish Balay       x_ptr+=len;
357827bd09bSSatish Balay     }
358827bd09bSSatish Balay 
359827bd09bSSatish Balay 
360827bd09bSSatish Balay     /* uu = X^T.u_l (comm portion) */
361b1c944f5SJed Brown     PCTFS_ssgl_radd  (uu, w, dim, stages);
362827bd09bSSatish Balay 
363827bd09bSSatish Balay     /* z = X.uu */
364ca8e9878SJed Brown     PCTFS_rvec_zero(z,n);
365827bd09bSSatish Balay     x_ptr=x;
366827bd09bSSatish Balay     iptr = col_indices;
367*db4deed7SKarl Rupp     for (k=0; k<i; k++) {
368827bd09bSSatish Balay       off = *iptr++;
369827bd09bSSatish Balay       len = *iptr++;
3700805154bSBarry Smith       dlen = PetscBLASIntCast(len);
37171a0148aSBarry Smith       BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
372827bd09bSSatish Balay       x_ptr+=len;
373827bd09bSSatish Balay     }
374827bd09bSSatish Balay 
375827bd09bSSatish Balay     /* compute v_l = v_l - z */
376ca8e9878SJed Brown     PCTFS_rvec_zero(v+a_n,a_m-a_n);
3770805154bSBarry Smith     dlen = PetscBLASIntCast(n);
37871a0148aSBarry Smith     BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);
379827bd09bSSatish Balay 
380827bd09bSSatish Balay     /* compute u_l = A.v_l */
381*db4deed7SKarl Rupp     if (a_n!=a_m) { PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); }
382ca8e9878SJed Brown     PCTFS_rvec_zero(u,n);
383827bd09bSSatish Balay     do_matvec(xxt_handle->mvi,v,u);
384827bd09bSSatish Balay 
385827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
3860805154bSBarry Smith     dlen = PetscBLASIntCast(n);
38771a0148aSBarry Smith     alpha = BLASdot_(&dlen,u,&i1,v,&i1);
388827bd09bSSatish Balay     /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
389b1c944f5SJed Brown     PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
390827bd09bSSatish Balay 
3918f1a2a5eSBarry Smith     alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha);
392827bd09bSSatish Balay 
393827bd09bSSatish Balay     /* check for small alpha                             */
394827bd09bSSatish Balay     /* LATER use this to detect and determine null space */
395c1235816SBarry Smith     if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);
396827bd09bSSatish Balay 
397827bd09bSSatish Balay     /* compute v_l = v_l/sqrt(alpha) */
398ca8e9878SJed Brown     PCTFS_rvec_scale(v,1.0/alpha,n);
399827bd09bSSatish Balay 
400827bd09bSSatish Balay     /* add newly generated column, v_l, to X */
401827bd09bSSatish Balay     flag = 1;
402827bd09bSSatish Balay     off=len=0;
403*db4deed7SKarl Rupp     for (k=0; k<n; k++) {
404*db4deed7SKarl Rupp       if (v[k]!=0.0) {
405827bd09bSSatish Balay         len=k;
406*db4deed7SKarl Rupp         if (flag) { off=k; flag=0; }
407827bd09bSSatish Balay       }
408827bd09bSSatish Balay     }
409827bd09bSSatish Balay 
410827bd09bSSatish Balay     len -= (off-1);
411827bd09bSSatish Balay 
412*db4deed7SKarl Rupp     if (len>0) {
413*db4deed7SKarl Rupp       if ((xxt_nnz+len)>xxt_max_nnz) {
414f1ed62a8SBarry Smith         ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr);
415827bd09bSSatish Balay         xxt_max_nnz *= 2;
416a501084fSBarry Smith         x_ptr = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar));
417ca8e9878SJed Brown         PCTFS_rvec_copy(x_ptr,x,xxt_nnz);
418a501084fSBarry Smith         free(x);
419827bd09bSSatish Balay         x = x_ptr;
420827bd09bSSatish Balay         x_ptr+=xxt_nnz;
421827bd09bSSatish Balay       }
422827bd09bSSatish Balay       xxt_nnz += len;
423ca8e9878SJed Brown       PCTFS_rvec_copy(x_ptr,v+off,len);
424827bd09bSSatish Balay 
425827bd09bSSatish Balay       /* keep track of number of zeros */
426*db4deed7SKarl Rupp       if (dim) {
427*db4deed7SKarl Rupp         for (k=0; k<len; k++) {
428*db4deed7SKarl Rupp           if (x_ptr[k]==0.0) { xxt_zero_nnz++; }
429827bd09bSSatish Balay         }
430*db4deed7SKarl Rupp       } else {
431*db4deed7SKarl Rupp         for (k=0; k<len; k++) {
432*db4deed7SKarl Rupp           if (x_ptr[k]==0.0) {xxt_zero_nnz_0++;}
433827bd09bSSatish Balay         }
434827bd09bSSatish Balay       }
435827bd09bSSatish Balay       col_indices[2*i] = off;
436827bd09bSSatish Balay       col_sz[i] = col_indices[2*i+1] = len;
437827bd09bSSatish Balay       col_vals[i] = x_ptr;
438827bd09bSSatish Balay     }
439*db4deed7SKarl Rupp     else {
440827bd09bSSatish Balay       col_indices[2*i] = 0;
441827bd09bSSatish Balay       col_sz[i] = col_indices[2*i+1] = 0;
442827bd09bSSatish Balay       col_vals[i] = x_ptr;
443827bd09bSSatish Balay     }
444827bd09bSSatish Balay   }
445827bd09bSSatish Balay 
446827bd09bSSatish Balay   /* close off stages for execution phase */
447*db4deed7SKarl Rupp   while (dim!=level) {
448827bd09bSSatish Balay     stages[dim++]=i;
44971a0148aSBarry Smith     ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr);
450827bd09bSSatish Balay   }
451827bd09bSSatish Balay   stages[dim]=i;
452827bd09bSSatish Balay 
453827bd09bSSatish Balay   xxt_handle->info->n=xxt_handle->mvi->n;
454827bd09bSSatish Balay   xxt_handle->info->m=m;
455827bd09bSSatish Balay   xxt_handle->info->nnz=xxt_nnz;
456827bd09bSSatish Balay   xxt_handle->info->max_nnz=xxt_max_nnz;
457827bd09bSSatish Balay   xxt_handle->info->msg_buf_sz=stages[level]-stages[0];
458a501084fSBarry Smith   xxt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
459a501084fSBarry Smith   xxt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
460827bd09bSSatish Balay   xxt_handle->info->x=x;
461827bd09bSSatish Balay   xxt_handle->info->col_vals=col_vals;
462827bd09bSSatish Balay   xxt_handle->info->col_sz=col_sz;
463827bd09bSSatish Balay   xxt_handle->info->col_indices=col_indices;
464827bd09bSSatish Balay   xxt_handle->info->stages=stages;
465827bd09bSSatish Balay   xxt_handle->info->nsolves=0;
466827bd09bSSatish Balay   xxt_handle->info->tot_solve_time=0.0;
467827bd09bSSatish Balay 
468a501084fSBarry Smith   free(segs);
469a501084fSBarry Smith   free(u);
470a501084fSBarry Smith   free(v);
471a501084fSBarry Smith   free(uu);
472a501084fSBarry Smith   free(z);
473a501084fSBarry Smith   free(w);
474827bd09bSSatish Balay 
475827bd09bSSatish Balay   return(0);
476827bd09bSSatish Balay }
477827bd09bSSatish Balay 
4787b1ae94cSBarry Smith /**************************************xxt.c***********************************/
4790924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle,  PetscScalar *uc)
480827bd09bSSatish Balay {
48152f87cdaSBarry Smith   PetscInt off, len, *iptr;
48252f87cdaSBarry Smith   PetscInt level       =xxt_handle->level;
48352f87cdaSBarry Smith   PetscInt n           =xxt_handle->info->n;
48452f87cdaSBarry Smith   PetscInt m           =xxt_handle->info->m;
48552f87cdaSBarry Smith   PetscInt *stages     =xxt_handle->info->stages;
48652f87cdaSBarry Smith   PetscInt *col_indices=xxt_handle->info->col_indices;
487a501084fSBarry Smith   PetscScalar *x_ptr, *uu_ptr;
488a501084fSBarry Smith   PetscScalar *solve_uu=xxt_handle->info->solve_uu;
489a501084fSBarry Smith   PetscScalar *solve_w =xxt_handle->info->solve_w;
490a501084fSBarry Smith   PetscScalar *x       =xxt_handle->info->x;
49171a0148aSBarry Smith   PetscBLASInt i1 = 1,dlen;
492827bd09bSSatish Balay 
4930924e98cSBarry Smith   PetscFunctionBegin;
494827bd09bSSatish Balay   uu_ptr=solve_uu;
495ca8e9878SJed Brown   PCTFS_rvec_zero(uu_ptr,m);
496827bd09bSSatish Balay 
497827bd09bSSatish Balay   /* x  = X.Y^T.b */
498827bd09bSSatish Balay   /* uu = Y^T.b */
499*db4deed7SKarl Rupp   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) {
500827bd09bSSatish Balay     off=*iptr++; len=*iptr++;
5010805154bSBarry Smith     dlen = PetscBLASIntCast(len);
50271a0148aSBarry Smith     *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1);
503827bd09bSSatish Balay   }
504827bd09bSSatish Balay 
505827bd09bSSatish Balay   /* comunication of beta */
506827bd09bSSatish Balay   uu_ptr=solve_uu;
507b1c944f5SJed Brown   if (level) { PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages); }
508827bd09bSSatish Balay 
509ca8e9878SJed Brown   PCTFS_rvec_zero(uc,n);
510827bd09bSSatish Balay 
511827bd09bSSatish Balay   /* x = X.uu */
512*db4deed7SKarl Rupp   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) {
513827bd09bSSatish Balay     off=*iptr++; len=*iptr++;
5140805154bSBarry Smith     dlen = PetscBLASIntCast(len);
51571a0148aSBarry Smith     BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
516827bd09bSSatish Balay   }
5170924e98cSBarry Smith   PetscFunctionReturn(0);
518827bd09bSSatish Balay }
519827bd09bSSatish Balay 
5207b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5210924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle)
522827bd09bSSatish Balay {
52352f87cdaSBarry Smith   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
524827bd09bSSatish Balay 
5250924e98cSBarry Smith   PetscFunctionBegin;
526e32f2f54SBarry Smith   if (!xxt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xxt_handle);
527827bd09bSSatish Balay 
528827bd09bSSatish Balay   vals[0]=vals[1]=xxt_handle->id;
529b1c944f5SJed Brown   PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
530e32f2f54SBarry 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);
5310924e98cSBarry Smith   PetscFunctionReturn(0);
532827bd09bSSatish Balay }
533827bd09bSSatish Balay 
5347b1ae94cSBarry Smith /**************************************xxt.c***********************************/
5350924e98cSBarry Smith static  PetscErrorCode det_separators(xxt_ADT xxt_handle)
536827bd09bSSatish Balay {
53752f87cdaSBarry Smith   PetscInt i, ct, id;
53852f87cdaSBarry Smith   PetscInt mask, edge, *iptr;
53952f87cdaSBarry Smith   PetscInt *dir, *used;
54052f87cdaSBarry Smith   PetscInt sum[4], w[4];
541a501084fSBarry Smith   PetscScalar rsum[4], rw[4];
54252f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
543a501084fSBarry Smith   PetscScalar *lhs, *rhs;
54452f87cdaSBarry Smith   PetscInt *nsep, *lnsep, *fo, nfo=0;
545ca8e9878SJed Brown   PCTFS_gs_ADT PCTFS_gs_handle=xxt_handle->mvi->PCTFS_gs_handle;
54652f87cdaSBarry Smith   PetscInt *local2global=xxt_handle->mvi->local2global;
54752f87cdaSBarry Smith   PetscInt  n=xxt_handle->mvi->n;
54852f87cdaSBarry Smith   PetscInt  m=xxt_handle->mvi->m;
54952f87cdaSBarry Smith   PetscInt level=xxt_handle->level;
550ab824b78SBarry Smith   PetscInt shared=0;
551827bd09bSSatish Balay 
5520924e98cSBarry Smith   PetscFunctionBegin;
55352f87cdaSBarry Smith   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
55452f87cdaSBarry Smith   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
55552f87cdaSBarry Smith   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
55652f87cdaSBarry Smith   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
55752f87cdaSBarry Smith   used = (PetscInt*)malloc(sizeof(PetscInt)*n);
558827bd09bSSatish Balay 
559ca8e9878SJed Brown   PCTFS_ivec_zero(dir  ,level+1);
560ca8e9878SJed Brown   PCTFS_ivec_zero(nsep ,level+1);
561ca8e9878SJed Brown   PCTFS_ivec_zero(lnsep,level+1);
562ca8e9878SJed Brown   PCTFS_ivec_set (fo   ,-1,n+1);
563ca8e9878SJed Brown   PCTFS_ivec_zero(used,n);
564827bd09bSSatish Balay 
5658cda6cd7SBarry Smith   lhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
5668cda6cd7SBarry Smith   rhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
567827bd09bSSatish Balay 
568827bd09bSSatish Balay   /* determine the # of unique dof */
569ca8e9878SJed Brown   PCTFS_rvec_zero(lhs,m);
570ca8e9878SJed Brown   PCTFS_rvec_set(lhs,1.0,n);
571ca8e9878SJed Brown   PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level);
572ca8e9878SJed Brown   PCTFS_rvec_zero(rsum,2);
573*db4deed7SKarl Rupp   for (ct=i=0;i<n;i++) {
574827bd09bSSatish Balay     if (lhs[i]!=0.0)
575827bd09bSSatish Balay       {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
576827bd09bSSatish Balay   }
577b1c944f5SJed Brown   PCTFS_grop_hc(rsum,rw,2,op,level);
578827bd09bSSatish Balay   rsum[0]+=0.1;
579827bd09bSSatish Balay   rsum[1]+=0.1;
580827bd09bSSatish Balay 
581*db4deed7SKarl Rupp   if (fabs(rsum[0]-rsum[1])>EPS) { shared=1; }
582827bd09bSSatish Balay 
58352f87cdaSBarry Smith   xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0];
58452f87cdaSBarry Smith   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0];
585827bd09bSSatish Balay 
586827bd09bSSatish Balay   /* determine separator sets top down */
587827bd09bSSatish Balay   if (shared)
588827bd09bSSatish Balay   {
589*db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
590*db4deed7SKarl Rupp 
591827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
592ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
593ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
594827bd09bSSatish Balay 
595827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
596ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
597ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
598827bd09bSSatish Balay 
599*db4deed7SKarl Rupp       for (i=0;i<n;i++) {
600*db4deed7SKarl Rupp         if (id< mask) {
601*db4deed7SKarl Rupp           if (lhs[i]!=0.0) { lhs[i]=1.0; }
602827bd09bSSatish Balay         }
603*db4deed7SKarl Rupp 
604*db4deed7SKarl Rupp         if (id>=mask) {
605*db4deed7SKarl Rupp           if (rhs[i]!=0.0) { rhs[i]=1.0; }
606827bd09bSSatish Balay         }
607827bd09bSSatish Balay       }
608827bd09bSSatish Balay 
609*db4deed7SKarl Rupp       if (id< mask) { PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1); }
610*db4deed7SKarl Rupp       else          { PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1); }
611827bd09bSSatish Balay 
612827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
613ca8e9878SJed Brown       PCTFS_rvec_zero(rsum,4);
614*db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
615*db4deed7SKarl Rupp 
616*db4deed7SKarl Rupp         if (!used[i]) {
617827bd09bSSatish Balay           /* number of unmarked dofs on node */
618827bd09bSSatish Balay           ct++;
619827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
620*db4deed7SKarl Rupp           if (id< mask) {
621*db4deed7SKarl Rupp             if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; }
622827bd09bSSatish Balay           }
623827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
624*db4deed7SKarl Rupp           if (id>=mask) {
625*db4deed7SKarl Rupp             if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; }
626827bd09bSSatish Balay           }
627827bd09bSSatish Balay         }
628*db4deed7SKarl Rupp 
629827bd09bSSatish Balay       }
630827bd09bSSatish Balay 
631827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
632827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
633827bd09bSSatish Balay       (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
634b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
635b1c944f5SJed Brown       PCTFS_grop_hc(rsum,rw,4,op,edge);
636827bd09bSSatish Balay       rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;
637827bd09bSSatish Balay 
638*db4deed7SKarl Rupp       if (id<mask) {
639827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
640*db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
641*db4deed7SKarl Rupp           if ((!used[i])&&(lhs[i]!=0.0)) {
642827bd09bSSatish Balay             ct++; nfo++;
643827bd09bSSatish Balay 
644c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
645827bd09bSSatish Balay 
646827bd09bSSatish Balay             *--iptr = local2global[i];
647827bd09bSSatish Balay             used[i]=edge;
648827bd09bSSatish Balay           }
649827bd09bSSatish Balay         }
650ca8e9878SJed Brown         if (ct>1) {PCTFS_ivec_sort(iptr,ct);}
651827bd09bSSatish Balay 
652827bd09bSSatish Balay         lnsep[edge]=ct;
65352f87cdaSBarry Smith         nsep[edge]=(PetscInt) rsum[0];
654827bd09bSSatish Balay         dir [edge]=LEFT;
655827bd09bSSatish Balay       }
656827bd09bSSatish Balay 
657*db4deed7SKarl Rupp       if (id>=mask) {
658827bd09bSSatish Balay         /* mark dofs I own that have signal and not in sep set */
659*db4deed7SKarl Rupp         for (ct=i=0;i<n;i++) {
660*db4deed7SKarl Rupp           if ((!used[i])&&(rhs[i]!=0.0)) {
661827bd09bSSatish Balay             ct++; nfo++;
662827bd09bSSatish Balay 
663c1235816SBarry Smith             if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n");
664827bd09bSSatish Balay 
665827bd09bSSatish Balay             *--iptr = local2global[i];
666827bd09bSSatish Balay             used[i]=edge;
667827bd09bSSatish Balay           }
668827bd09bSSatish Balay         }
669ca8e9878SJed Brown         if (ct>1) { PCTFS_ivec_sort(iptr,ct); }
670827bd09bSSatish Balay 
671827bd09bSSatish Balay         lnsep[edge]=ct;
67252f87cdaSBarry Smith         nsep[edge]= (PetscInt) rsum[1];
673827bd09bSSatish Balay         dir [edge]=RIGHT;
674827bd09bSSatish Balay       }
675827bd09bSSatish Balay 
676827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
677827bd09bSSatish Balay       /* do we need full set of separators for this?                */
678827bd09bSSatish Balay 
679827bd09bSSatish Balay       /* fold rhs hc into lower */
680*db4deed7SKarl Rupp       if (id>=mask) { id-=mask; }
681827bd09bSSatish Balay     }
682*db4deed7SKarl Rupp   } else {
683*db4deed7SKarl Rupp     for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
684827bd09bSSatish Balay       /* set rsh of hc, fire, and collect lhs responses */
685ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
686ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
687827bd09bSSatish Balay 
688827bd09bSSatish Balay       /* set lsh of hc, fire, and collect rhs responses */
689ca8e9878SJed Brown       (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
690ca8e9878SJed Brown       PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
691827bd09bSSatish Balay 
692827bd09bSSatish Balay       /* count number of dofs I own that have signal and not in sep set */
693*db4deed7SKarl Rupp       for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
694*db4deed7SKarl Rupp         if (!used[i]) {
695827bd09bSSatish Balay           /* number of unmarked dofs on node */
696827bd09bSSatish Balay           ct++;
697827bd09bSSatish Balay           /* number of dofs to be marked on lhs hc */
698827bd09bSSatish Balay           if ((id< mask)&&(lhs[i]!=0.0)) { sum[0]++; }
699827bd09bSSatish Balay           /* number of dofs to be marked on rhs hc */
700827bd09bSSatish Balay           if ((id>=mask)&&(rhs[i]!=0.0)) { sum[1]++; }
701827bd09bSSatish Balay         }
702827bd09bSSatish Balay       }
703827bd09bSSatish Balay 
704827bd09bSSatish Balay       /* go for load balance - choose half with most unmarked dofs, bias LHS */
705827bd09bSSatish Balay       (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
706b1c944f5SJed Brown       PCTFS_giop_hc(sum,w,4,op,edge);
707827bd09bSSatish Balay 
708827bd09bSSatish Balay       /* lhs hc wins */
709*db4deed7SKarl Rupp       if (sum[2]>=sum[3]) {
710*db4deed7SKarl Rupp         if (id<mask) {
711827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
712*db4deed7SKarl Rupp           for (ct=i=0;i<n;i++) {
713*db4deed7SKarl Rupp             if ((!used[i])&&(lhs[i]!=0.0)) {
714827bd09bSSatish Balay               ct++; nfo++;
715827bd09bSSatish Balay               *--iptr = local2global[i];
716827bd09bSSatish Balay               used[i]=edge;
717827bd09bSSatish Balay             }
718827bd09bSSatish Balay           }
719ca8e9878SJed Brown           if (ct>1) { PCTFS_ivec_sort(iptr,ct); }
720827bd09bSSatish Balay           lnsep[edge]=ct;
721827bd09bSSatish Balay         }
722827bd09bSSatish Balay         nsep[edge]=sum[0];
723827bd09bSSatish Balay         dir [edge]=LEFT;
724*db4deed7SKarl Rupp       } else { /* rhs hc wins */
725*db4deed7SKarl Rupp         if (id>=mask) {
726827bd09bSSatish Balay           /* mark dofs I own that have signal and not in sep set */
727*db4deed7SKarl Rupp           for (ct=i=0;i<n;i++) {
728*db4deed7SKarl Rupp             if ((!used[i])&&(rhs[i]!=0.0)) {
729827bd09bSSatish Balay               ct++; nfo++;
730827bd09bSSatish Balay               *--iptr = local2global[i];
731827bd09bSSatish Balay               used[i]=edge;
732827bd09bSSatish Balay             }
733827bd09bSSatish Balay           }
734ca8e9878SJed Brown           if (ct>1) {PCTFS_ivec_sort(iptr,ct);}
735827bd09bSSatish Balay           lnsep[edge]=ct;
736827bd09bSSatish Balay         }
737827bd09bSSatish Balay         nsep[edge]=sum[1];
738827bd09bSSatish Balay         dir [edge]=RIGHT;
739827bd09bSSatish Balay       }
740827bd09bSSatish Balay       /* LATER or we can recur on these to order seps at this level */
741827bd09bSSatish Balay       /* do we need full set of separators for this?                */
742827bd09bSSatish Balay 
743827bd09bSSatish Balay       /* fold rhs hc into lower */
744*db4deed7SKarl Rupp       if (id>=mask) { id-=mask; }
745827bd09bSSatish Balay     }
746827bd09bSSatish Balay   }
747827bd09bSSatish Balay 
748827bd09bSSatish Balay 
749827bd09bSSatish Balay   /* level 0 is on processor case - so mark the remainder */
750*db4deed7SKarl Rupp   for (ct=i=0;i<n;i++) {
751*db4deed7SKarl Rupp     if (!used[i]) {
752827bd09bSSatish Balay       ct++; nfo++;
753827bd09bSSatish Balay       *--iptr = local2global[i];
754827bd09bSSatish Balay       used[i]=edge;
755827bd09bSSatish Balay     }
756827bd09bSSatish Balay   }
757ca8e9878SJed Brown   if (ct>1) { PCTFS_ivec_sort(iptr,ct); }
758827bd09bSSatish Balay   lnsep[edge]=ct;
759827bd09bSSatish Balay   nsep [edge]=ct;
760827bd09bSSatish Balay   dir  [edge]=LEFT;
761827bd09bSSatish Balay 
762827bd09bSSatish Balay   xxt_handle->info->nsep=nsep;
763827bd09bSSatish Balay   xxt_handle->info->lnsep=lnsep;
764827bd09bSSatish Balay   xxt_handle->info->fo=fo;
765827bd09bSSatish Balay   xxt_handle->info->nfo=nfo;
766827bd09bSSatish Balay 
767a501084fSBarry Smith   free(dir);
768a501084fSBarry Smith   free(lhs);
769a501084fSBarry Smith   free(rhs);
770a501084fSBarry Smith   free(used);
7710924e98cSBarry Smith   PetscFunctionReturn(0);
772827bd09bSSatish Balay }
773827bd09bSSatish Balay 
7747b1ae94cSBarry Smith /**************************************xxt.c***********************************/
77552f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
776827bd09bSSatish Balay {
777827bd09bSSatish Balay   mv_info *mvi;
778827bd09bSSatish Balay 
779827bd09bSSatish Balay 
780a501084fSBarry Smith   mvi = (mv_info*)malloc(sizeof(mv_info));
781827bd09bSSatish Balay   mvi->n=n;
782827bd09bSSatish Balay   mvi->m=m;
783827bd09bSSatish Balay   mvi->n_global=-1;
784827bd09bSSatish Balay   mvi->m_global=-1;
78552f87cdaSBarry Smith   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
786ca8e9878SJed Brown   PCTFS_ivec_copy(mvi->local2global,local2global,m);
787827bd09bSSatish Balay   mvi->local2global[m] = INT_MAX;
788a501084fSBarry Smith   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
789827bd09bSSatish Balay   mvi->grid_data=grid_data;
790827bd09bSSatish Balay 
791827bd09bSSatish Balay   /* set xxt communication handle to perform restricted matvec */
792ca8e9878SJed Brown   mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
793827bd09bSSatish Balay 
794827bd09bSSatish Balay   return(mvi);
795827bd09bSSatish Balay }
796827bd09bSSatish Balay 
7977b1ae94cSBarry Smith /**************************************xxt.c***********************************/
7980924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
799827bd09bSSatish Balay {
8000924e98cSBarry Smith   PetscFunctionBegin;
801827bd09bSSatish Balay   A->matvec((mv_info*)A->grid_data,v,u);
8020924e98cSBarry Smith   PetscFunctionReturn(0);
803827bd09bSSatish Balay }
804827bd09bSSatish Balay 
805827bd09bSSatish Balay 
806827bd09bSSatish Balay 
807