xref: /petsc/src/ksp/pc/impls/tfs/gs.c (revision ca8e9878d215fbd0d6b1ef802734bc6b1f5fa99c)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /***********************************gs.c***************************************
3827bd09bSSatish Balay 
4827bd09bSSatish Balay Author: Henry M. Tufo III
5827bd09bSSatish Balay 
6827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
7827bd09bSSatish Balay 
8827bd09bSSatish Balay snail-mail:
9827bd09bSSatish Balay Division of Applied Mathematics
10827bd09bSSatish Balay Brown University
11827bd09bSSatish Balay Providence, RI 02912
12827bd09bSSatish Balay 
13827bd09bSSatish Balay Last Modification:
14827bd09bSSatish Balay 6.21.97
15827bd09bSSatish Balay ************************************gs.c**************************************/
16827bd09bSSatish Balay 
17827bd09bSSatish Balay /***********************************gs.c***************************************
18827bd09bSSatish Balay File Description:
19827bd09bSSatish Balay -----------------
20827bd09bSSatish Balay 
21827bd09bSSatish Balay ************************************gs.c**************************************/
22827bd09bSSatish Balay 
23c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
2439945688SSatish Balay 
25827bd09bSSatish Balay /* default length of number of items via tree - doubles if exceeded */
26827bd09bSSatish Balay #define TREE_BUF_SZ 2048;
27827bd09bSSatish Balay #define GS_VEC_SZ   1
28827bd09bSSatish Balay 
29827bd09bSSatish Balay 
30827bd09bSSatish Balay 
31827bd09bSSatish Balay /***********************************gs.c***************************************
32827bd09bSSatish Balay Type: struct gather_scatter_id
33827bd09bSSatish Balay ------------------------------
34827bd09bSSatish Balay 
35827bd09bSSatish Balay ************************************gs.c**************************************/
36827bd09bSSatish Balay typedef struct gather_scatter_id {
3752f87cdaSBarry Smith   PetscInt id;
3852f87cdaSBarry Smith   PetscInt nel_min;
3952f87cdaSBarry Smith   PetscInt nel_max;
4052f87cdaSBarry Smith   PetscInt nel_sum;
4152f87cdaSBarry Smith   PetscInt negl;
4252f87cdaSBarry Smith   PetscInt gl_max;
4352f87cdaSBarry Smith   PetscInt gl_min;
4452f87cdaSBarry Smith   PetscInt repeats;
4552f87cdaSBarry Smith   PetscInt ordered;
4652f87cdaSBarry Smith   PetscInt positive;
47a501084fSBarry Smith   PetscScalar *vals;
48827bd09bSSatish Balay 
49827bd09bSSatish Balay   /* bit mask info */
5052f87cdaSBarry Smith   PetscInt *my_proc_mask;
5152f87cdaSBarry Smith   PetscInt mask_sz;
5252f87cdaSBarry Smith   PetscInt *ngh_buf;
5352f87cdaSBarry Smith   PetscInt ngh_buf_sz;
5452f87cdaSBarry Smith   PetscInt *nghs;
5552f87cdaSBarry Smith   PetscInt num_nghs;
5652f87cdaSBarry Smith   PetscInt max_nghs;
5752f87cdaSBarry Smith   PetscInt *pw_nghs;
5852f87cdaSBarry Smith   PetscInt num_pw_nghs;
5952f87cdaSBarry Smith   PetscInt *tree_nghs;
6052f87cdaSBarry Smith   PetscInt num_tree_nghs;
61827bd09bSSatish Balay 
6252f87cdaSBarry Smith   PetscInt num_loads;
63827bd09bSSatish Balay 
64827bd09bSSatish Balay   /* repeats == true -> local info */
6552f87cdaSBarry Smith   PetscInt nel;         /* number of unique elememts */
6652f87cdaSBarry Smith   PetscInt *elms;       /* of size nel */
6752f87cdaSBarry Smith   PetscInt nel_total;
6852f87cdaSBarry Smith   PetscInt *local_elms; /* of size nel_total */
6952f87cdaSBarry Smith   PetscInt *companion;  /* of size nel_total */
70827bd09bSSatish Balay 
71827bd09bSSatish Balay   /* local info */
7252f87cdaSBarry Smith   PetscInt num_local_total;
7352f87cdaSBarry Smith   PetscInt local_strength;
7452f87cdaSBarry Smith   PetscInt num_local;
7552f87cdaSBarry Smith   PetscInt *num_local_reduce;
7652f87cdaSBarry Smith   PetscInt **local_reduce;
7752f87cdaSBarry Smith   PetscInt num_local_gop;
7852f87cdaSBarry Smith   PetscInt *num_gop_local_reduce;
7952f87cdaSBarry Smith   PetscInt **gop_local_reduce;
80827bd09bSSatish Balay 
81827bd09bSSatish Balay   /* pairwise info */
8252f87cdaSBarry Smith   PetscInt level;
8352f87cdaSBarry Smith   PetscInt num_pairs;
8452f87cdaSBarry Smith   PetscInt max_pairs;
8552f87cdaSBarry Smith   PetscInt loc_node_pairs;
8652f87cdaSBarry Smith   PetscInt max_node_pairs;
8752f87cdaSBarry Smith   PetscInt min_node_pairs;
8852f87cdaSBarry Smith   PetscInt avg_node_pairs;
8952f87cdaSBarry Smith   PetscInt *pair_list;
9052f87cdaSBarry Smith   PetscInt *msg_sizes;
9152f87cdaSBarry Smith   PetscInt **node_list;
9252f87cdaSBarry Smith   PetscInt len_pw_list;
9352f87cdaSBarry Smith   PetscInt *pw_elm_list;
94a501084fSBarry Smith   PetscScalar *pw_vals;
95827bd09bSSatish Balay 
96827bd09bSSatish Balay   MPI_Request *msg_ids_in;
97827bd09bSSatish Balay   MPI_Request *msg_ids_out;
98827bd09bSSatish Balay 
99a501084fSBarry Smith   PetscScalar *out;
100a501084fSBarry Smith   PetscScalar *in;
10152f87cdaSBarry Smith   PetscInt msg_total;
102827bd09bSSatish Balay 
103827bd09bSSatish Balay   /* tree - crystal accumulator info */
10452f87cdaSBarry Smith   PetscInt max_left_over;
10552f87cdaSBarry Smith   PetscInt *pre;
10652f87cdaSBarry Smith   PetscInt *in_num;
10752f87cdaSBarry Smith   PetscInt *out_num;
10852f87cdaSBarry Smith   PetscInt **in_list;
10952f87cdaSBarry Smith   PetscInt **out_list;
110827bd09bSSatish Balay 
111827bd09bSSatish Balay   /* new tree work*/
11252f87cdaSBarry Smith   PetscInt  tree_nel;
11352f87cdaSBarry Smith   PetscInt *tree_elms;
114a501084fSBarry Smith   PetscScalar *tree_buf;
115a501084fSBarry Smith   PetscScalar *tree_work;
116827bd09bSSatish Balay 
11752f87cdaSBarry Smith   PetscInt  tree_map_sz;
11852f87cdaSBarry Smith   PetscInt *tree_map_in;
11952f87cdaSBarry Smith   PetscInt *tree_map_out;
120827bd09bSSatish Balay 
121827bd09bSSatish Balay   /* current memory status */
12252f87cdaSBarry Smith   PetscInt gl_bss_min;
12352f87cdaSBarry Smith   PetscInt gl_perm_min;
124827bd09bSSatish Balay 
125*ca8e9878SJed Brown   /* max segment size for PCTFS_gs_gop_vec() */
12652f87cdaSBarry Smith   PetscInt vec_sz;
127827bd09bSSatish Balay 
128827bd09bSSatish Balay   /* hack to make paul happy */
129*ca8e9878SJed Brown   MPI_Comm PCTFS_gs_comm;
130827bd09bSSatish Balay 
131*ca8e9878SJed Brown } PCTFS_gs_id;
132827bd09bSSatish Balay 
133*ca8e9878SJed Brown static PCTFS_gs_id *gsi_check_args(PetscInt *elms, PetscInt nel, PetscInt level);
134*ca8e9878SJed Brown static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs);
135*ca8e9878SJed Brown static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs);
136*ca8e9878SJed Brown static PetscErrorCode set_pairwise(PCTFS_gs_id *gs);
137*ca8e9878SJed Brown static PCTFS_gs_id * gsi_new(void);
138*ca8e9878SJed Brown static PetscErrorCode set_tree(PCTFS_gs_id *gs);
139827bd09bSSatish Balay 
140827bd09bSSatish Balay /* same for all but vector flavor */
141*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_out(PCTFS_gs_id *gs, PetscScalar *vals);
142827bd09bSSatish Balay /* vector flavor */
143*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_out(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step);
144827bd09bSSatish Balay 
145*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step);
146*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step);
147*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step);
148*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step);
149*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_tree_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step);
150827bd09bSSatish Balay 
151827bd09bSSatish Balay 
152*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_plus(PCTFS_gs_id *gs, PetscScalar *vals);
153*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals);
154827bd09bSSatish Balay 
155*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim);
156*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim);
157*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim);
158827bd09bSSatish Balay 
159827bd09bSSatish Balay /* global vars */
160827bd09bSSatish Balay /* from comm.c module */
161827bd09bSSatish Balay 
16252f87cdaSBarry Smith static PetscInt num_gs_ids = 0;
163827bd09bSSatish Balay 
164827bd09bSSatish Balay /* should make this dynamic ... later */
16552f87cdaSBarry Smith static PetscInt msg_buf=MAX_MSG_BUF;
16652f87cdaSBarry Smith static PetscInt vec_sz=GS_VEC_SZ;
16752f87cdaSBarry Smith static PetscInt *tree_buf=NULL;
16852f87cdaSBarry Smith static PetscInt tree_buf_sz=0;
16952f87cdaSBarry Smith static PetscInt ntree=0;
170827bd09bSSatish Balay 
171f1ed62a8SBarry Smith /***************************************************************************/
172*ca8e9878SJed Brown PetscErrorCode PCTFS_gs_init_vec_sz(PetscInt size)
173827bd09bSSatish Balay {
1743fdc5746SBarry Smith   PetscFunctionBegin;
175827bd09bSSatish Balay   vec_sz = size;
1763fdc5746SBarry Smith   PetscFunctionReturn(0);
177827bd09bSSatish Balay }
178827bd09bSSatish Balay 
179f1ed62a8SBarry Smith /******************************************************************************/
180*ca8e9878SJed Brown PetscErrorCode PCTFS_gs_init_msg_buf_sz(PetscInt buf_size)
181827bd09bSSatish Balay {
1823fdc5746SBarry Smith   PetscFunctionBegin;
183827bd09bSSatish Balay   msg_buf = buf_size;
1843fdc5746SBarry Smith   PetscFunctionReturn(0);
185827bd09bSSatish Balay }
186827bd09bSSatish Balay 
187f1ed62a8SBarry Smith /******************************************************************************/
188*ca8e9878SJed Brown PCTFS_gs_id *PCTFS_gs_init( PetscInt *elms, PetscInt nel, PetscInt level)
189827bd09bSSatish Balay {
190*ca8e9878SJed Brown    PCTFS_gs_id *gs;
191*ca8e9878SJed Brown   MPI_Group PCTFS_gs_group;
192*ca8e9878SJed Brown   MPI_Comm  PCTFS_gs_comm;
193f1ed62a8SBarry Smith   PetscErrorCode ierr;
194827bd09bSSatish Balay 
1953fdc5746SBarry Smith   PetscFunctionBegin;
196827bd09bSSatish Balay   /* ensure that communication package has been initialized */
197b1c944f5SJed Brown   PCTFS_comm_init();
198827bd09bSSatish Balay 
199827bd09bSSatish Balay 
200827bd09bSSatish Balay   /* determines if we have enough dynamic/semi-static memory */
201827bd09bSSatish Balay   /* checks input, allocs and sets gd_id template            */
202827bd09bSSatish Balay   gs = gsi_check_args(elms,nel,level);
203827bd09bSSatish Balay 
204827bd09bSSatish Balay   /* only bit mask version up and working for the moment    */
205827bd09bSSatish Balay   /* LATER :: get int list version working for sparse pblms */
206f1ed62a8SBarry Smith   ierr = gsi_via_bit_mask(gs);CHKERRABORT(PETSC_COMM_WORLD,ierr);
207827bd09bSSatish Balay 
208827bd09bSSatish Balay 
209*ca8e9878SJed Brown   ierr = MPI_Comm_group(MPI_COMM_WORLD,&PCTFS_gs_group);CHKERRABORT(PETSC_COMM_WORLD,ierr);
210*ca8e9878SJed Brown   ierr = MPI_Comm_create(MPI_COMM_WORLD,PCTFS_gs_group,&PCTFS_gs_comm);CHKERRABORT(PETSC_COMM_WORLD,ierr);
211*ca8e9878SJed Brown   gs->PCTFS_gs_comm=PCTFS_gs_comm;
212827bd09bSSatish Balay 
213827bd09bSSatish Balay   return(gs);
214827bd09bSSatish Balay }
215827bd09bSSatish Balay 
216f1ed62a8SBarry Smith /******************************************************************************/
217*ca8e9878SJed Brown static PCTFS_gs_id *gsi_new(void)
218827bd09bSSatish Balay {
219f1ed62a8SBarry Smith   PetscErrorCode ierr;
220*ca8e9878SJed Brown   PCTFS_gs_id *gs;
221*ca8e9878SJed Brown   gs = (PCTFS_gs_id *) malloc(sizeof(PCTFS_gs_id));
222*ca8e9878SJed Brown   ierr = PetscMemzero(gs,sizeof(PCTFS_gs_id));CHKERRABORT(PETSC_COMM_WORLD,ierr);
223827bd09bSSatish Balay   return(gs);
224827bd09bSSatish Balay }
225827bd09bSSatish Balay 
226f1ed62a8SBarry Smith /******************************************************************************/
227*ca8e9878SJed Brown static PCTFS_gs_id * gsi_check_args(PetscInt *in_elms, PetscInt nel, PetscInt level)
228827bd09bSSatish Balay {
22952f87cdaSBarry Smith    PetscInt i, j, k, t2;
23052f87cdaSBarry Smith   PetscInt *companion, *elms, *unique, *iptr;
23152f87cdaSBarry Smith   PetscInt num_local=0, *num_to_reduce, **local_reduce;
23252f87cdaSBarry Smith   PetscInt oprs[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_MIN,GL_B_AND};
23352f87cdaSBarry Smith   PetscInt vals[sizeof(oprs)/sizeof(oprs[0])-1];
23452f87cdaSBarry Smith   PetscInt work[sizeof(oprs)/sizeof(oprs[0])-1];
235*ca8e9878SJed Brown   PCTFS_gs_id *gs;
236d1528f56SBarry Smith   PetscErrorCode ierr;
237827bd09bSSatish Balay 
238827bd09bSSatish Balay 
239c1235816SBarry Smith   if (!in_elms) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"elms point to nothing!!!\n");
240c1235816SBarry Smith   if (nel<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"can't have fewer than 0 elms!!!\n");
241827bd09bSSatish Balay 
242827bd09bSSatish Balay   if (nel==0)
243f1ed62a8SBarry Smith     {ierr = PetscInfo(0,"I don't have any elements!!!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);}
244827bd09bSSatish Balay 
245827bd09bSSatish Balay   /* get space for gs template */
246827bd09bSSatish Balay   gs = gsi_new();
247827bd09bSSatish Balay   gs->id = ++num_gs_ids;
248827bd09bSSatish Balay 
249827bd09bSSatish Balay   /* hmt 6.4.99                                            */
250827bd09bSSatish Balay   /* caller can set global ids that don't participate to 0 */
251*ca8e9878SJed Brown   /* PCTFS_gs_init ignores all zeros in elm list                 */
252827bd09bSSatish Balay   /* negative global ids are still invalid                 */
253827bd09bSSatish Balay   for (i=j=0;i<nel;i++)
254827bd09bSSatish Balay     {if (in_elms[i]!=0) {j++;}}
255827bd09bSSatish Balay 
256827bd09bSSatish Balay   k=nel; nel=j;
257827bd09bSSatish Balay 
258827bd09bSSatish Balay   /* copy over in_elms list and create inverse map */
25952f87cdaSBarry Smith   elms = (PetscInt*) malloc((nel+1)*sizeof(PetscInt));
26052f87cdaSBarry Smith   companion = (PetscInt*) malloc(nel*sizeof(PetscInt));
2611d7d0905SBarry Smith 
262827bd09bSSatish Balay   for (i=j=0;i<k;i++)
263827bd09bSSatish Balay     {
264827bd09bSSatish Balay       if (in_elms[i]!=0)
265827bd09bSSatish Balay         {elms[j] = in_elms[i]; companion[j++] = i;}
266827bd09bSSatish Balay     }
267827bd09bSSatish Balay 
268c1235816SBarry Smith   if (j!=nel) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"nel j mismatch!\n");
269827bd09bSSatish Balay 
270827bd09bSSatish Balay   /* pre-pass ... check to see if sorted */
271827bd09bSSatish Balay   elms[nel] = INT_MAX;
272827bd09bSSatish Balay   iptr = elms;
273827bd09bSSatish Balay   unique = elms+1;
274827bd09bSSatish Balay   j=0;
275827bd09bSSatish Balay   while (*iptr!=INT_MAX)
276827bd09bSSatish Balay     {
277827bd09bSSatish Balay       if (*iptr++>*unique++)
278827bd09bSSatish Balay         {j=1; break;}
279827bd09bSSatish Balay     }
280827bd09bSSatish Balay 
281827bd09bSSatish Balay   /* set up inverse map */
282827bd09bSSatish Balay   if (j)
283827bd09bSSatish Balay     {
284f1ed62a8SBarry Smith       ierr = PetscInfo(0,"gsi_check_args() :: elm list *not* sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);
285*ca8e9878SJed Brown       ierr = PCTFS_SMI_sort((void*)elms, (void*)companion, nel, SORT_INTEGER);CHKERRABORT(PETSC_COMM_WORLD,ierr);
286827bd09bSSatish Balay     }
287827bd09bSSatish Balay   else
288f1ed62a8SBarry Smith     {ierr = PetscInfo(0,"gsi_check_args() :: elm list sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);}
289827bd09bSSatish Balay   elms[nel] = INT_MIN;
290827bd09bSSatish Balay 
291827bd09bSSatish Balay   /* first pass */
292827bd09bSSatish Balay   /* determine number of unique elements, check pd */
293827bd09bSSatish Balay   for (i=k=0;i<nel;i+=j)
294827bd09bSSatish Balay     {
295827bd09bSSatish Balay       t2 = elms[i];
296827bd09bSSatish Balay       j=++i;
297827bd09bSSatish Balay 
298827bd09bSSatish Balay       /* clump 'em for now */
299827bd09bSSatish Balay       while (elms[j]==t2) {j++;}
300827bd09bSSatish Balay 
301827bd09bSSatish Balay       /* how many together and num local */
302827bd09bSSatish Balay       if (j-=i)
303827bd09bSSatish Balay         {num_local++; k+=j;}
304827bd09bSSatish Balay     }
305827bd09bSSatish Balay 
306827bd09bSSatish Balay   /* how many unique elements? */
307827bd09bSSatish Balay   gs->repeats=k;
308827bd09bSSatish Balay   gs->nel = nel-k;
309827bd09bSSatish Balay 
310827bd09bSSatish Balay 
311827bd09bSSatish Balay   /* number of repeats? */
312827bd09bSSatish Balay   gs->num_local = num_local;
313827bd09bSSatish Balay   num_local+=2;
31452f87cdaSBarry Smith   gs->local_reduce=local_reduce=(PetscInt **)malloc(num_local*sizeof(PetscInt*));
31552f87cdaSBarry Smith   gs->num_local_reduce=num_to_reduce=(PetscInt*) malloc(num_local*sizeof(PetscInt));
316827bd09bSSatish Balay 
31752f87cdaSBarry Smith   unique = (PetscInt*) malloc((gs->nel+1)*sizeof(PetscInt));
318827bd09bSSatish Balay   gs->elms = unique;
319827bd09bSSatish Balay   gs->nel_total = nel;
320827bd09bSSatish Balay   gs->local_elms = elms;
321827bd09bSSatish Balay   gs->companion = companion;
322827bd09bSSatish Balay 
323827bd09bSSatish Balay   /* compess map as well as keep track of local ops */
324827bd09bSSatish Balay   for (num_local=i=j=0;i<gs->nel;i++)
325827bd09bSSatish Balay     {
326827bd09bSSatish Balay       k=j;
327827bd09bSSatish Balay       t2 = unique[i] = elms[j];
328827bd09bSSatish Balay       companion[i] = companion[j];
329827bd09bSSatish Balay 
330827bd09bSSatish Balay       while (elms[j]==t2) {j++;}
331827bd09bSSatish Balay 
332827bd09bSSatish Balay       if ((t2=(j-k))>1)
333827bd09bSSatish Balay         {
334827bd09bSSatish Balay           /* number together */
335827bd09bSSatish Balay           num_to_reduce[num_local] = t2++;
33652f87cdaSBarry Smith           iptr = local_reduce[num_local++] = (PetscInt*)malloc(t2*sizeof(PetscInt));
337827bd09bSSatish Balay 
338827bd09bSSatish Balay           /* to use binary searching don't remap until we check intersection */
339827bd09bSSatish Balay           *iptr++ = i;
340827bd09bSSatish Balay 
341827bd09bSSatish Balay           /* note that we're skipping the first one */
342827bd09bSSatish Balay           while (++k<j)
343827bd09bSSatish Balay             {*(iptr++) = companion[k];}
344827bd09bSSatish Balay           *iptr = -1;
345827bd09bSSatish Balay         }
346827bd09bSSatish Balay     }
347827bd09bSSatish Balay 
348827bd09bSSatish Balay   /* sentinel for ngh_buf */
349827bd09bSSatish Balay   unique[gs->nel]=INT_MAX;
350827bd09bSSatish Balay 
351827bd09bSSatish Balay   /* for two partition sort hack */
352827bd09bSSatish Balay   num_to_reduce[num_local] = 0;
353827bd09bSSatish Balay   local_reduce[num_local] = NULL;
354827bd09bSSatish Balay   num_to_reduce[++num_local] = 0;
355827bd09bSSatish Balay   local_reduce[num_local] = NULL;
356827bd09bSSatish Balay 
357827bd09bSSatish Balay   /* load 'em up */
358827bd09bSSatish Balay   /* note one extra to hold NON_UNIFORM flag!!! */
359827bd09bSSatish Balay   vals[2] = vals[1] = vals[0] = nel;
360827bd09bSSatish Balay   if (gs->nel>0)
361827bd09bSSatish Balay     {
3621d7d0905SBarry Smith        vals[3] = unique[0];
3631d7d0905SBarry Smith        vals[4] = unique[gs->nel-1];
364827bd09bSSatish Balay     }
365827bd09bSSatish Balay   else
366827bd09bSSatish Balay     {
3671d7d0905SBarry Smith        vals[3] = INT_MAX;
3681d7d0905SBarry Smith        vals[4] = INT_MIN;
369827bd09bSSatish Balay     }
370827bd09bSSatish Balay   vals[5] = level;
371827bd09bSSatish Balay   vals[6] = num_gs_ids;
372827bd09bSSatish Balay 
373827bd09bSSatish Balay   /* GLOBAL: send 'em out */
374b1c944f5SJed Brown   ierr = PCTFS_giop(vals,work,sizeof(oprs)/sizeof(oprs[0])-1,oprs);CHKERRABORT(PETSC_COMM_WORLD,ierr);
375827bd09bSSatish Balay 
376827bd09bSSatish Balay   /* must be semi-pos def - only pairwise depends on this */
377827bd09bSSatish Balay   /* LATER - remove this restriction */
378c1235816SBarry Smith   if (vals[3]<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system not semi-pos def \n");
379c1235816SBarry Smith   if (vals[4]==INT_MAX) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system ub too large !\n");
380827bd09bSSatish Balay 
381827bd09bSSatish Balay   gs->nel_min = vals[0];
382827bd09bSSatish Balay   gs->nel_max = vals[1];
383827bd09bSSatish Balay   gs->nel_sum = vals[2];
384827bd09bSSatish Balay   gs->gl_min  = vals[3];
385827bd09bSSatish Balay   gs->gl_max  = vals[4];
386827bd09bSSatish Balay   gs->negl    = vals[4]-vals[3]+1;
387827bd09bSSatish Balay 
388c1235816SBarry Smith   if (gs->negl<=0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system empty or neg :: %d\n");
389827bd09bSSatish Balay 
390827bd09bSSatish Balay   /* LATER :: add level == -1 -> program selects level */
391827bd09bSSatish Balay   if (vals[5]<0)
392827bd09bSSatish Balay     {vals[5]=0;}
393b1c944f5SJed Brown   else if (vals[5]>PCTFS_num_nodes)
394b1c944f5SJed Brown     {vals[5]=PCTFS_num_nodes;}
395827bd09bSSatish Balay   gs->level = vals[5];
396827bd09bSSatish Balay 
397827bd09bSSatish Balay   return(gs);
398827bd09bSSatish Balay }
399827bd09bSSatish Balay 
400f1ed62a8SBarry Smith /******************************************************************************/
401*ca8e9878SJed Brown static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs)
402827bd09bSSatish Balay {
40352f87cdaSBarry Smith    PetscInt i, nel, *elms;
40452f87cdaSBarry Smith   PetscInt t1;
40552f87cdaSBarry Smith   PetscInt **reduce;
40652f87cdaSBarry Smith   PetscInt *map;
407f1ed62a8SBarry Smith   PetscErrorCode ierr;
408827bd09bSSatish Balay 
409f1ed62a8SBarry Smith   PetscFunctionBegin;
410*ca8e9878SJed Brown   /* totally local removes ... PCTFS_ct_bits == 0 */
411827bd09bSSatish Balay   get_ngh_buf(gs);
412827bd09bSSatish Balay 
413827bd09bSSatish Balay   if (gs->level)
414827bd09bSSatish Balay     {set_pairwise(gs);}
415827bd09bSSatish Balay 
416827bd09bSSatish Balay   if (gs->max_left_over)
417827bd09bSSatish Balay     {set_tree(gs);}
418827bd09bSSatish Balay 
419827bd09bSSatish Balay   /* intersection local and pairwise/tree? */
420827bd09bSSatish Balay   gs->num_local_total = gs->num_local;
421827bd09bSSatish Balay   gs->gop_local_reduce = gs->local_reduce;
422827bd09bSSatish Balay   gs->num_gop_local_reduce = gs->num_local_reduce;
423827bd09bSSatish Balay 
424827bd09bSSatish Balay   map = gs->companion;
425827bd09bSSatish Balay 
426827bd09bSSatish Balay   /* is there any local compression */
427d890fc11SSatish Balay   if (!gs->num_local) {
428827bd09bSSatish Balay     gs->local_strength = NONE;
429827bd09bSSatish Balay     gs->num_local_gop = 0;
430d890fc11SSatish Balay   } else {
431827bd09bSSatish Balay       /* ok find intersection */
432827bd09bSSatish Balay       map = gs->companion;
433827bd09bSSatish Balay       reduce = gs->local_reduce;
434827bd09bSSatish Balay       for (i=0, t1=0; i<gs->num_local; i++, reduce++)
435827bd09bSSatish Balay         {
436*ca8e9878SJed Brown           if ((PCTFS_ivec_binary_search(**reduce,gs->pw_elm_list,gs->len_pw_list)>=0)
437827bd09bSSatish Balay               ||
438*ca8e9878SJed Brown               PCTFS_ivec_binary_search(**reduce,gs->tree_map_in,gs->tree_map_sz)>=0)
439827bd09bSSatish Balay             {
440827bd09bSSatish Balay               t1++;
441e32f2f54SBarry Smith               if (gs->num_local_reduce[i]<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nobody in list?");
442827bd09bSSatish Balay               gs->num_local_reduce[i] *= -1;
443827bd09bSSatish Balay             }
444827bd09bSSatish Balay            **reduce=map[**reduce];
445827bd09bSSatish Balay         }
446827bd09bSSatish Balay 
447827bd09bSSatish Balay       /* intersection is empty */
448827bd09bSSatish Balay       if (!t1)
449827bd09bSSatish Balay         {
450827bd09bSSatish Balay           gs->local_strength = FULL;
451827bd09bSSatish Balay           gs->num_local_gop = 0;
452827bd09bSSatish Balay         }
453827bd09bSSatish Balay       /* intersection not empty */
454827bd09bSSatish Balay       else
455827bd09bSSatish Balay         {
456827bd09bSSatish Balay           gs->local_strength = PARTIAL;
457*ca8e9878SJed Brown           ierr = PCTFS_SMI_sort((void*)gs->num_local_reduce, (void*)gs->local_reduce, gs->num_local + 1, SORT_INT_PTR);CHKERRQ(ierr);
458827bd09bSSatish Balay 
459827bd09bSSatish Balay           gs->num_local_gop = t1;
460827bd09bSSatish Balay           gs->num_local_total =  gs->num_local;
461827bd09bSSatish Balay           gs->num_local    -= t1;
462827bd09bSSatish Balay           gs->gop_local_reduce = gs->local_reduce;
463827bd09bSSatish Balay           gs->num_gop_local_reduce = gs->num_local_reduce;
464827bd09bSSatish Balay 
465827bd09bSSatish Balay           for (i=0; i<t1; i++)
466827bd09bSSatish Balay             {
467e32f2f54SBarry Smith               if (gs->num_gop_local_reduce[i]>=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"they aren't negative?");
468827bd09bSSatish Balay               gs->num_gop_local_reduce[i] *= -1;
469827bd09bSSatish Balay               gs->local_reduce++;
470827bd09bSSatish Balay               gs->num_local_reduce++;
471827bd09bSSatish Balay             }
472827bd09bSSatish Balay           gs->local_reduce++;
473827bd09bSSatish Balay           gs->num_local_reduce++;
474827bd09bSSatish Balay         }
475827bd09bSSatish Balay     }
476827bd09bSSatish Balay 
477827bd09bSSatish Balay   elms = gs->pw_elm_list;
478827bd09bSSatish Balay   nel  = gs->len_pw_list;
479827bd09bSSatish Balay   for (i=0; i<nel; i++)
480827bd09bSSatish Balay     {elms[i] = map[elms[i]];}
481827bd09bSSatish Balay 
482827bd09bSSatish Balay   elms = gs->tree_map_in;
483827bd09bSSatish Balay   nel  = gs->tree_map_sz;
484827bd09bSSatish Balay   for (i=0; i<nel; i++)
485827bd09bSSatish Balay     {elms[i] = map[elms[i]];}
486827bd09bSSatish Balay 
487827bd09bSSatish Balay   /* clean up */
488a501084fSBarry Smith   free((void*) gs->local_elms);
489a501084fSBarry Smith   free((void*) gs->companion);
490a501084fSBarry Smith   free((void*) gs->elms);
491a501084fSBarry Smith   free((void*) gs->ngh_buf);
492827bd09bSSatish Balay   gs->local_elms = gs->companion = gs->elms = gs->ngh_buf = NULL;
4933fdc5746SBarry Smith   PetscFunctionReturn(0);
494827bd09bSSatish Balay }
495827bd09bSSatish Balay 
496f1ed62a8SBarry Smith /******************************************************************************/
49752f87cdaSBarry Smith static PetscErrorCode place_in_tree( PetscInt elm)
498827bd09bSSatish Balay {
49952f87cdaSBarry Smith    PetscInt *tp, n;
500827bd09bSSatish Balay 
5013fdc5746SBarry Smith   PetscFunctionBegin;
502827bd09bSSatish Balay   if (ntree==tree_buf_sz)
503827bd09bSSatish Balay     {
504827bd09bSSatish Balay       if (tree_buf_sz)
505827bd09bSSatish Balay         {
506827bd09bSSatish Balay           tp = tree_buf;
507827bd09bSSatish Balay           n = tree_buf_sz;
508827bd09bSSatish Balay           tree_buf_sz<<=1;
50952f87cdaSBarry Smith           tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt));
510*ca8e9878SJed Brown           PCTFS_ivec_copy(tree_buf,tp,n);
511a501084fSBarry Smith           free(tp);
512827bd09bSSatish Balay         }
513827bd09bSSatish Balay       else
514827bd09bSSatish Balay         {
515827bd09bSSatish Balay           tree_buf_sz = TREE_BUF_SZ;
51652f87cdaSBarry Smith           tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt));
517827bd09bSSatish Balay         }
518827bd09bSSatish Balay     }
519827bd09bSSatish Balay 
520827bd09bSSatish Balay   tree_buf[ntree++] = elm;
5213fdc5746SBarry Smith   PetscFunctionReturn(0);
522827bd09bSSatish Balay }
523827bd09bSSatish Balay 
524f1ed62a8SBarry Smith /******************************************************************************/
525*ca8e9878SJed Brown static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs)
526827bd09bSSatish Balay {
52752f87cdaSBarry Smith    PetscInt i, j, npw=0, ntree_map=0;
52852f87cdaSBarry Smith   PetscInt p_mask_size, ngh_buf_size, buf_size;
52952f87cdaSBarry Smith   PetscInt *p_mask, *sh_proc_mask, *pw_sh_proc_mask;
53052f87cdaSBarry Smith   PetscInt *ngh_buf, *buf1, *buf2;
53152f87cdaSBarry Smith   PetscInt offset, per_load, num_loads, or_ct, start, end;
53252f87cdaSBarry Smith   PetscInt *ptr1, *ptr2, i_start, negl, nel, *elms;
53352f87cdaSBarry Smith   PetscInt oper=GL_B_OR;
53452f87cdaSBarry Smith   PetscInt *ptr3, *t_mask, level, ct1, ct2;
535f1ed62a8SBarry Smith   PetscErrorCode ierr;
536827bd09bSSatish Balay 
5373fdc5746SBarry Smith   PetscFunctionBegin;
538827bd09bSSatish Balay   /* to make life easier */
539827bd09bSSatish Balay   nel   = gs->nel;
540827bd09bSSatish Balay   elms  = gs->elms;
541827bd09bSSatish Balay   level = gs->level;
542827bd09bSSatish Balay 
543b1c944f5SJed Brown   /* det #bytes needed for processor bit masks and init w/mask cor. to PCTFS_my_id */
544*ca8e9878SJed Brown   p_mask = (PetscInt*) malloc(p_mask_size=PCTFS_len_bit_mask(PCTFS_num_nodes));
545*ca8e9878SJed Brown   ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr);
546827bd09bSSatish Balay 
547827bd09bSSatish Balay   /* allocate space for masks and info bufs */
54852f87cdaSBarry Smith   gs->nghs = sh_proc_mask = (PetscInt*) malloc(p_mask_size);
54952f87cdaSBarry Smith   gs->pw_nghs = pw_sh_proc_mask = (PetscInt*) malloc(p_mask_size);
550827bd09bSSatish Balay   gs->ngh_buf_sz = ngh_buf_size = p_mask_size*nel;
55152f87cdaSBarry Smith   t_mask = (PetscInt*) malloc(p_mask_size);
55252f87cdaSBarry Smith   gs->ngh_buf = ngh_buf = (PetscInt*) malloc(ngh_buf_size);
553827bd09bSSatish Balay 
554827bd09bSSatish Balay   /* comm buffer size ... memory usage bounded by ~2*msg_buf */
555827bd09bSSatish Balay   /* had thought I could exploit rendezvous threshold */
556827bd09bSSatish Balay 
557827bd09bSSatish Balay   /* default is one pass */
558827bd09bSSatish Balay   per_load = negl  = gs->negl;
559827bd09bSSatish Balay   gs->num_loads = num_loads = 1;
560827bd09bSSatish Balay   i=p_mask_size*negl;
561827bd09bSSatish Balay 
562827bd09bSSatish Balay   /* possible overflow on buffer size */
563827bd09bSSatish Balay   /* overflow hack                    */
564827bd09bSSatish Balay   if (i<0) {i=INT_MAX;}
565827bd09bSSatish Balay 
56639945688SSatish Balay   buf_size = PetscMin(msg_buf,i);
567827bd09bSSatish Balay 
568827bd09bSSatish Balay   /* can we do it? */
569e32f2f54SBarry Smith   if (p_mask_size>buf_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"get_ngh_buf() :: buf<pms :: %d>%d\n",p_mask_size,buf_size);
570827bd09bSSatish Balay 
571b1c944f5SJed Brown   /* get PCTFS_giop buf space ... make *only* one malloc */
57252f87cdaSBarry Smith   buf1 = (PetscInt*) malloc(buf_size<<1);
573827bd09bSSatish Balay 
574827bd09bSSatish Balay   /* more than one gior exchange needed? */
575827bd09bSSatish Balay   if (buf_size!=i)
576827bd09bSSatish Balay     {
577827bd09bSSatish Balay       per_load = buf_size/p_mask_size;
578827bd09bSSatish Balay       buf_size = per_load*p_mask_size;
579827bd09bSSatish Balay       gs->num_loads = num_loads = negl/per_load + (negl%per_load>0);
580827bd09bSSatish Balay     }
581827bd09bSSatish Balay 
582827bd09bSSatish Balay 
583827bd09bSSatish Balay   /* convert buf sizes from #bytes to #ints - 32 bit only! */
584a501084fSBarry Smith   p_mask_size/=sizeof(PetscInt); ngh_buf_size/=sizeof(PetscInt); buf_size/=sizeof(PetscInt);
585827bd09bSSatish Balay 
586b1c944f5SJed Brown   /* find PCTFS_giop work space */
587827bd09bSSatish Balay   buf2 = buf1+buf_size;
588827bd09bSSatish Balay 
589827bd09bSSatish Balay   /* hold #ints needed for processor masks */
590827bd09bSSatish Balay   gs->mask_sz=p_mask_size;
591827bd09bSSatish Balay 
592827bd09bSSatish Balay   /* init buffers */
593*ca8e9878SJed Brown   ierr = PCTFS_ivec_zero(sh_proc_mask,p_mask_size);CHKERRQ(ierr);
594*ca8e9878SJed Brown   ierr = PCTFS_ivec_zero(pw_sh_proc_mask,p_mask_size);CHKERRQ(ierr);
595*ca8e9878SJed Brown   ierr = PCTFS_ivec_zero(ngh_buf,ngh_buf_size);CHKERRQ(ierr);
596827bd09bSSatish Balay 
597827bd09bSSatish Balay   /* HACK reset tree info */
598827bd09bSSatish Balay   tree_buf=NULL;
599827bd09bSSatish Balay   tree_buf_sz=ntree=0;
600827bd09bSSatish Balay 
601827bd09bSSatish Balay   /* ok do it */
602827bd09bSSatish Balay   for (ptr1=ngh_buf,ptr2=elms,end=gs->gl_min,or_ct=i=0; or_ct<num_loads; or_ct++)
603827bd09bSSatish Balay     {
604827bd09bSSatish Balay       /* identity for bitwise or is 000...000 */
605*ca8e9878SJed Brown       PCTFS_ivec_zero(buf1,buf_size);
606827bd09bSSatish Balay 
607827bd09bSSatish Balay       /* load msg buffer */
608827bd09bSSatish Balay       for (start=end,end+=per_load,i_start=i; (offset=*ptr2)<end; i++, ptr2++)
609827bd09bSSatish Balay         {
610827bd09bSSatish Balay           offset = (offset-start)*p_mask_size;
611*ca8e9878SJed Brown           PCTFS_ivec_copy(buf1+offset,p_mask,p_mask_size);
612827bd09bSSatish Balay         }
613827bd09bSSatish Balay 
614827bd09bSSatish Balay       /* GLOBAL: pass buffer */
615b1c944f5SJed Brown       ierr = PCTFS_giop(buf1,buf2,buf_size,&oper);CHKERRQ(ierr);
616827bd09bSSatish Balay 
617827bd09bSSatish Balay 
618827bd09bSSatish Balay       /* unload buffer into ngh_buf */
619827bd09bSSatish Balay       ptr2=(elms+i_start);
620827bd09bSSatish Balay       for(ptr3=buf1,j=start; j<end; ptr3+=p_mask_size,j++)
621827bd09bSSatish Balay         {
622827bd09bSSatish Balay           /* I own it ... may have to pairwise it */
623827bd09bSSatish Balay           if (j==*ptr2)
624827bd09bSSatish Balay             {
625827bd09bSSatish Balay               /* do i share it w/anyone? */
626*ca8e9878SJed Brown               ct1 = PCTFS_ct_bits((char *)ptr3,p_mask_size*sizeof(PetscInt));
627827bd09bSSatish Balay               /* guess not */
628827bd09bSSatish Balay               if (ct1<2)
629827bd09bSSatish Balay                 {ptr2++; ptr1+=p_mask_size; continue;}
630827bd09bSSatish Balay 
631827bd09bSSatish Balay               /* i do ... so keep info and turn off my bit */
632*ca8e9878SJed Brown               PCTFS_ivec_copy(ptr1,ptr3,p_mask_size);
633*ca8e9878SJed Brown               ierr = PCTFS_ivec_xor(ptr1,p_mask,p_mask_size);CHKERRQ(ierr);
634*ca8e9878SJed Brown               ierr = PCTFS_ivec_or(sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr);
635827bd09bSSatish Balay 
636827bd09bSSatish Balay               /* is it to be done pairwise? */
637827bd09bSSatish Balay               if (--ct1<=level)
638827bd09bSSatish Balay                 {
639827bd09bSSatish Balay                   npw++;
640827bd09bSSatish Balay 
641827bd09bSSatish Balay                   /* turn on high bit to indicate pw need to process */
642827bd09bSSatish Balay                   *ptr2++ |= TOP_BIT;
643*ca8e9878SJed Brown                   ierr = PCTFS_ivec_or(pw_sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr);
644827bd09bSSatish Balay                   ptr1+=p_mask_size;
645827bd09bSSatish Balay                   continue;
646827bd09bSSatish Balay                 }
647827bd09bSSatish Balay 
648827bd09bSSatish Balay               /* get set for next and note that I have a tree contribution */
649827bd09bSSatish Balay               /* could save exact elm index for tree here -> save a search */
650827bd09bSSatish Balay               ptr2++; ptr1+=p_mask_size; ntree_map++;
651827bd09bSSatish Balay             }
652827bd09bSSatish Balay           /* i don't but still might be involved in tree */
653827bd09bSSatish Balay           else
654827bd09bSSatish Balay             {
655827bd09bSSatish Balay 
656827bd09bSSatish Balay               /* shared by how many? */
657*ca8e9878SJed Brown               ct1 = PCTFS_ct_bits((char *)ptr3,p_mask_size*sizeof(PetscInt));
658827bd09bSSatish Balay 
659827bd09bSSatish Balay               /* none! */
660f1ed62a8SBarry Smith               if (ct1<2) continue;
661827bd09bSSatish Balay 
662827bd09bSSatish Balay               /* is it going to be done pairwise? but not by me of course!*/
663f1ed62a8SBarry Smith               if (--ct1<=level) continue;
664827bd09bSSatish Balay             }
665827bd09bSSatish Balay           /* LATER we're going to have to process it NOW */
666827bd09bSSatish Balay           /* nope ... tree it */
667f1ed62a8SBarry Smith           ierr = place_in_tree(j);CHKERRQ(ierr);
668827bd09bSSatish Balay         }
669827bd09bSSatish Balay     }
670827bd09bSSatish Balay 
671a501084fSBarry Smith   free((void*)t_mask);
672a501084fSBarry Smith   free((void*)buf1);
673827bd09bSSatish Balay 
674827bd09bSSatish Balay   gs->len_pw_list=npw;
675*ca8e9878SJed Brown   gs->num_nghs = PCTFS_ct_bits((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt));
676827bd09bSSatish Balay 
677827bd09bSSatish Balay   /* expand from bit mask list to int list and save ngh list */
67852f87cdaSBarry Smith   gs->nghs = (PetscInt*) malloc(gs->num_nghs * sizeof(PetscInt));
679*ca8e9878SJed Brown   PCTFS_bm_to_proc((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt),gs->nghs);
680827bd09bSSatish Balay 
681*ca8e9878SJed Brown   gs->num_pw_nghs = PCTFS_ct_bits((char *)pw_sh_proc_mask,p_mask_size*sizeof(PetscInt));
682827bd09bSSatish Balay 
683827bd09bSSatish Balay   oper = GL_MAX;
684827bd09bSSatish Balay   ct1 = gs->num_nghs;
685b1c944f5SJed Brown   ierr = PCTFS_giop(&ct1,&ct2,1,&oper);CHKERRQ(ierr);
686827bd09bSSatish Balay   gs->max_nghs = ct1;
687827bd09bSSatish Balay 
688827bd09bSSatish Balay   gs->tree_map_sz  = ntree_map;
689827bd09bSSatish Balay   gs->max_left_over=ntree;
690827bd09bSSatish Balay 
691a501084fSBarry Smith   free((void*)p_mask);
692a501084fSBarry Smith   free((void*)sh_proc_mask);
6933fdc5746SBarry Smith   PetscFunctionReturn(0);
694827bd09bSSatish Balay }
695827bd09bSSatish Balay 
696f1ed62a8SBarry Smith /******************************************************************************/
697*ca8e9878SJed Brown static PetscErrorCode set_pairwise(PCTFS_gs_id *gs)
698827bd09bSSatish Balay {
69952f87cdaSBarry Smith    PetscInt i, j;
70052f87cdaSBarry Smith   PetscInt p_mask_size;
70152f87cdaSBarry Smith   PetscInt *p_mask, *sh_proc_mask, *tmp_proc_mask;
70252f87cdaSBarry Smith   PetscInt *ngh_buf, *buf2;
70352f87cdaSBarry Smith   PetscInt offset;
70452f87cdaSBarry Smith   PetscInt *msg_list, *msg_size, **msg_nodes, nprs;
70552f87cdaSBarry Smith   PetscInt *pairwise_elm_list, len_pair_list=0;
70652f87cdaSBarry Smith   PetscInt *iptr, t1, i_start, nel, *elms;
70752f87cdaSBarry Smith   PetscInt ct;
708f1ed62a8SBarry Smith   PetscErrorCode ierr;
709827bd09bSSatish Balay 
7103fdc5746SBarry Smith   PetscFunctionBegin;
711827bd09bSSatish Balay   /* to make life easier */
712827bd09bSSatish Balay   nel  = gs->nel;
713827bd09bSSatish Balay   elms = gs->elms;
714827bd09bSSatish Balay   ngh_buf = gs->ngh_buf;
715827bd09bSSatish Balay   sh_proc_mask  = gs->pw_nghs;
716827bd09bSSatish Balay 
717827bd09bSSatish Balay   /* need a few temp masks */
718*ca8e9878SJed Brown   p_mask_size   = PCTFS_len_bit_mask(PCTFS_num_nodes);
71952f87cdaSBarry Smith   p_mask        = (PetscInt*) malloc(p_mask_size);
72052f87cdaSBarry Smith   tmp_proc_mask = (PetscInt*) malloc(p_mask_size);
721827bd09bSSatish Balay 
722b1c944f5SJed Brown   /* set mask to my PCTFS_my_id's bit mask */
723*ca8e9878SJed Brown   ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr);
724827bd09bSSatish Balay 
725a501084fSBarry Smith   p_mask_size /= sizeof(PetscInt);
726827bd09bSSatish Balay 
727827bd09bSSatish Balay   len_pair_list=gs->len_pw_list;
72852f87cdaSBarry Smith   gs->pw_elm_list=pairwise_elm_list=(PetscInt*)malloc((len_pair_list+1)*sizeof(PetscInt));
729827bd09bSSatish Balay 
730827bd09bSSatish Balay   /* how many processors (nghs) do we have to exchange with? */
731*ca8e9878SJed Brown   nprs=gs->num_pairs=PCTFS_ct_bits((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt));
732827bd09bSSatish Balay 
733827bd09bSSatish Balay 
734*ca8e9878SJed Brown   /* allocate space for PCTFS_gs_gop() info */
73552f87cdaSBarry Smith   gs->pair_list = msg_list = (PetscInt *)  malloc(sizeof(PetscInt)*nprs);
73652f87cdaSBarry Smith   gs->msg_sizes = msg_size  = (PetscInt *)  malloc(sizeof(PetscInt)*nprs);
73752f87cdaSBarry Smith   gs->node_list = msg_nodes = (PetscInt **) malloc(sizeof(PetscInt*)*(nprs+1));
738827bd09bSSatish Balay 
739827bd09bSSatish Balay   /* init msg_size list */
740*ca8e9878SJed Brown   ierr = PCTFS_ivec_zero(msg_size,nprs);CHKERRQ(ierr);
741827bd09bSSatish Balay 
742827bd09bSSatish Balay   /* expand from bit mask list to int list */
743*ca8e9878SJed Brown   ierr = PCTFS_bm_to_proc((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt),msg_list);CHKERRQ(ierr);
744827bd09bSSatish Balay 
745827bd09bSSatish Balay   /* keep list of elements being handled pairwise */
746827bd09bSSatish Balay   for (i=j=0;i<nel;i++)
747827bd09bSSatish Balay     {
748827bd09bSSatish Balay       if (elms[i] & TOP_BIT)
749827bd09bSSatish Balay         {elms[i] ^= TOP_BIT; pairwise_elm_list[j++] = i;}
750827bd09bSSatish Balay     }
751827bd09bSSatish Balay   pairwise_elm_list[j] = -1;
752827bd09bSSatish Balay 
753a501084fSBarry Smith   gs->msg_ids_out = (MPI_Request *)  malloc(sizeof(MPI_Request)*(nprs+1));
754827bd09bSSatish Balay   gs->msg_ids_out[nprs] = MPI_REQUEST_NULL;
755a501084fSBarry Smith   gs->msg_ids_in = (MPI_Request *)  malloc(sizeof(MPI_Request)*(nprs+1));
756827bd09bSSatish Balay   gs->msg_ids_in[nprs] = MPI_REQUEST_NULL;
757a501084fSBarry Smith   gs->pw_vals = (PetscScalar *) malloc(sizeof(PetscScalar)*len_pair_list*vec_sz);
758827bd09bSSatish Balay 
759827bd09bSSatish Balay   /* find who goes to each processor */
760827bd09bSSatish Balay   for (i_start=i=0;i<nprs;i++)
761827bd09bSSatish Balay     {
762827bd09bSSatish Balay       /* processor i's mask */
763*ca8e9878SJed Brown       ierr = PCTFS_set_bit_mask(p_mask,p_mask_size*sizeof(PetscInt),msg_list[i]);CHKERRQ(ierr);
764827bd09bSSatish Balay 
765827bd09bSSatish Balay       /* det # going to processor i */
766827bd09bSSatish Balay       for (ct=j=0;j<len_pair_list;j++)
767827bd09bSSatish Balay         {
768827bd09bSSatish Balay           buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size);
769*ca8e9878SJed Brown           ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr);
770*ca8e9878SJed Brown           if (PCTFS_ct_bits((char *)tmp_proc_mask,p_mask_size*sizeof(PetscInt)))
771827bd09bSSatish Balay             {ct++;}
772827bd09bSSatish Balay         }
773827bd09bSSatish Balay       msg_size[i] = ct;
77439945688SSatish Balay       i_start = PetscMax(i_start,ct);
775827bd09bSSatish Balay 
776827bd09bSSatish Balay       /*space to hold nodes in message to first neighbor */
77752f87cdaSBarry Smith       msg_nodes[i] = iptr = (PetscInt*) malloc(sizeof(PetscInt)*(ct+1));
778827bd09bSSatish Balay 
779827bd09bSSatish Balay       for (j=0;j<len_pair_list;j++)
780827bd09bSSatish Balay         {
781827bd09bSSatish Balay           buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size);
782*ca8e9878SJed Brown           ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr);
783*ca8e9878SJed Brown           if (PCTFS_ct_bits((char *)tmp_proc_mask,p_mask_size*sizeof(PetscInt)))
784827bd09bSSatish Balay             {*iptr++ = j;}
785827bd09bSSatish Balay         }
786827bd09bSSatish Balay       *iptr = -1;
787827bd09bSSatish Balay     }
788827bd09bSSatish Balay   msg_nodes[nprs] = NULL;
789827bd09bSSatish Balay 
790827bd09bSSatish Balay   j=gs->loc_node_pairs=i_start;
791827bd09bSSatish Balay   t1 = GL_MAX;
792b1c944f5SJed Brown   ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr);
793827bd09bSSatish Balay   gs->max_node_pairs = i_start;
794827bd09bSSatish Balay 
795827bd09bSSatish Balay   i_start=j;
796827bd09bSSatish Balay   t1 = GL_MIN;
797b1c944f5SJed Brown   ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr);
798827bd09bSSatish Balay   gs->min_node_pairs = i_start;
799827bd09bSSatish Balay 
800827bd09bSSatish Balay   i_start=j;
801827bd09bSSatish Balay   t1 = GL_ADD;
802b1c944f5SJed Brown   ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr);
803b1c944f5SJed Brown   gs->avg_node_pairs = i_start/PCTFS_num_nodes + 1;
804827bd09bSSatish Balay 
805827bd09bSSatish Balay   i_start=nprs;
806827bd09bSSatish Balay   t1 = GL_MAX;
807b1c944f5SJed Brown   PCTFS_giop(&i_start,&offset,1,&t1);
808827bd09bSSatish Balay   gs->max_pairs = i_start;
809827bd09bSSatish Balay 
810827bd09bSSatish Balay 
811827bd09bSSatish Balay   /* remap pairwise in tail of gsi_via_bit_mask() */
812*ca8e9878SJed Brown   gs->msg_total = PCTFS_ivec_sum(gs->msg_sizes,nprs);
813a501084fSBarry Smith   gs->out = (PetscScalar *) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz);
814a501084fSBarry Smith   gs->in  = (PetscScalar *) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz);
815827bd09bSSatish Balay 
816827bd09bSSatish Balay   /* reset malloc pool */
817a501084fSBarry Smith   free((void*)p_mask);
818a501084fSBarry Smith   free((void*)tmp_proc_mask);
8193fdc5746SBarry Smith   PetscFunctionReturn(0);
820827bd09bSSatish Balay }
821827bd09bSSatish Balay 
822f1ed62a8SBarry Smith /* to do pruned tree just save ngh buf copy for each one and decode here!
823827bd09bSSatish Balay ******************************************************************************/
824*ca8e9878SJed Brown static PetscErrorCode set_tree(PCTFS_gs_id *gs)
825827bd09bSSatish Balay {
82652f87cdaSBarry Smith   PetscInt i, j, n, nel;
82752f87cdaSBarry Smith   PetscInt *iptr_in, *iptr_out, *tree_elms, *elms;
828827bd09bSSatish Balay 
8293fdc5746SBarry Smith   PetscFunctionBegin;
830827bd09bSSatish Balay   /* local work ptrs */
831827bd09bSSatish Balay   elms = gs->elms;
832827bd09bSSatish Balay   nel     = gs->nel;
833827bd09bSSatish Balay 
834827bd09bSSatish Balay   /* how many via tree */
835827bd09bSSatish Balay   gs->tree_nel  = n = ntree;
836827bd09bSSatish Balay   gs->tree_elms = tree_elms = iptr_in = tree_buf;
837a501084fSBarry Smith   gs->tree_buf  = (PetscScalar *) malloc(sizeof(PetscScalar)*n*vec_sz);
838a501084fSBarry Smith   gs->tree_work = (PetscScalar *) malloc(sizeof(PetscScalar)*n*vec_sz);
839827bd09bSSatish Balay   j=gs->tree_map_sz;
84052f87cdaSBarry Smith   gs->tree_map_in = iptr_in  = (PetscInt*) malloc(sizeof(PetscInt)*(j+1));
84152f87cdaSBarry Smith   gs->tree_map_out = iptr_out = (PetscInt*) malloc(sizeof(PetscInt)*(j+1));
842827bd09bSSatish Balay 
843827bd09bSSatish Balay   /* search the longer of the two lists */
844827bd09bSSatish Balay   /* note ... could save this info in get_ngh_buf and save searches */
845827bd09bSSatish Balay   if (n<=nel)
846827bd09bSSatish Balay     {
847827bd09bSSatish Balay       /* bijective fct w/remap - search elm list */
848827bd09bSSatish Balay       for (i=0; i<n; i++)
849827bd09bSSatish Balay         {
850*ca8e9878SJed Brown           if ((j=PCTFS_ivec_binary_search(*tree_elms++,elms,nel))>=0)
851827bd09bSSatish Balay             {*iptr_in++ = j; *iptr_out++ = i;}
852827bd09bSSatish Balay         }
853827bd09bSSatish Balay     }
854827bd09bSSatish Balay   else
855827bd09bSSatish Balay     {
856827bd09bSSatish Balay       for (i=0; i<nel; i++)
857827bd09bSSatish Balay         {
858*ca8e9878SJed Brown           if ((j=PCTFS_ivec_binary_search(*elms++,tree_elms,n))>=0)
859827bd09bSSatish Balay             {*iptr_in++ = i; *iptr_out++ = j;}
860827bd09bSSatish Balay         }
861827bd09bSSatish Balay     }
862827bd09bSSatish Balay 
863827bd09bSSatish Balay   /* sentinel */
864827bd09bSSatish Balay   *iptr_in = *iptr_out = -1;
8653fdc5746SBarry Smith   PetscFunctionReturn(0);
866827bd09bSSatish Balay }
867827bd09bSSatish Balay 
868f1ed62a8SBarry Smith /******************************************************************************/
869*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_out( PCTFS_gs_id *gs,  PetscScalar *vals)
870827bd09bSSatish Balay {
87152f87cdaSBarry Smith   PetscInt *num, *map, **reduce;
872a501084fSBarry Smith   PetscScalar tmp;
873827bd09bSSatish Balay 
8743fdc5746SBarry Smith   PetscFunctionBegin;
875827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
876827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
877827bd09bSSatish Balay   while ((map = *reduce++))
878827bd09bSSatish Balay     {
879827bd09bSSatish Balay       /* wall */
880827bd09bSSatish Balay       if (*num == 2)
881827bd09bSSatish Balay         {
882827bd09bSSatish Balay           num ++;
883827bd09bSSatish Balay           vals[map[1]] = vals[map[0]];
884827bd09bSSatish Balay         }
885827bd09bSSatish Balay       /* corner shared by three elements */
886827bd09bSSatish Balay       else if (*num == 3)
887827bd09bSSatish Balay         {
888827bd09bSSatish Balay           num ++;
889827bd09bSSatish Balay           vals[map[2]] = vals[map[1]] = vals[map[0]];
890827bd09bSSatish Balay         }
891827bd09bSSatish Balay       /* corner shared by four elements */
892827bd09bSSatish Balay       else if (*num == 4)
893827bd09bSSatish Balay         {
894827bd09bSSatish Balay           num ++;
895827bd09bSSatish Balay           vals[map[3]] = vals[map[2]] = vals[map[1]] = vals[map[0]];
896827bd09bSSatish Balay         }
897827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
898827bd09bSSatish Balay       else
899827bd09bSSatish Balay         {
900827bd09bSSatish Balay           num++;
901827bd09bSSatish Balay           tmp = *(vals + *map++);
902827bd09bSSatish Balay           while (*map >= 0)
903827bd09bSSatish Balay             {*(vals + *map++) = tmp;}
904827bd09bSSatish Balay         }
905827bd09bSSatish Balay     }
9063fdc5746SBarry Smith   PetscFunctionReturn(0);
907827bd09bSSatish Balay }
908827bd09bSSatish Balay 
9097b1ae94cSBarry Smith /******************************************************************************/
910*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_plus( PCTFS_gs_id *gs,  PetscScalar *vals)
911827bd09bSSatish Balay {
91252f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
913a501084fSBarry Smith    PetscScalar tmp;
914827bd09bSSatish Balay 
9153fdc5746SBarry Smith   PetscFunctionBegin;
916827bd09bSSatish Balay   num    = gs->num_local_reduce;
917827bd09bSSatish Balay   reduce = gs->local_reduce;
918827bd09bSSatish Balay   while ((map = *reduce))
919827bd09bSSatish Balay     {
920827bd09bSSatish Balay       /* wall */
921827bd09bSSatish Balay       if (*num == 2)
922827bd09bSSatish Balay         {
923827bd09bSSatish Balay           num ++; reduce++;
924827bd09bSSatish Balay           vals[map[1]] = vals[map[0]] += vals[map[1]];
925827bd09bSSatish Balay         }
926827bd09bSSatish Balay       /* corner shared by three elements */
927827bd09bSSatish Balay       else if (*num == 3)
928827bd09bSSatish Balay         {
929827bd09bSSatish Balay           num ++; reduce++;
930827bd09bSSatish Balay           vals[map[2]]=vals[map[1]]=vals[map[0]]+=(vals[map[1]]+vals[map[2]]);
931827bd09bSSatish Balay         }
932827bd09bSSatish Balay       /* corner shared by four elements */
933827bd09bSSatish Balay       else if (*num == 4)
934827bd09bSSatish Balay         {
935827bd09bSSatish Balay           num ++; reduce++;
936827bd09bSSatish Balay           vals[map[1]]=vals[map[2]]=vals[map[3]]=vals[map[0]] +=
937827bd09bSSatish Balay                                  (vals[map[1]] + vals[map[2]] + vals[map[3]]);
938827bd09bSSatish Balay         }
939827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
940827bd09bSSatish Balay       else
941827bd09bSSatish Balay         {
942827bd09bSSatish Balay           num ++;
943827bd09bSSatish Balay           tmp = 0.0;
944827bd09bSSatish Balay           while (*map >= 0)
945827bd09bSSatish Balay             {tmp += *(vals + *map++);}
946827bd09bSSatish Balay 
947827bd09bSSatish Balay           map = *reduce++;
948827bd09bSSatish Balay           while (*map >= 0)
949827bd09bSSatish Balay             {*(vals + *map++) = tmp;}
950827bd09bSSatish Balay         }
951827bd09bSSatish Balay     }
9523fdc5746SBarry Smith   PetscFunctionReturn(0);
953827bd09bSSatish Balay }
954827bd09bSSatish Balay 
9557b1ae94cSBarry Smith /******************************************************************************/
956*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_in_plus( PCTFS_gs_id *gs,  PetscScalar *vals)
957827bd09bSSatish Balay {
95852f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
959a501084fSBarry Smith    PetscScalar *base;
960827bd09bSSatish Balay 
9613fdc5746SBarry Smith   PetscFunctionBegin;
962827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
963827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
964827bd09bSSatish Balay   while ((map = *reduce++))
965827bd09bSSatish Balay     {
966827bd09bSSatish Balay       /* wall */
967827bd09bSSatish Balay       if (*num == 2)
968827bd09bSSatish Balay         {
969827bd09bSSatish Balay           num ++;
970827bd09bSSatish Balay           vals[map[0]] += vals[map[1]];
971827bd09bSSatish Balay         }
972827bd09bSSatish Balay       /* corner shared by three elements */
973827bd09bSSatish Balay       else if (*num == 3)
974827bd09bSSatish Balay         {
975827bd09bSSatish Balay           num ++;
976827bd09bSSatish Balay           vals[map[0]] += (vals[map[1]] + vals[map[2]]);
977827bd09bSSatish Balay         }
978827bd09bSSatish Balay       /* corner shared by four elements */
979827bd09bSSatish Balay       else if (*num == 4)
980827bd09bSSatish Balay         {
981827bd09bSSatish Balay           num ++;
982827bd09bSSatish Balay           vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]);
983827bd09bSSatish Balay         }
984827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
985827bd09bSSatish Balay       else
986827bd09bSSatish Balay         {
987827bd09bSSatish Balay           num++;
988827bd09bSSatish Balay           base = vals + *map++;
989827bd09bSSatish Balay           while (*map >= 0)
990827bd09bSSatish Balay             {*base += *(vals + *map++);}
991827bd09bSSatish Balay         }
992827bd09bSSatish Balay     }
9933fdc5746SBarry Smith   PetscFunctionReturn(0);
994827bd09bSSatish Balay }
995827bd09bSSatish Balay 
9967b1ae94cSBarry Smith /******************************************************************************/
997*ca8e9878SJed Brown PetscErrorCode PCTFS_gs_free( PCTFS_gs_id *gs)
998827bd09bSSatish Balay {
99952f87cdaSBarry Smith    PetscInt i;
1000827bd09bSSatish Balay 
10013fdc5746SBarry Smith   PetscFunctionBegin;
1002a501084fSBarry Smith   if (gs->nghs) {free((void*) gs->nghs);}
1003a501084fSBarry Smith   if (gs->pw_nghs) {free((void*) gs->pw_nghs);}
1004827bd09bSSatish Balay 
1005827bd09bSSatish Balay   /* tree */
1006827bd09bSSatish Balay   if (gs->max_left_over)
1007827bd09bSSatish Balay     {
1008a501084fSBarry Smith       if (gs->tree_elms) {free((void*) gs->tree_elms);}
1009a501084fSBarry Smith       if (gs->tree_buf) {free((void*) gs->tree_buf);}
1010a501084fSBarry Smith       if (gs->tree_work) {free((void*) gs->tree_work);}
1011a501084fSBarry Smith       if (gs->tree_map_in) {free((void*) gs->tree_map_in);}
1012a501084fSBarry Smith       if (gs->tree_map_out) {free((void*) gs->tree_map_out);}
1013827bd09bSSatish Balay     }
1014827bd09bSSatish Balay 
1015827bd09bSSatish Balay   /* pairwise info */
1016827bd09bSSatish Balay   if (gs->num_pairs)
1017827bd09bSSatish Balay     {
1018827bd09bSSatish Balay       /* should be NULL already */
1019a501084fSBarry Smith       if (gs->ngh_buf) {free((void*) gs->ngh_buf);}
1020a501084fSBarry Smith       if (gs->elms) {free((void*) gs->elms);}
1021a501084fSBarry Smith       if (gs->local_elms) {free((void*) gs->local_elms);}
1022a501084fSBarry Smith       if (gs->companion) {free((void*) gs->companion);}
1023827bd09bSSatish Balay 
1024827bd09bSSatish Balay       /* only set if pairwise */
1025a501084fSBarry Smith       if (gs->vals) {free((void*) gs->vals);}
1026a501084fSBarry Smith       if (gs->in) {free((void*) gs->in);}
1027a501084fSBarry Smith       if (gs->out) {free((void*) gs->out);}
1028a501084fSBarry Smith       if (gs->msg_ids_in) {free((void*) gs->msg_ids_in);}
1029a501084fSBarry Smith       if (gs->msg_ids_out) {free((void*) gs->msg_ids_out);}
1030a501084fSBarry Smith       if (gs->pw_vals) {free((void*) gs->pw_vals);}
1031a501084fSBarry Smith       if (gs->pw_elm_list) {free((void*) gs->pw_elm_list);}
1032827bd09bSSatish Balay       if (gs->node_list)
1033827bd09bSSatish Balay         {
1034827bd09bSSatish Balay           for (i=0;i<gs->num_pairs;i++)
1035a501084fSBarry Smith             {if (gs->node_list[i]) {free((void*) gs->node_list[i]);}}
1036a501084fSBarry Smith           free((void*) gs->node_list);
1037827bd09bSSatish Balay         }
1038a501084fSBarry Smith       if (gs->msg_sizes) {free((void*) gs->msg_sizes);}
1039a501084fSBarry Smith       if (gs->pair_list) {free((void*) gs->pair_list);}
1040827bd09bSSatish Balay     }
1041827bd09bSSatish Balay 
1042827bd09bSSatish Balay   /* local info */
1043827bd09bSSatish Balay   if (gs->num_local_total>=0)
1044827bd09bSSatish Balay     {
1045827bd09bSSatish Balay       for (i=0;i<gs->num_local_total+1;i++)
1046827bd09bSSatish Balay         /*      for (i=0;i<gs->num_local_total;i++) */
1047827bd09bSSatish Balay         {
1048827bd09bSSatish Balay           if (gs->num_gop_local_reduce[i])
1049a501084fSBarry Smith             {free((void*) gs->gop_local_reduce[i]);}
1050827bd09bSSatish Balay         }
1051827bd09bSSatish Balay     }
1052827bd09bSSatish Balay 
1053827bd09bSSatish Balay   /* if intersection tree/pairwise and local isn't empty */
1054a501084fSBarry Smith   if (gs->gop_local_reduce) {free((void*) gs->gop_local_reduce);}
1055a501084fSBarry Smith   if (gs->num_gop_local_reduce) {free((void*) gs->num_gop_local_reduce);}
1056827bd09bSSatish Balay 
1057a501084fSBarry Smith   free((void*) gs);
10583fdc5746SBarry Smith   PetscFunctionReturn(0);
1059827bd09bSSatish Balay }
1060827bd09bSSatish Balay 
10617b1ae94cSBarry Smith /******************************************************************************/
1062*ca8e9878SJed Brown PetscErrorCode PCTFS_gs_gop_vec( PCTFS_gs_id *gs,  PetscScalar *vals,  const char *op,  PetscInt step)
1063827bd09bSSatish Balay {
1064d1528f56SBarry Smith   PetscErrorCode ierr;
1065d1528f56SBarry Smith 
10663fdc5746SBarry Smith   PetscFunctionBegin;
1067827bd09bSSatish Balay   switch (*op) {
1068827bd09bSSatish Balay   case '+':
1069*ca8e9878SJed Brown     PCTFS_gs_gop_vec_plus(gs,vals,step);
1070827bd09bSSatish Balay     break;
1071827bd09bSSatish Balay   default:
1072*ca8e9878SJed Brown     ierr = PetscInfo1(0,"PCTFS_gs_gop_vec() :: %c is not a valid op",op[0]);CHKERRQ(ierr);
1073*ca8e9878SJed Brown     ierr = PetscInfo(0,"PCTFS_gs_gop_vec() :: default :: plus");CHKERRQ(ierr);
1074*ca8e9878SJed Brown     PCTFS_gs_gop_vec_plus(gs,vals,step);
1075827bd09bSSatish Balay     break;
1076827bd09bSSatish Balay   }
10773fdc5746SBarry Smith   PetscFunctionReturn(0);
1078827bd09bSSatish Balay }
1079827bd09bSSatish Balay 
10807b1ae94cSBarry Smith /******************************************************************************/
1081*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_plus( PCTFS_gs_id *gs,  PetscScalar *vals,  PetscInt step)
1082827bd09bSSatish Balay {
10833fdc5746SBarry Smith   PetscFunctionBegin;
1084*ca8e9878SJed Brown   if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_gs_gop_vec() passed NULL gs handle!!!");
1085827bd09bSSatish Balay 
1086827bd09bSSatish Balay   /* local only operations!!! */
1087827bd09bSSatish Balay   if (gs->num_local)
1088*ca8e9878SJed Brown     {PCTFS_gs_gop_vec_local_plus(gs,vals,step);}
1089827bd09bSSatish Balay 
1090827bd09bSSatish Balay   /* if intersection tree/pairwise and local isn't empty */
1091827bd09bSSatish Balay   if (gs->num_local_gop)
1092827bd09bSSatish Balay     {
1093*ca8e9878SJed Brown       PCTFS_gs_gop_vec_local_in_plus(gs,vals,step);
1094827bd09bSSatish Balay 
1095827bd09bSSatish Balay       /* pairwise */
1096827bd09bSSatish Balay       if (gs->num_pairs)
1097*ca8e9878SJed Brown         {PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step);}
1098827bd09bSSatish Balay 
1099827bd09bSSatish Balay       /* tree */
1100827bd09bSSatish Balay       else if (gs->max_left_over)
1101*ca8e9878SJed Brown         {PCTFS_gs_gop_vec_tree_plus(gs,vals,step);}
1102827bd09bSSatish Balay 
1103*ca8e9878SJed Brown       PCTFS_gs_gop_vec_local_out(gs,vals,step);
1104827bd09bSSatish Balay     }
1105827bd09bSSatish Balay   /* if intersection tree/pairwise and local is empty */
1106827bd09bSSatish Balay   else
1107827bd09bSSatish Balay     {
1108827bd09bSSatish Balay       /* pairwise */
1109827bd09bSSatish Balay       if (gs->num_pairs)
1110*ca8e9878SJed Brown         {PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step);}
1111827bd09bSSatish Balay 
1112827bd09bSSatish Balay       /* tree */
1113827bd09bSSatish Balay       else if (gs->max_left_over)
1114*ca8e9878SJed Brown         {PCTFS_gs_gop_vec_tree_plus(gs,vals,step);}
1115827bd09bSSatish Balay     }
11163fdc5746SBarry Smith   PetscFunctionReturn(0);
1117827bd09bSSatish Balay }
1118827bd09bSSatish Balay 
11197b1ae94cSBarry Smith /******************************************************************************/
1120*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_plus( PCTFS_gs_id *gs,  PetscScalar *vals, PetscInt step)
1121827bd09bSSatish Balay {
112252f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
1123a501084fSBarry Smith    PetscScalar *base;
1124827bd09bSSatish Balay 
11253fdc5746SBarry Smith   PetscFunctionBegin;
1126827bd09bSSatish Balay   num    = gs->num_local_reduce;
1127827bd09bSSatish Balay   reduce = gs->local_reduce;
1128827bd09bSSatish Balay   while ((map = *reduce))
1129827bd09bSSatish Balay     {
1130827bd09bSSatish Balay       base = vals + map[0] * step;
1131827bd09bSSatish Balay 
1132827bd09bSSatish Balay       /* wall */
1133827bd09bSSatish Balay       if (*num == 2)
1134827bd09bSSatish Balay         {
1135827bd09bSSatish Balay           num++; reduce++;
1136*ca8e9878SJed Brown           PCTFS_rvec_add (base,vals+map[1]*step,step);
1137*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[1]*step,base,step);
1138827bd09bSSatish Balay         }
1139827bd09bSSatish Balay       /* corner shared by three elements */
1140827bd09bSSatish Balay       else if (*num == 3)
1141827bd09bSSatish Balay         {
1142827bd09bSSatish Balay           num++; reduce++;
1143*ca8e9878SJed Brown           PCTFS_rvec_add (base,vals+map[1]*step,step);
1144*ca8e9878SJed Brown           PCTFS_rvec_add (base,vals+map[2]*step,step);
1145*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[2]*step,base,step);
1146*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[1]*step,base,step);
1147827bd09bSSatish Balay         }
1148827bd09bSSatish Balay       /* corner shared by four elements */
1149827bd09bSSatish Balay       else if (*num == 4)
1150827bd09bSSatish Balay         {
1151827bd09bSSatish Balay           num++; reduce++;
1152*ca8e9878SJed Brown           PCTFS_rvec_add (base,vals+map[1]*step,step);
1153*ca8e9878SJed Brown           PCTFS_rvec_add (base,vals+map[2]*step,step);
1154*ca8e9878SJed Brown           PCTFS_rvec_add (base,vals+map[3]*step,step);
1155*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[3]*step,base,step);
1156*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[2]*step,base,step);
1157*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[1]*step,base,step);
1158827bd09bSSatish Balay         }
1159827bd09bSSatish Balay       /* general case ... odd geoms ... 3D */
1160827bd09bSSatish Balay       else
1161827bd09bSSatish Balay         {
1162827bd09bSSatish Balay           num++;
1163827bd09bSSatish Balay           while (*++map >= 0)
1164*ca8e9878SJed Brown             {PCTFS_rvec_add (base,vals+*map*step,step);}
1165827bd09bSSatish Balay 
1166827bd09bSSatish Balay           map = *reduce;
1167827bd09bSSatish Balay           while (*++map >= 0)
1168*ca8e9878SJed Brown             {PCTFS_rvec_copy(vals+*map*step,base,step);}
1169827bd09bSSatish Balay 
1170827bd09bSSatish Balay           reduce++;
1171827bd09bSSatish Balay         }
1172827bd09bSSatish Balay     }
11733fdc5746SBarry Smith   PetscFunctionReturn(0);
1174827bd09bSSatish Balay }
1175827bd09bSSatish Balay 
11767b1ae94cSBarry Smith /******************************************************************************/
1177*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus( PCTFS_gs_id *gs,  PetscScalar *vals, PetscInt step)
1178827bd09bSSatish Balay {
117952f87cdaSBarry Smith    PetscInt  *num, *map, **reduce;
1180a501084fSBarry Smith    PetscScalar *base;
11813fdc5746SBarry Smith   PetscFunctionBegin;
1182827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
1183827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
1184827bd09bSSatish Balay   while ((map = *reduce++))
1185827bd09bSSatish Balay     {
1186827bd09bSSatish Balay       base = vals + map[0] * step;
1187827bd09bSSatish Balay 
1188827bd09bSSatish Balay       /* wall */
1189827bd09bSSatish Balay       if (*num == 2)
1190827bd09bSSatish Balay         {
1191827bd09bSSatish Balay           num ++;
1192*ca8e9878SJed Brown           PCTFS_rvec_add(base,vals+map[1]*step,step);
1193827bd09bSSatish Balay         }
1194827bd09bSSatish Balay       /* corner shared by three elements */
1195827bd09bSSatish Balay       else if (*num == 3)
1196827bd09bSSatish Balay         {
1197827bd09bSSatish Balay           num ++;
1198*ca8e9878SJed Brown           PCTFS_rvec_add(base,vals+map[1]*step,step);
1199*ca8e9878SJed Brown           PCTFS_rvec_add(base,vals+map[2]*step,step);
1200827bd09bSSatish Balay         }
1201827bd09bSSatish Balay       /* corner shared by four elements */
1202827bd09bSSatish Balay       else if (*num == 4)
1203827bd09bSSatish Balay         {
1204827bd09bSSatish Balay           num ++;
1205*ca8e9878SJed Brown           PCTFS_rvec_add(base,vals+map[1]*step,step);
1206*ca8e9878SJed Brown           PCTFS_rvec_add(base,vals+map[2]*step,step);
1207*ca8e9878SJed Brown           PCTFS_rvec_add(base,vals+map[3]*step,step);
1208827bd09bSSatish Balay         }
1209827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
1210827bd09bSSatish Balay       else
1211827bd09bSSatish Balay         {
1212827bd09bSSatish Balay           num++;
1213827bd09bSSatish Balay           while (*++map >= 0)
1214*ca8e9878SJed Brown             {PCTFS_rvec_add(base,vals+*map*step,step);}
1215827bd09bSSatish Balay         }
1216827bd09bSSatish Balay     }
12173fdc5746SBarry Smith   PetscFunctionReturn(0);
1218827bd09bSSatish Balay }
1219827bd09bSSatish Balay 
12207b1ae94cSBarry Smith /******************************************************************************/
1221*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_out( PCTFS_gs_id *gs,  PetscScalar *vals, PetscInt step)
1222827bd09bSSatish Balay {
122352f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
1224a501084fSBarry Smith    PetscScalar *base;
1225827bd09bSSatish Balay 
12263fdc5746SBarry Smith   PetscFunctionBegin;
1227827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
1228827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
1229827bd09bSSatish Balay   while ((map = *reduce++))
1230827bd09bSSatish Balay     {
1231827bd09bSSatish Balay       base = vals + map[0] * step;
1232827bd09bSSatish Balay 
1233827bd09bSSatish Balay       /* wall */
1234827bd09bSSatish Balay       if (*num == 2)
1235827bd09bSSatish Balay         {
1236827bd09bSSatish Balay           num ++;
1237*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[1]*step,base,step);
1238827bd09bSSatish Balay         }
1239827bd09bSSatish Balay       /* corner shared by three elements */
1240827bd09bSSatish Balay       else if (*num == 3)
1241827bd09bSSatish Balay         {
1242827bd09bSSatish Balay           num ++;
1243*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[1]*step,base,step);
1244*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[2]*step,base,step);
1245827bd09bSSatish Balay         }
1246827bd09bSSatish Balay       /* corner shared by four elements */
1247827bd09bSSatish Balay       else if (*num == 4)
1248827bd09bSSatish Balay         {
1249827bd09bSSatish Balay           num ++;
1250*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[1]*step,base,step);
1251*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[2]*step,base,step);
1252*ca8e9878SJed Brown           PCTFS_rvec_copy(vals+map[3]*step,base,step);
1253827bd09bSSatish Balay         }
1254827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
1255827bd09bSSatish Balay       else
1256827bd09bSSatish Balay         {
1257827bd09bSSatish Balay           num++;
1258827bd09bSSatish Balay           while (*++map >= 0)
1259*ca8e9878SJed Brown             {PCTFS_rvec_copy(vals+*map*step,base,step);}
1260827bd09bSSatish Balay         }
1261827bd09bSSatish Balay     }
12623fdc5746SBarry Smith   PetscFunctionReturn(0);
1263827bd09bSSatish Balay }
1264827bd09bSSatish Balay 
12657b1ae94cSBarry Smith /******************************************************************************/
1266*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus( PCTFS_gs_id *gs,  PetscScalar *in_vals, PetscInt step)
1267827bd09bSSatish Balay {
1268a501084fSBarry Smith   PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2;
126952f87cdaSBarry Smith   PetscInt *iptr, *msg_list, *msg_size, **msg_nodes;
127052f87cdaSBarry Smith   PetscInt *pw, *list, *size, **nodes;
1271827bd09bSSatish Balay   MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out;
1272827bd09bSSatish Balay   MPI_Status status;
12730805154bSBarry Smith   PetscBLASInt i1 = 1,dstep;
12743fdc5746SBarry Smith   PetscErrorCode ierr;
1275827bd09bSSatish Balay 
12763fdc5746SBarry Smith   PetscFunctionBegin;
1277a501084fSBarry Smith   /* strip and load s */
1278827bd09bSSatish Balay   msg_list =list         = gs->pair_list;
1279827bd09bSSatish Balay   msg_size =size         = gs->msg_sizes;
1280827bd09bSSatish Balay   msg_nodes=nodes        = gs->node_list;
1281827bd09bSSatish Balay   iptr=pw                = gs->pw_elm_list;
1282827bd09bSSatish Balay   dptr1=dptr3            = gs->pw_vals;
1283827bd09bSSatish Balay   msg_ids_in  = ids_in   = gs->msg_ids_in;
1284827bd09bSSatish Balay   msg_ids_out = ids_out  = gs->msg_ids_out;
1285827bd09bSSatish Balay   dptr2                  = gs->out;
1286827bd09bSSatish Balay   in1=in2                = gs->in;
1287827bd09bSSatish Balay 
1288827bd09bSSatish Balay   /* post the receives */
1289827bd09bSSatish Balay   /*  msg_nodes=nodes; */
1290827bd09bSSatish Balay   do
1291827bd09bSSatish Balay     {
1292827bd09bSSatish Balay       /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the
1293827bd09bSSatish Balay          second one *list and do list++ afterwards */
1294*ca8e9878SJed Brown       ierr = MPI_Irecv(in1, *size *step, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr);
12959182e22cSBarry Smith       list++;msg_ids_in++;
1296827bd09bSSatish Balay       in1 += *size++ *step;
1297827bd09bSSatish Balay     }
1298827bd09bSSatish Balay   while (*++msg_nodes);
1299827bd09bSSatish Balay   msg_nodes=nodes;
1300827bd09bSSatish Balay 
1301827bd09bSSatish Balay   /* load gs values into in out gs buffers */
1302827bd09bSSatish Balay   while (*iptr >= 0)
1303827bd09bSSatish Balay     {
1304*ca8e9878SJed Brown       PCTFS_rvec_copy(dptr3,in_vals + *iptr*step,step);
1305827bd09bSSatish Balay       dptr3+=step;
1306827bd09bSSatish Balay       iptr++;
1307827bd09bSSatish Balay     }
1308827bd09bSSatish Balay 
1309827bd09bSSatish Balay   /* load out buffers and post the sends */
1310827bd09bSSatish Balay   while ((iptr = *msg_nodes++))
1311827bd09bSSatish Balay     {
1312827bd09bSSatish Balay       dptr3 = dptr2;
1313827bd09bSSatish Balay       while (*iptr >= 0)
1314827bd09bSSatish Balay         {
1315*ca8e9878SJed Brown           PCTFS_rvec_copy(dptr2,dptr1 + *iptr*step,step);
1316827bd09bSSatish Balay           dptr2+=step;
1317827bd09bSSatish Balay           iptr++;
1318827bd09bSSatish Balay         }
1319*ca8e9878SJed Brown       ierr = MPI_Isend(dptr3, *msg_size *step, MPIU_SCALAR, *msg_list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr);
13209182e22cSBarry Smith       msg_size++; msg_list++;msg_ids_out++;
1321827bd09bSSatish Balay     }
1322827bd09bSSatish Balay 
1323827bd09bSSatish Balay   /* tree */
1324827bd09bSSatish Balay   if (gs->max_left_over)
1325*ca8e9878SJed Brown     {PCTFS_gs_gop_vec_tree_plus(gs,in_vals,step);}
1326827bd09bSSatish Balay 
1327827bd09bSSatish Balay   /* process the received data */
1328827bd09bSSatish Balay   msg_nodes=nodes;
1329a501084fSBarry Smith   while ((iptr = *nodes++)){
1330a501084fSBarry Smith     PetscScalar d1 = 1.0;
1331827bd09bSSatish Balay       /* Should I check the return value of MPI_Wait() or status? */
1332827bd09bSSatish Balay       /* Can this loop be replaced by a call to MPI_Waitall()? */
13339182e22cSBarry Smith       ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr);
13349182e22cSBarry Smith       ids_in++;
1335a501084fSBarry Smith       while (*iptr >= 0) {
13360805154bSBarry Smith 	dstep = PetscBLASIntCast(step);
13374a0de3f6SBarry Smith         BLASaxpy_(&dstep,&d1,in2,&i1,dptr1 + *iptr*step,&i1);
1338827bd09bSSatish Balay 	in2+=step;
1339827bd09bSSatish Balay 	iptr++;
1340827bd09bSSatish Balay       }
1341827bd09bSSatish Balay   }
1342827bd09bSSatish Balay 
1343827bd09bSSatish Balay   /* replace vals */
1344827bd09bSSatish Balay   while (*pw >= 0)
1345827bd09bSSatish Balay     {
1346*ca8e9878SJed Brown       PCTFS_rvec_copy(in_vals + *pw*step,dptr1,step);
1347827bd09bSSatish Balay       dptr1+=step;
1348827bd09bSSatish Balay       pw++;
1349827bd09bSSatish Balay     }
1350827bd09bSSatish Balay 
1351827bd09bSSatish Balay   /* clear isend message handles */
1352827bd09bSSatish Balay   /* This changed for clarity though it could be the same */
1353827bd09bSSatish Balay   while (*msg_nodes++)
1354827bd09bSSatish Balay     /* Should I check the return value of MPI_Wait() or status? */
1355827bd09bSSatish Balay     /* Can this loop be replaced by a call to MPI_Waitall()? */
13569182e22cSBarry Smith     {ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr);ids_out++;}
13579182e22cSBarry Smith 
1358827bd09bSSatish Balay 
13593fdc5746SBarry Smith   PetscFunctionReturn(0);
1360827bd09bSSatish Balay }
1361827bd09bSSatish Balay 
13627b1ae94cSBarry Smith /******************************************************************************/
1363*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_tree_plus( PCTFS_gs_id *gs,  PetscScalar *vals,  PetscInt step)
1364827bd09bSSatish Balay {
136552f87cdaSBarry Smith   PetscInt size, *in, *out;
1366a501084fSBarry Smith   PetscScalar *buf, *work;
136752f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
1368a501084fSBarry Smith   PetscBLASInt i1 = 1;
1369827bd09bSSatish Balay 
13703fdc5746SBarry Smith   PetscFunctionBegin;
1371827bd09bSSatish Balay   /* copy over to local variables */
1372827bd09bSSatish Balay   in   = gs->tree_map_in;
1373827bd09bSSatish Balay   out  = gs->tree_map_out;
1374827bd09bSSatish Balay   buf  = gs->tree_buf;
1375827bd09bSSatish Balay   work = gs->tree_work;
1376827bd09bSSatish Balay   size = gs->tree_nel*step;
1377827bd09bSSatish Balay 
1378827bd09bSSatish Balay   /* zero out collection buffer */
1379*ca8e9878SJed Brown   PCTFS_rvec_zero(buf,size);
1380827bd09bSSatish Balay 
1381827bd09bSSatish Balay 
1382827bd09bSSatish Balay   /* copy over my contributions */
1383827bd09bSSatish Balay   while (*in >= 0)
1384827bd09bSSatish Balay     {
13850805154bSBarry Smith       PetscBLASInt dstep = PetscBLASIntCast(step);
13866e4f4d19SBarry Smith       BLAScopy_(&dstep,vals + *in++*step,&i1,buf + *out++*step,&i1);
1387827bd09bSSatish Balay     }
1388827bd09bSSatish Balay 
1389827bd09bSSatish Balay   /* perform fan in/out on full buffer */
1390b1c944f5SJed Brown   /* must change PCTFS_grop to handle the blas */
1391b1c944f5SJed Brown   PCTFS_grop(buf,work,size,op);
1392827bd09bSSatish Balay 
1393827bd09bSSatish Balay   /* reset */
1394827bd09bSSatish Balay   in   = gs->tree_map_in;
1395827bd09bSSatish Balay   out  = gs->tree_map_out;
1396827bd09bSSatish Balay 
1397827bd09bSSatish Balay   /* get the portion of the results I need */
1398827bd09bSSatish Balay   while (*in >= 0)
1399827bd09bSSatish Balay     {
14000805154bSBarry Smith       PetscBLASInt dstep = PetscBLASIntCast(step);
14016e4f4d19SBarry Smith       BLAScopy_(&dstep,buf + *out++*step,&i1,vals + *in++*step,&i1);
1402827bd09bSSatish Balay     }
14033fdc5746SBarry Smith   PetscFunctionReturn(0);
1404827bd09bSSatish Balay }
1405827bd09bSSatish Balay 
14067b1ae94cSBarry Smith /******************************************************************************/
1407*ca8e9878SJed Brown PetscErrorCode PCTFS_gs_gop_hc( PCTFS_gs_id *gs,  PetscScalar *vals,  const char *op,  PetscInt dim)
1408827bd09bSSatish Balay {
1409d1528f56SBarry Smith   PetscErrorCode ierr;
1410d1528f56SBarry Smith 
14113fdc5746SBarry Smith   PetscFunctionBegin;
1412827bd09bSSatish Balay   switch (*op) {
1413827bd09bSSatish Balay   case '+':
1414*ca8e9878SJed Brown     PCTFS_gs_gop_plus_hc(gs,vals,dim);
1415827bd09bSSatish Balay     break;
1416827bd09bSSatish Balay   default:
1417*ca8e9878SJed Brown     ierr = PetscInfo1(0,"PCTFS_gs_gop_hc() :: %c is not a valid op",op[0]);CHKERRQ(ierr);
1418*ca8e9878SJed Brown     ierr = PetscInfo(0,"PCTFS_gs_gop_hc() :: default :: plus\n");CHKERRQ(ierr);
1419*ca8e9878SJed Brown     PCTFS_gs_gop_plus_hc(gs,vals,dim);
1420827bd09bSSatish Balay     break;
1421827bd09bSSatish Balay   }
14223fdc5746SBarry Smith   PetscFunctionReturn(0);
1423827bd09bSSatish Balay }
1424827bd09bSSatish Balay 
14257b1ae94cSBarry Smith /******************************************************************************/
1426*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_plus_hc( PCTFS_gs_id *gs,  PetscScalar *vals, PetscInt dim)
1427827bd09bSSatish Balay {
14283fdc5746SBarry Smith   PetscFunctionBegin;
1429827bd09bSSatish Balay   /* if there's nothing to do return */
1430827bd09bSSatish Balay   if (dim<=0)
14313fdc5746SBarry Smith     {  PetscFunctionReturn(0);}
1432827bd09bSSatish Balay 
1433827bd09bSSatish Balay   /* can't do more dimensions then exist */
1434b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
1435827bd09bSSatish Balay 
1436827bd09bSSatish Balay   /* local only operations!!! */
1437827bd09bSSatish Balay   if (gs->num_local)
1438*ca8e9878SJed Brown     {PCTFS_gs_gop_local_plus(gs,vals);}
1439827bd09bSSatish Balay 
1440827bd09bSSatish Balay   /* if intersection tree/pairwise and local isn't empty */
1441827bd09bSSatish Balay   if (gs->num_local_gop)
1442827bd09bSSatish Balay     {
1443*ca8e9878SJed Brown       PCTFS_gs_gop_local_in_plus(gs,vals);
1444827bd09bSSatish Balay 
1445827bd09bSSatish Balay       /* pairwise will do tree inside ... */
1446827bd09bSSatish Balay       if (gs->num_pairs)
1447*ca8e9878SJed Brown         {PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim);}
1448827bd09bSSatish Balay 
1449827bd09bSSatish Balay       /* tree only */
1450827bd09bSSatish Balay       else if (gs->max_left_over)
1451*ca8e9878SJed Brown         {PCTFS_gs_gop_tree_plus_hc(gs,vals,dim);}
1452827bd09bSSatish Balay 
1453*ca8e9878SJed Brown       PCTFS_gs_gop_local_out(gs,vals);
1454827bd09bSSatish Balay     }
1455827bd09bSSatish Balay   /* if intersection tree/pairwise and local is empty */
1456827bd09bSSatish Balay   else
1457827bd09bSSatish Balay     {
1458827bd09bSSatish Balay       /* pairwise will do tree inside */
1459827bd09bSSatish Balay       if (gs->num_pairs)
1460*ca8e9878SJed Brown         {PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim);}
1461827bd09bSSatish Balay 
1462827bd09bSSatish Balay       /* tree */
1463827bd09bSSatish Balay       else if (gs->max_left_over)
1464*ca8e9878SJed Brown         {PCTFS_gs_gop_tree_plus_hc(gs,vals,dim);}
1465827bd09bSSatish Balay     }
14663fdc5746SBarry Smith   PetscFunctionReturn(0);
1467827bd09bSSatish Balay }
1468827bd09bSSatish Balay 
14697b1ae94cSBarry Smith /******************************************************************************/
1470*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc( PCTFS_gs_id *gs,  PetscScalar *in_vals, PetscInt dim)
1471827bd09bSSatish Balay {
1472a501084fSBarry Smith    PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2;
147352f87cdaSBarry Smith    PetscInt *iptr, *msg_list, *msg_size, **msg_nodes;
147452f87cdaSBarry Smith    PetscInt *pw, *list, *size, **nodes;
1475827bd09bSSatish Balay   MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out;
1476827bd09bSSatish Balay   MPI_Status status;
147752f87cdaSBarry Smith   PetscInt i, mask=1;
14783fdc5746SBarry Smith   PetscErrorCode ierr;
1479827bd09bSSatish Balay 
14803fdc5746SBarry Smith   PetscFunctionBegin;
1481827bd09bSSatish Balay   for (i=1; i<dim; i++)
1482827bd09bSSatish Balay     {mask<<=1; mask++;}
1483827bd09bSSatish Balay 
1484827bd09bSSatish Balay 
1485a501084fSBarry Smith   /* strip and load s */
1486827bd09bSSatish Balay   msg_list =list         = gs->pair_list;
1487827bd09bSSatish Balay   msg_size =size         = gs->msg_sizes;
1488827bd09bSSatish Balay   msg_nodes=nodes        = gs->node_list;
1489827bd09bSSatish Balay   iptr=pw                = gs->pw_elm_list;
1490827bd09bSSatish Balay   dptr1=dptr3            = gs->pw_vals;
1491827bd09bSSatish Balay   msg_ids_in  = ids_in   = gs->msg_ids_in;
1492827bd09bSSatish Balay   msg_ids_out = ids_out  = gs->msg_ids_out;
1493827bd09bSSatish Balay   dptr2                  = gs->out;
1494827bd09bSSatish Balay   in1=in2                = gs->in;
1495827bd09bSSatish Balay 
1496827bd09bSSatish Balay   /* post the receives */
1497827bd09bSSatish Balay   /*  msg_nodes=nodes; */
1498827bd09bSSatish Balay   do
1499827bd09bSSatish Balay     {
1500827bd09bSSatish Balay       /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the
1501827bd09bSSatish Balay          second one *list and do list++ afterwards */
1502b1c944f5SJed Brown       if ((PCTFS_my_id|mask)==(*list|mask))
1503827bd09bSSatish Balay         {
1504*ca8e9878SJed Brown           ierr = MPI_Irecv(in1, *size, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr);
15059182e22cSBarry Smith           list++; msg_ids_in++;in1 += *size++;
1506827bd09bSSatish Balay         }
1507827bd09bSSatish Balay       else
1508827bd09bSSatish Balay         {list++; size++;}
1509827bd09bSSatish Balay     }
1510827bd09bSSatish Balay   while (*++msg_nodes);
1511827bd09bSSatish Balay 
1512827bd09bSSatish Balay   /* load gs values into in out gs buffers */
1513827bd09bSSatish Balay   while (*iptr >= 0)
1514827bd09bSSatish Balay     {*dptr3++ = *(in_vals + *iptr++);}
1515827bd09bSSatish Balay 
1516827bd09bSSatish Balay   /* load out buffers and post the sends */
1517827bd09bSSatish Balay   msg_nodes=nodes;
1518827bd09bSSatish Balay   list = msg_list;
1519827bd09bSSatish Balay   while ((iptr = *msg_nodes++))
1520827bd09bSSatish Balay     {
1521b1c944f5SJed Brown       if ((PCTFS_my_id|mask)==(*list|mask))
1522827bd09bSSatish Balay         {
1523827bd09bSSatish Balay           dptr3 = dptr2;
1524827bd09bSSatish Balay           while (*iptr >= 0)
1525827bd09bSSatish Balay             {*dptr2++ = *(dptr1 + *iptr++);}
1526827bd09bSSatish Balay           /* CHECK PERSISTENT COMMS MODE FOR ALL THIS STUFF */
1527827bd09bSSatish Balay           /* is msg_ids_out++ correct? */
1528*ca8e9878SJed Brown           ierr = MPI_Isend(dptr3, *msg_size, MPIU_SCALAR, *list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr);
15299182e22cSBarry Smith           msg_size++;list++;msg_ids_out++;
1530827bd09bSSatish Balay         }
1531827bd09bSSatish Balay       else
1532827bd09bSSatish Balay         {list++; msg_size++;}
1533827bd09bSSatish Balay     }
1534827bd09bSSatish Balay 
1535827bd09bSSatish Balay   /* do the tree while we're waiting */
1536827bd09bSSatish Balay   if (gs->max_left_over)
1537*ca8e9878SJed Brown     {PCTFS_gs_gop_tree_plus_hc(gs,in_vals,dim);}
1538827bd09bSSatish Balay 
1539827bd09bSSatish Balay   /* process the received data */
1540827bd09bSSatish Balay   msg_nodes=nodes;
1541827bd09bSSatish Balay   list = msg_list;
1542827bd09bSSatish Balay   while ((iptr = *nodes++))
1543827bd09bSSatish Balay     {
1544b1c944f5SJed Brown       if ((PCTFS_my_id|mask)==(*list|mask))
1545827bd09bSSatish Balay         {
1546827bd09bSSatish Balay           /* Should I check the return value of MPI_Wait() or status? */
1547827bd09bSSatish Balay           /* Can this loop be replaced by a call to MPI_Waitall()? */
15489182e22cSBarry Smith           ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr);
15499182e22cSBarry Smith           ids_in++;
1550827bd09bSSatish Balay           while (*iptr >= 0)
1551827bd09bSSatish Balay             {*(dptr1 + *iptr++) += *in2++;}
1552827bd09bSSatish Balay         }
1553827bd09bSSatish Balay       list++;
1554827bd09bSSatish Balay     }
1555827bd09bSSatish Balay 
1556827bd09bSSatish Balay   /* replace vals */
1557827bd09bSSatish Balay   while (*pw >= 0)
1558827bd09bSSatish Balay     {*(in_vals + *pw++) = *dptr1++;}
1559827bd09bSSatish Balay 
1560827bd09bSSatish Balay   /* clear isend message handles */
1561827bd09bSSatish Balay   /* This changed for clarity though it could be the same */
1562827bd09bSSatish Balay   while (*msg_nodes++)
1563827bd09bSSatish Balay     {
1564b1c944f5SJed Brown       if ((PCTFS_my_id|mask)==(*msg_list|mask))
1565827bd09bSSatish Balay         {
1566827bd09bSSatish Balay           /* Should I check the return value of MPI_Wait() or status? */
1567827bd09bSSatish Balay           /* Can this loop be replaced by a call to MPI_Waitall()? */
15689182e22cSBarry Smith           ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr);
15699182e22cSBarry Smith           ids_out++;
1570827bd09bSSatish Balay         }
1571827bd09bSSatish Balay       msg_list++;
1572827bd09bSSatish Balay     }
1573827bd09bSSatish Balay 
15743fdc5746SBarry Smith   PetscFunctionReturn(0);
1575827bd09bSSatish Balay }
1576827bd09bSSatish Balay 
15777b1ae94cSBarry Smith /******************************************************************************/
1578*ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim)
1579827bd09bSSatish Balay {
158052f87cdaSBarry Smith   PetscInt size;
158152f87cdaSBarry Smith   PetscInt *in, *out;
1582a501084fSBarry Smith   PetscScalar *buf, *work;
158352f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
1584827bd09bSSatish Balay 
15853fdc5746SBarry Smith   PetscFunctionBegin;
1586827bd09bSSatish Balay   in   = gs->tree_map_in;
1587827bd09bSSatish Balay   out  = gs->tree_map_out;
1588827bd09bSSatish Balay   buf  = gs->tree_buf;
1589827bd09bSSatish Balay   work = gs->tree_work;
1590827bd09bSSatish Balay   size = gs->tree_nel;
1591827bd09bSSatish Balay 
1592*ca8e9878SJed Brown   PCTFS_rvec_zero(buf,size);
1593827bd09bSSatish Balay 
1594827bd09bSSatish Balay   while (*in >= 0)
1595827bd09bSSatish Balay     {*(buf + *out++) = *(vals + *in++);}
1596827bd09bSSatish Balay 
1597827bd09bSSatish Balay   in   = gs->tree_map_in;
1598827bd09bSSatish Balay   out  = gs->tree_map_out;
1599827bd09bSSatish Balay 
1600b1c944f5SJed Brown   PCTFS_grop_hc(buf,work,size,op,dim);
1601827bd09bSSatish Balay 
1602827bd09bSSatish Balay   while (*in >= 0)
1603827bd09bSSatish Balay     {*(vals + *in++) = *(buf + *out++);}
16043fdc5746SBarry Smith   PetscFunctionReturn(0);
1605827bd09bSSatish Balay }
1606827bd09bSSatish Balay 
1607827bd09bSSatish Balay 
1608827bd09bSSatish Balay 
1609