xref: /petsc/src/ksp/pc/impls/tfs/gs.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2827bd09bSSatish Balay 
3827bd09bSSatish Balay /***********************************gs.c***************************************
4827bd09bSSatish Balay 
5827bd09bSSatish Balay Author: Henry M. Tufo III
6827bd09bSSatish Balay 
7827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
8827bd09bSSatish Balay 
9827bd09bSSatish Balay snail-mail:
10827bd09bSSatish Balay Division of Applied Mathematics
11827bd09bSSatish Balay Brown University
12827bd09bSSatish Balay Providence, RI 02912
13827bd09bSSatish Balay 
14827bd09bSSatish Balay Last Modification:
15827bd09bSSatish Balay 6.21.97
16827bd09bSSatish Balay ************************************gs.c**************************************/
17827bd09bSSatish Balay 
18827bd09bSSatish Balay /***********************************gs.c***************************************
19827bd09bSSatish Balay File Description:
20827bd09bSSatish Balay -----------------
21827bd09bSSatish Balay 
22827bd09bSSatish Balay ************************************gs.c**************************************/
23827bd09bSSatish Balay 
24*7c4f633dSBarry Smith #include "../src/ksp/pc/impls/tfs/tfs.h"
2539945688SSatish Balay 
26827bd09bSSatish Balay /* default length of number of items via tree - doubles if exceeded */
27827bd09bSSatish Balay #define TREE_BUF_SZ 2048;
28827bd09bSSatish Balay #define GS_VEC_SZ   1
29827bd09bSSatish Balay 
30827bd09bSSatish Balay 
31827bd09bSSatish Balay 
32827bd09bSSatish Balay /***********************************gs.c***************************************
33827bd09bSSatish Balay Type: struct gather_scatter_id
34827bd09bSSatish Balay ------------------------------
35827bd09bSSatish Balay 
36827bd09bSSatish Balay ************************************gs.c**************************************/
37827bd09bSSatish Balay typedef struct gather_scatter_id {
3852f87cdaSBarry Smith   PetscInt id;
3952f87cdaSBarry Smith   PetscInt nel_min;
4052f87cdaSBarry Smith   PetscInt nel_max;
4152f87cdaSBarry Smith   PetscInt nel_sum;
4252f87cdaSBarry Smith   PetscInt negl;
4352f87cdaSBarry Smith   PetscInt gl_max;
4452f87cdaSBarry Smith   PetscInt gl_min;
4552f87cdaSBarry Smith   PetscInt repeats;
4652f87cdaSBarry Smith   PetscInt ordered;
4752f87cdaSBarry Smith   PetscInt positive;
48a501084fSBarry Smith   PetscScalar *vals;
49827bd09bSSatish Balay 
50827bd09bSSatish Balay   /* bit mask info */
5152f87cdaSBarry Smith   PetscInt *my_proc_mask;
5252f87cdaSBarry Smith   PetscInt mask_sz;
5352f87cdaSBarry Smith   PetscInt *ngh_buf;
5452f87cdaSBarry Smith   PetscInt ngh_buf_sz;
5552f87cdaSBarry Smith   PetscInt *nghs;
5652f87cdaSBarry Smith   PetscInt num_nghs;
5752f87cdaSBarry Smith   PetscInt max_nghs;
5852f87cdaSBarry Smith   PetscInt *pw_nghs;
5952f87cdaSBarry Smith   PetscInt num_pw_nghs;
6052f87cdaSBarry Smith   PetscInt *tree_nghs;
6152f87cdaSBarry Smith   PetscInt num_tree_nghs;
62827bd09bSSatish Balay 
6352f87cdaSBarry Smith   PetscInt num_loads;
64827bd09bSSatish Balay 
65827bd09bSSatish Balay   /* repeats == true -> local info */
6652f87cdaSBarry Smith   PetscInt nel;         /* number of unique elememts */
6752f87cdaSBarry Smith   PetscInt *elms;       /* of size nel */
6852f87cdaSBarry Smith   PetscInt nel_total;
6952f87cdaSBarry Smith   PetscInt *local_elms; /* of size nel_total */
7052f87cdaSBarry Smith   PetscInt *companion;  /* of size nel_total */
71827bd09bSSatish Balay 
72827bd09bSSatish Balay   /* local info */
7352f87cdaSBarry Smith   PetscInt num_local_total;
7452f87cdaSBarry Smith   PetscInt local_strength;
7552f87cdaSBarry Smith   PetscInt num_local;
7652f87cdaSBarry Smith   PetscInt *num_local_reduce;
7752f87cdaSBarry Smith   PetscInt **local_reduce;
7852f87cdaSBarry Smith   PetscInt num_local_gop;
7952f87cdaSBarry Smith   PetscInt *num_gop_local_reduce;
8052f87cdaSBarry Smith   PetscInt **gop_local_reduce;
81827bd09bSSatish Balay 
82827bd09bSSatish Balay   /* pairwise info */
8352f87cdaSBarry Smith   PetscInt level;
8452f87cdaSBarry Smith   PetscInt num_pairs;
8552f87cdaSBarry Smith   PetscInt max_pairs;
8652f87cdaSBarry Smith   PetscInt loc_node_pairs;
8752f87cdaSBarry Smith   PetscInt max_node_pairs;
8852f87cdaSBarry Smith   PetscInt min_node_pairs;
8952f87cdaSBarry Smith   PetscInt avg_node_pairs;
9052f87cdaSBarry Smith   PetscInt *pair_list;
9152f87cdaSBarry Smith   PetscInt *msg_sizes;
9252f87cdaSBarry Smith   PetscInt **node_list;
9352f87cdaSBarry Smith   PetscInt len_pw_list;
9452f87cdaSBarry Smith   PetscInt *pw_elm_list;
95a501084fSBarry Smith   PetscScalar *pw_vals;
96827bd09bSSatish Balay 
97827bd09bSSatish Balay   MPI_Request *msg_ids_in;
98827bd09bSSatish Balay   MPI_Request *msg_ids_out;
99827bd09bSSatish Balay 
100a501084fSBarry Smith   PetscScalar *out;
101a501084fSBarry Smith   PetscScalar *in;
10252f87cdaSBarry Smith   PetscInt msg_total;
103827bd09bSSatish Balay 
104827bd09bSSatish Balay   /* tree - crystal accumulator info */
10552f87cdaSBarry Smith   PetscInt max_left_over;
10652f87cdaSBarry Smith   PetscInt *pre;
10752f87cdaSBarry Smith   PetscInt *in_num;
10852f87cdaSBarry Smith   PetscInt *out_num;
10952f87cdaSBarry Smith   PetscInt **in_list;
11052f87cdaSBarry Smith   PetscInt **out_list;
111827bd09bSSatish Balay 
112827bd09bSSatish Balay   /* new tree work*/
11352f87cdaSBarry Smith   PetscInt  tree_nel;
11452f87cdaSBarry Smith   PetscInt *tree_elms;
115a501084fSBarry Smith   PetscScalar *tree_buf;
116a501084fSBarry Smith   PetscScalar *tree_work;
117827bd09bSSatish Balay 
11852f87cdaSBarry Smith   PetscInt  tree_map_sz;
11952f87cdaSBarry Smith   PetscInt *tree_map_in;
12052f87cdaSBarry Smith   PetscInt *tree_map_out;
121827bd09bSSatish Balay 
122827bd09bSSatish Balay   /* current memory status */
12352f87cdaSBarry Smith   PetscInt gl_bss_min;
12452f87cdaSBarry Smith   PetscInt gl_perm_min;
125827bd09bSSatish Balay 
126827bd09bSSatish Balay   /* max segment size for gs_gop_vec() */
12752f87cdaSBarry Smith   PetscInt vec_sz;
128827bd09bSSatish Balay 
129827bd09bSSatish Balay   /* hack to make paul happy */
130827bd09bSSatish Balay   MPI_Comm gs_comm;
131827bd09bSSatish Balay 
132827bd09bSSatish Balay } gs_id;
133827bd09bSSatish Balay 
13452f87cdaSBarry Smith static gs_id *gsi_check_args(PetscInt *elms, PetscInt nel, PetscInt level);
1353fdc5746SBarry Smith static PetscErrorCode gsi_via_bit_mask(gs_id *gs);
1363fdc5746SBarry Smith static PetscErrorCode get_ngh_buf(gs_id *gs);
1373fdc5746SBarry Smith static PetscErrorCode set_pairwise(gs_id *gs);
138827bd09bSSatish Balay static gs_id * gsi_new(void);
1393fdc5746SBarry Smith static PetscErrorCode set_tree(gs_id *gs);
140827bd09bSSatish Balay 
141827bd09bSSatish Balay /* same for all but vector flavor */
1423fdc5746SBarry Smith static PetscErrorCode gs_gop_local_out(gs_id *gs, PetscScalar *vals);
143827bd09bSSatish Balay /* vector flavor */
14452f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_local_out(gs_id *gs, PetscScalar *vals, PetscInt step);
145827bd09bSSatish Balay 
14652f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_plus(gs_id *gs, PetscScalar *in_vals, PetscInt step);
14752f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_pairwise_plus(gs_id *gs, PetscScalar *in_vals, PetscInt step);
14852f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_local_plus(gs_id *gs, PetscScalar *vals, PetscInt step);
14952f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_local_in_plus(gs_id *gs, PetscScalar *vals, PetscInt step);
15052f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_tree_plus(gs_id *gs, PetscScalar *vals, PetscInt step);
151827bd09bSSatish Balay 
152827bd09bSSatish Balay 
1533fdc5746SBarry Smith static PetscErrorCode gs_gop_local_plus(gs_id *gs, PetscScalar *vals);
1543fdc5746SBarry Smith static PetscErrorCode gs_gop_local_in_plus(gs_id *gs, PetscScalar *vals);
155827bd09bSSatish Balay 
15652f87cdaSBarry Smith static PetscErrorCode gs_gop_plus_hc(gs_id *gs, PetscScalar *in_vals, PetscInt dim);
15752f87cdaSBarry Smith static PetscErrorCode gs_gop_pairwise_plus_hc(gs_id *gs, PetscScalar *in_vals, PetscInt dim);
15852f87cdaSBarry Smith static PetscErrorCode gs_gop_tree_plus_hc(gs_id *gs, PetscScalar *vals, PetscInt dim);
159827bd09bSSatish Balay 
160827bd09bSSatish Balay /* global vars */
161827bd09bSSatish Balay /* from comm.c module */
162827bd09bSSatish Balay 
16352f87cdaSBarry Smith static PetscInt num_gs_ids = 0;
164827bd09bSSatish Balay 
165827bd09bSSatish Balay /* should make this dynamic ... later */
16652f87cdaSBarry Smith static PetscInt msg_buf=MAX_MSG_BUF;
16752f87cdaSBarry Smith static PetscInt vec_sz=GS_VEC_SZ;
16852f87cdaSBarry Smith static PetscInt *tree_buf=NULL;
16952f87cdaSBarry Smith static PetscInt tree_buf_sz=0;
17052f87cdaSBarry Smith static PetscInt ntree=0;
171827bd09bSSatish Balay 
172f1ed62a8SBarry Smith /***************************************************************************/
17352f87cdaSBarry Smith PetscErrorCode gs_init_vec_sz(PetscInt size)
174827bd09bSSatish Balay {
1753fdc5746SBarry Smith   PetscFunctionBegin;
176827bd09bSSatish Balay   vec_sz = size;
1773fdc5746SBarry Smith   PetscFunctionReturn(0);
178827bd09bSSatish Balay }
179827bd09bSSatish Balay 
180f1ed62a8SBarry Smith /******************************************************************************/
18152f87cdaSBarry Smith PetscErrorCode gs_init_msg_buf_sz(PetscInt buf_size)
182827bd09bSSatish Balay {
1833fdc5746SBarry Smith   PetscFunctionBegin;
184827bd09bSSatish Balay   msg_buf = buf_size;
1853fdc5746SBarry Smith   PetscFunctionReturn(0);
186827bd09bSSatish Balay }
187827bd09bSSatish Balay 
188f1ed62a8SBarry Smith /******************************************************************************/
18952f87cdaSBarry Smith gs_id *gs_init( PetscInt *elms, PetscInt nel, PetscInt level)
190827bd09bSSatish Balay {
191a501084fSBarry Smith    gs_id *gs;
192827bd09bSSatish Balay   MPI_Group gs_group;
193827bd09bSSatish Balay   MPI_Comm  gs_comm;
194f1ed62a8SBarry Smith   PetscErrorCode ierr;
195827bd09bSSatish Balay 
1963fdc5746SBarry Smith   PetscFunctionBegin;
197827bd09bSSatish Balay   /* ensure that communication package has been initialized */
198827bd09bSSatish Balay   comm_init();
199827bd09bSSatish Balay 
200827bd09bSSatish Balay 
201827bd09bSSatish Balay   /* determines if we have enough dynamic/semi-static memory */
202827bd09bSSatish Balay   /* checks input, allocs and sets gd_id template            */
203827bd09bSSatish Balay   gs = gsi_check_args(elms,nel,level);
204827bd09bSSatish Balay 
205827bd09bSSatish Balay   /* only bit mask version up and working for the moment    */
206827bd09bSSatish Balay   /* LATER :: get int list version working for sparse pblms */
207f1ed62a8SBarry Smith   ierr = gsi_via_bit_mask(gs);CHKERRABORT(PETSC_COMM_WORLD,ierr);
208827bd09bSSatish Balay 
209827bd09bSSatish Balay 
210f1ed62a8SBarry Smith   ierr = MPI_Comm_group(MPI_COMM_WORLD,&gs_group);CHKERRABORT(PETSC_COMM_WORLD,ierr);
211f1ed62a8SBarry Smith   ierr = MPI_Comm_create(MPI_COMM_WORLD,gs_group,&gs_comm);CHKERRABORT(PETSC_COMM_WORLD,ierr);
212827bd09bSSatish Balay   gs->gs_comm=gs_comm;
213827bd09bSSatish Balay 
214827bd09bSSatish Balay   return(gs);
215827bd09bSSatish Balay }
216827bd09bSSatish Balay 
217f1ed62a8SBarry Smith /******************************************************************************/
2180924e98cSBarry Smith static gs_id *gsi_new(void)
219827bd09bSSatish Balay {
220f1ed62a8SBarry Smith   PetscErrorCode ierr;
221827bd09bSSatish Balay   gs_id *gs;
222330ea6edSBarry Smith   gs = (gs_id *) malloc(sizeof(gs_id));
223f1ed62a8SBarry Smith   ierr = PetscMemzero(gs,sizeof(gs_id));CHKERRABORT(PETSC_COMM_WORLD,ierr);
224827bd09bSSatish Balay   return(gs);
225827bd09bSSatish Balay }
226827bd09bSSatish Balay 
227f1ed62a8SBarry Smith /******************************************************************************/
22852f87cdaSBarry Smith static gs_id * gsi_check_args(PetscInt *in_elms, PetscInt nel, PetscInt level)
229827bd09bSSatish Balay {
23052f87cdaSBarry Smith    PetscInt i, j, k, t2;
23152f87cdaSBarry Smith   PetscInt *companion, *elms, *unique, *iptr;
23252f87cdaSBarry Smith   PetscInt num_local=0, *num_to_reduce, **local_reduce;
23352f87cdaSBarry Smith   PetscInt oprs[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_MIN,GL_B_AND};
23452f87cdaSBarry Smith   PetscInt vals[sizeof(oprs)/sizeof(oprs[0])-1];
23552f87cdaSBarry Smith   PetscInt work[sizeof(oprs)/sizeof(oprs[0])-1];
236827bd09bSSatish Balay   gs_id *gs;
237d1528f56SBarry Smith   PetscErrorCode ierr;
238827bd09bSSatish Balay 
239827bd09bSSatish Balay 
240827bd09bSSatish Balay   if (!in_elms)
241388eb383SBarry Smith     {SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"elms point to nothing!!!\n");}
242827bd09bSSatish Balay 
243827bd09bSSatish Balay   if (nel<0)
244388eb383SBarry Smith     {SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"can't have fewer than 0 elms!!!\n");}
245827bd09bSSatish Balay 
246827bd09bSSatish Balay   if (nel==0)
247f1ed62a8SBarry Smith     {ierr = PetscInfo(0,"I don't have any elements!!!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);}
248827bd09bSSatish Balay 
249827bd09bSSatish Balay   /* get space for gs template */
250827bd09bSSatish Balay   gs = gsi_new();
251827bd09bSSatish Balay   gs->id = ++num_gs_ids;
252827bd09bSSatish Balay 
253827bd09bSSatish Balay   /* hmt 6.4.99                                            */
254827bd09bSSatish Balay   /* caller can set global ids that don't participate to 0 */
255827bd09bSSatish Balay   /* gs_init ignores all zeros in elm list                 */
256827bd09bSSatish Balay   /* negative global ids are still invalid                 */
257827bd09bSSatish Balay   for (i=j=0;i<nel;i++)
258827bd09bSSatish Balay     {if (in_elms[i]!=0) {j++;}}
259827bd09bSSatish Balay 
260827bd09bSSatish Balay   k=nel; nel=j;
261827bd09bSSatish Balay 
262827bd09bSSatish Balay   /* copy over in_elms list and create inverse map */
26352f87cdaSBarry Smith   elms = (PetscInt*) malloc((nel+1)*sizeof(PetscInt));
26452f87cdaSBarry Smith   companion = (PetscInt*) malloc(nel*sizeof(PetscInt));
2651d7d0905SBarry Smith 
266827bd09bSSatish Balay   for (i=j=0;i<k;i++)
267827bd09bSSatish Balay     {
268827bd09bSSatish Balay       if (in_elms[i]!=0)
269827bd09bSSatish Balay         {elms[j] = in_elms[i]; companion[j++] = i;}
270827bd09bSSatish Balay     }
271827bd09bSSatish Balay 
272827bd09bSSatish Balay   if (j!=nel)
273388eb383SBarry Smith     {SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"nel j mismatch!\n");}
274827bd09bSSatish Balay 
275827bd09bSSatish Balay   /* pre-pass ... check to see if sorted */
276827bd09bSSatish Balay   elms[nel] = INT_MAX;
277827bd09bSSatish Balay   iptr = elms;
278827bd09bSSatish Balay   unique = elms+1;
279827bd09bSSatish Balay   j=0;
280827bd09bSSatish Balay   while (*iptr!=INT_MAX)
281827bd09bSSatish Balay     {
282827bd09bSSatish Balay       if (*iptr++>*unique++)
283827bd09bSSatish Balay         {j=1; break;}
284827bd09bSSatish Balay     }
285827bd09bSSatish Balay 
286827bd09bSSatish Balay   /* set up inverse map */
287827bd09bSSatish Balay   if (j)
288827bd09bSSatish Balay     {
289f1ed62a8SBarry Smith       ierr = PetscInfo(0,"gsi_check_args() :: elm list *not* sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);
290f1ed62a8SBarry Smith       ierr = SMI_sort((void*)elms, (void*)companion, nel, SORT_INTEGER);CHKERRABORT(PETSC_COMM_WORLD,ierr);
291827bd09bSSatish Balay     }
292827bd09bSSatish Balay   else
293f1ed62a8SBarry Smith     {ierr = PetscInfo(0,"gsi_check_args() :: elm list sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);}
294827bd09bSSatish Balay   elms[nel] = INT_MIN;
295827bd09bSSatish Balay 
296827bd09bSSatish Balay   /* first pass */
297827bd09bSSatish Balay   /* determine number of unique elements, check pd */
298827bd09bSSatish Balay   for (i=k=0;i<nel;i+=j)
299827bd09bSSatish Balay     {
300827bd09bSSatish Balay       t2 = elms[i];
301827bd09bSSatish Balay       j=++i;
302827bd09bSSatish Balay 
303827bd09bSSatish Balay       /* clump 'em for now */
304827bd09bSSatish Balay       while (elms[j]==t2) {j++;}
305827bd09bSSatish Balay 
306827bd09bSSatish Balay       /* how many together and num local */
307827bd09bSSatish Balay       if (j-=i)
308827bd09bSSatish Balay         {num_local++; k+=j;}
309827bd09bSSatish Balay     }
310827bd09bSSatish Balay 
311827bd09bSSatish Balay   /* how many unique elements? */
312827bd09bSSatish Balay   gs->repeats=k;
313827bd09bSSatish Balay   gs->nel = nel-k;
314827bd09bSSatish Balay 
315827bd09bSSatish Balay 
316827bd09bSSatish Balay   /* number of repeats? */
317827bd09bSSatish Balay   gs->num_local = num_local;
318827bd09bSSatish Balay   num_local+=2;
31952f87cdaSBarry Smith   gs->local_reduce=local_reduce=(PetscInt **)malloc(num_local*sizeof(PetscInt*));
32052f87cdaSBarry Smith   gs->num_local_reduce=num_to_reduce=(PetscInt*) malloc(num_local*sizeof(PetscInt));
321827bd09bSSatish Balay 
32252f87cdaSBarry Smith   unique = (PetscInt*) malloc((gs->nel+1)*sizeof(PetscInt));
323827bd09bSSatish Balay   gs->elms = unique;
324827bd09bSSatish Balay   gs->nel_total = nel;
325827bd09bSSatish Balay   gs->local_elms = elms;
326827bd09bSSatish Balay   gs->companion = companion;
327827bd09bSSatish Balay 
328827bd09bSSatish Balay   /* compess map as well as keep track of local ops */
329827bd09bSSatish Balay   for (num_local=i=j=0;i<gs->nel;i++)
330827bd09bSSatish Balay     {
331827bd09bSSatish Balay       k=j;
332827bd09bSSatish Balay       t2 = unique[i] = elms[j];
333827bd09bSSatish Balay       companion[i] = companion[j];
334827bd09bSSatish Balay 
335827bd09bSSatish Balay       while (elms[j]==t2) {j++;}
336827bd09bSSatish Balay 
337827bd09bSSatish Balay       if ((t2=(j-k))>1)
338827bd09bSSatish Balay         {
339827bd09bSSatish Balay           /* number together */
340827bd09bSSatish Balay           num_to_reduce[num_local] = t2++;
34152f87cdaSBarry Smith           iptr = local_reduce[num_local++] = (PetscInt*)malloc(t2*sizeof(PetscInt));
342827bd09bSSatish Balay 
343827bd09bSSatish Balay           /* to use binary searching don't remap until we check intersection */
344827bd09bSSatish Balay           *iptr++ = i;
345827bd09bSSatish Balay 
346827bd09bSSatish Balay           /* note that we're skipping the first one */
347827bd09bSSatish Balay           while (++k<j)
348827bd09bSSatish Balay             {*(iptr++) = companion[k];}
349827bd09bSSatish Balay           *iptr = -1;
350827bd09bSSatish Balay         }
351827bd09bSSatish Balay     }
352827bd09bSSatish Balay 
353827bd09bSSatish Balay   /* sentinel for ngh_buf */
354827bd09bSSatish Balay   unique[gs->nel]=INT_MAX;
355827bd09bSSatish Balay 
356827bd09bSSatish Balay   /* for two partition sort hack */
357827bd09bSSatish Balay   num_to_reduce[num_local] = 0;
358827bd09bSSatish Balay   local_reduce[num_local] = NULL;
359827bd09bSSatish Balay   num_to_reduce[++num_local] = 0;
360827bd09bSSatish Balay   local_reduce[num_local] = NULL;
361827bd09bSSatish Balay 
362827bd09bSSatish Balay   /* load 'em up */
363827bd09bSSatish Balay   /* note one extra to hold NON_UNIFORM flag!!! */
364827bd09bSSatish Balay   vals[2] = vals[1] = vals[0] = nel;
365827bd09bSSatish Balay   if (gs->nel>0)
366827bd09bSSatish Balay     {
3671d7d0905SBarry Smith        vals[3] = unique[0];
3681d7d0905SBarry Smith        vals[4] = unique[gs->nel-1];
369827bd09bSSatish Balay     }
370827bd09bSSatish Balay   else
371827bd09bSSatish Balay     {
3721d7d0905SBarry Smith        vals[3] = INT_MAX;
3731d7d0905SBarry Smith        vals[4] = INT_MIN;
374827bd09bSSatish Balay     }
375827bd09bSSatish Balay   vals[5] = level;
376827bd09bSSatish Balay   vals[6] = num_gs_ids;
377827bd09bSSatish Balay 
378827bd09bSSatish Balay   /* GLOBAL: send 'em out */
379f1ed62a8SBarry Smith   ierr = giop(vals,work,sizeof(oprs)/sizeof(oprs[0])-1,oprs);CHKERRABORT(PETSC_COMM_WORLD,ierr);
380827bd09bSSatish Balay 
381827bd09bSSatish Balay   /* must be semi-pos def - only pairwise depends on this */
382827bd09bSSatish Balay   /* LATER - remove this restriction */
383827bd09bSSatish Balay   if (vals[3]<0)
384388eb383SBarry Smith     {SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system not semi-pos def \n");}
385827bd09bSSatish Balay 
386827bd09bSSatish Balay   if (vals[4]==INT_MAX)
387388eb383SBarry Smith     {SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system ub too large !\n");}
388827bd09bSSatish Balay 
389827bd09bSSatish Balay   gs->nel_min = vals[0];
390827bd09bSSatish Balay   gs->nel_max = vals[1];
391827bd09bSSatish Balay   gs->nel_sum = vals[2];
392827bd09bSSatish Balay   gs->gl_min  = vals[3];
393827bd09bSSatish Balay   gs->gl_max  = vals[4];
394827bd09bSSatish Balay   gs->negl    = vals[4]-vals[3]+1;
395827bd09bSSatish Balay 
396827bd09bSSatish Balay   if (gs->negl<=0)
397388eb383SBarry Smith     {SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system empty or neg :: %d\n");}
398827bd09bSSatish Balay 
399827bd09bSSatish Balay   /* LATER :: add level == -1 -> program selects level */
400827bd09bSSatish Balay   if (vals[5]<0)
401827bd09bSSatish Balay     {vals[5]=0;}
402827bd09bSSatish Balay   else if (vals[5]>num_nodes)
403827bd09bSSatish Balay     {vals[5]=num_nodes;}
404827bd09bSSatish Balay   gs->level = vals[5];
405827bd09bSSatish Balay 
406827bd09bSSatish Balay   return(gs);
407827bd09bSSatish Balay }
408827bd09bSSatish Balay 
409f1ed62a8SBarry Smith /******************************************************************************/
4100924e98cSBarry Smith static PetscErrorCode gsi_via_bit_mask(gs_id *gs)
411827bd09bSSatish Balay {
41252f87cdaSBarry Smith    PetscInt i, nel, *elms;
41352f87cdaSBarry Smith   PetscInt t1;
41452f87cdaSBarry Smith   PetscInt **reduce;
41552f87cdaSBarry Smith   PetscInt *map;
416f1ed62a8SBarry Smith   PetscErrorCode ierr;
417827bd09bSSatish Balay 
418f1ed62a8SBarry Smith   PetscFunctionBegin;
419827bd09bSSatish Balay   /* totally local removes ... ct_bits == 0 */
420827bd09bSSatish Balay   get_ngh_buf(gs);
421827bd09bSSatish Balay 
422827bd09bSSatish Balay   if (gs->level)
423827bd09bSSatish Balay     {set_pairwise(gs);}
424827bd09bSSatish Balay 
425827bd09bSSatish Balay   if (gs->max_left_over)
426827bd09bSSatish Balay     {set_tree(gs);}
427827bd09bSSatish Balay 
428827bd09bSSatish Balay   /* intersection local and pairwise/tree? */
429827bd09bSSatish Balay   gs->num_local_total = gs->num_local;
430827bd09bSSatish Balay   gs->gop_local_reduce = gs->local_reduce;
431827bd09bSSatish Balay   gs->num_gop_local_reduce = gs->num_local_reduce;
432827bd09bSSatish Balay 
433827bd09bSSatish Balay   map = gs->companion;
434827bd09bSSatish Balay 
435827bd09bSSatish Balay   /* is there any local compression */
436d890fc11SSatish Balay   if (!gs->num_local) {
437827bd09bSSatish Balay     gs->local_strength = NONE;
438827bd09bSSatish Balay     gs->num_local_gop = 0;
439d890fc11SSatish Balay   } else {
440827bd09bSSatish Balay       /* ok find intersection */
441827bd09bSSatish Balay       map = gs->companion;
442827bd09bSSatish Balay       reduce = gs->local_reduce;
443827bd09bSSatish Balay       for (i=0, t1=0; i<gs->num_local; i++, reduce++)
444827bd09bSSatish Balay         {
445827bd09bSSatish Balay           if ((ivec_binary_search(**reduce,gs->pw_elm_list,gs->len_pw_list)>=0)
446827bd09bSSatish Balay               ||
447827bd09bSSatish Balay               ivec_binary_search(**reduce,gs->tree_map_in,gs->tree_map_sz)>=0)
448827bd09bSSatish Balay             {
449827bd09bSSatish Balay               t1++;
450f1ed62a8SBarry Smith               if (gs->num_local_reduce[i]<=0) SETERRQ(PETSC_ERR_PLIB,"nobody in list?");
451827bd09bSSatish Balay               gs->num_local_reduce[i] *= -1;
452827bd09bSSatish Balay             }
453827bd09bSSatish Balay            **reduce=map[**reduce];
454827bd09bSSatish Balay         }
455827bd09bSSatish Balay 
456827bd09bSSatish Balay       /* intersection is empty */
457827bd09bSSatish Balay       if (!t1)
458827bd09bSSatish Balay         {
459827bd09bSSatish Balay           gs->local_strength = FULL;
460827bd09bSSatish Balay           gs->num_local_gop = 0;
461827bd09bSSatish Balay         }
462827bd09bSSatish Balay       /* intersection not empty */
463827bd09bSSatish Balay       else
464827bd09bSSatish Balay         {
465827bd09bSSatish Balay           gs->local_strength = PARTIAL;
466f1ed62a8SBarry Smith           ierr = SMI_sort((void*)gs->num_local_reduce, (void*)gs->local_reduce, gs->num_local + 1, SORT_INT_PTR);CHKERRQ(ierr);
467827bd09bSSatish Balay 
468827bd09bSSatish Balay           gs->num_local_gop = t1;
469827bd09bSSatish Balay           gs->num_local_total =  gs->num_local;
470827bd09bSSatish Balay           gs->num_local    -= t1;
471827bd09bSSatish Balay           gs->gop_local_reduce = gs->local_reduce;
472827bd09bSSatish Balay           gs->num_gop_local_reduce = gs->num_local_reduce;
473827bd09bSSatish Balay 
474827bd09bSSatish Balay           for (i=0; i<t1; i++)
475827bd09bSSatish Balay             {
476f1ed62a8SBarry Smith               if (gs->num_gop_local_reduce[i]>=0) SETERRQ(PETSC_ERR_PLIB,"they aren't negative?");
477827bd09bSSatish Balay               gs->num_gop_local_reduce[i] *= -1;
478827bd09bSSatish Balay               gs->local_reduce++;
479827bd09bSSatish Balay               gs->num_local_reduce++;
480827bd09bSSatish Balay             }
481827bd09bSSatish Balay           gs->local_reduce++;
482827bd09bSSatish Balay           gs->num_local_reduce++;
483827bd09bSSatish Balay         }
484827bd09bSSatish Balay     }
485827bd09bSSatish Balay 
486827bd09bSSatish Balay   elms = gs->pw_elm_list;
487827bd09bSSatish Balay   nel  = gs->len_pw_list;
488827bd09bSSatish Balay   for (i=0; i<nel; i++)
489827bd09bSSatish Balay     {elms[i] = map[elms[i]];}
490827bd09bSSatish Balay 
491827bd09bSSatish Balay   elms = gs->tree_map_in;
492827bd09bSSatish Balay   nel  = gs->tree_map_sz;
493827bd09bSSatish Balay   for (i=0; i<nel; i++)
494827bd09bSSatish Balay     {elms[i] = map[elms[i]];}
495827bd09bSSatish Balay 
496827bd09bSSatish Balay   /* clean up */
497a501084fSBarry Smith   free((void*) gs->local_elms);
498a501084fSBarry Smith   free((void*) gs->companion);
499a501084fSBarry Smith   free((void*) gs->elms);
500a501084fSBarry Smith   free((void*) gs->ngh_buf);
501827bd09bSSatish Balay   gs->local_elms = gs->companion = gs->elms = gs->ngh_buf = NULL;
5023fdc5746SBarry Smith   PetscFunctionReturn(0);
503827bd09bSSatish Balay }
504827bd09bSSatish Balay 
505f1ed62a8SBarry Smith /******************************************************************************/
50652f87cdaSBarry Smith static PetscErrorCode place_in_tree( PetscInt elm)
507827bd09bSSatish Balay {
50852f87cdaSBarry Smith    PetscInt *tp, n;
509827bd09bSSatish Balay 
5103fdc5746SBarry Smith   PetscFunctionBegin;
511827bd09bSSatish Balay   if (ntree==tree_buf_sz)
512827bd09bSSatish Balay     {
513827bd09bSSatish Balay       if (tree_buf_sz)
514827bd09bSSatish Balay         {
515827bd09bSSatish Balay           tp = tree_buf;
516827bd09bSSatish Balay           n = tree_buf_sz;
517827bd09bSSatish Balay           tree_buf_sz<<=1;
51852f87cdaSBarry Smith           tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt));
519827bd09bSSatish Balay           ivec_copy(tree_buf,tp,n);
520a501084fSBarry Smith           free(tp);
521827bd09bSSatish Balay         }
522827bd09bSSatish Balay       else
523827bd09bSSatish Balay         {
524827bd09bSSatish Balay           tree_buf_sz = TREE_BUF_SZ;
52552f87cdaSBarry Smith           tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt));
526827bd09bSSatish Balay         }
527827bd09bSSatish Balay     }
528827bd09bSSatish Balay 
529827bd09bSSatish Balay   tree_buf[ntree++] = elm;
5303fdc5746SBarry Smith   PetscFunctionReturn(0);
531827bd09bSSatish Balay }
532827bd09bSSatish Balay 
533f1ed62a8SBarry Smith /******************************************************************************/
5340924e98cSBarry Smith static PetscErrorCode get_ngh_buf(gs_id *gs)
535827bd09bSSatish Balay {
53652f87cdaSBarry Smith    PetscInt i, j, npw=0, ntree_map=0;
53752f87cdaSBarry Smith   PetscInt p_mask_size, ngh_buf_size, buf_size;
53852f87cdaSBarry Smith   PetscInt *p_mask, *sh_proc_mask, *pw_sh_proc_mask;
53952f87cdaSBarry Smith   PetscInt *ngh_buf, *buf1, *buf2;
54052f87cdaSBarry Smith   PetscInt offset, per_load, num_loads, or_ct, start, end;
54152f87cdaSBarry Smith   PetscInt *ptr1, *ptr2, i_start, negl, nel, *elms;
54252f87cdaSBarry Smith   PetscInt oper=GL_B_OR;
54352f87cdaSBarry Smith   PetscInt *ptr3, *t_mask, level, ct1, ct2;
544f1ed62a8SBarry Smith   PetscErrorCode ierr;
545827bd09bSSatish Balay 
5463fdc5746SBarry Smith   PetscFunctionBegin;
547827bd09bSSatish Balay   /* to make life easier */
548827bd09bSSatish Balay   nel   = gs->nel;
549827bd09bSSatish Balay   elms  = gs->elms;
550827bd09bSSatish Balay   level = gs->level;
551827bd09bSSatish Balay 
552827bd09bSSatish Balay   /* det #bytes needed for processor bit masks and init w/mask cor. to my_id */
55352f87cdaSBarry Smith   p_mask = (PetscInt*) malloc(p_mask_size=len_bit_mask(num_nodes));
554f1ed62a8SBarry Smith   ierr = set_bit_mask(p_mask,p_mask_size,my_id);CHKERRQ(ierr);
555827bd09bSSatish Balay 
556827bd09bSSatish Balay   /* allocate space for masks and info bufs */
55752f87cdaSBarry Smith   gs->nghs = sh_proc_mask = (PetscInt*) malloc(p_mask_size);
55852f87cdaSBarry Smith   gs->pw_nghs = pw_sh_proc_mask = (PetscInt*) malloc(p_mask_size);
559827bd09bSSatish Balay   gs->ngh_buf_sz = ngh_buf_size = p_mask_size*nel;
56052f87cdaSBarry Smith   t_mask = (PetscInt*) malloc(p_mask_size);
56152f87cdaSBarry Smith   gs->ngh_buf = ngh_buf = (PetscInt*) malloc(ngh_buf_size);
562827bd09bSSatish Balay 
563827bd09bSSatish Balay   /* comm buffer size ... memory usage bounded by ~2*msg_buf */
564827bd09bSSatish Balay   /* had thought I could exploit rendezvous threshold */
565827bd09bSSatish Balay 
566827bd09bSSatish Balay   /* default is one pass */
567827bd09bSSatish Balay   per_load = negl  = gs->negl;
568827bd09bSSatish Balay   gs->num_loads = num_loads = 1;
569827bd09bSSatish Balay   i=p_mask_size*negl;
570827bd09bSSatish Balay 
571827bd09bSSatish Balay   /* possible overflow on buffer size */
572827bd09bSSatish Balay   /* overflow hack                    */
573827bd09bSSatish Balay   if (i<0) {i=INT_MAX;}
574827bd09bSSatish Balay 
57539945688SSatish Balay   buf_size = PetscMin(msg_buf,i);
576827bd09bSSatish Balay 
577827bd09bSSatish Balay   /* can we do it? */
578f1ed62a8SBarry Smith   if (p_mask_size>buf_size) SETERRQ2(PETSC_ERR_PLIB,"get_ngh_buf() :: buf<pms :: %d>%d\n",p_mask_size,buf_size);
579827bd09bSSatish Balay 
580827bd09bSSatish Balay   /* get giop buf space ... make *only* one malloc */
58152f87cdaSBarry Smith   buf1 = (PetscInt*) malloc(buf_size<<1);
582827bd09bSSatish Balay 
583827bd09bSSatish Balay   /* more than one gior exchange needed? */
584827bd09bSSatish Balay   if (buf_size!=i)
585827bd09bSSatish Balay     {
586827bd09bSSatish Balay       per_load = buf_size/p_mask_size;
587827bd09bSSatish Balay       buf_size = per_load*p_mask_size;
588827bd09bSSatish Balay       gs->num_loads = num_loads = negl/per_load + (negl%per_load>0);
589827bd09bSSatish Balay     }
590827bd09bSSatish Balay 
591827bd09bSSatish Balay 
592827bd09bSSatish Balay   /* convert buf sizes from #bytes to #ints - 32 bit only! */
593a501084fSBarry Smith   p_mask_size/=sizeof(PetscInt); ngh_buf_size/=sizeof(PetscInt); buf_size/=sizeof(PetscInt);
594827bd09bSSatish Balay 
595827bd09bSSatish Balay   /* find giop work space */
596827bd09bSSatish Balay   buf2 = buf1+buf_size;
597827bd09bSSatish Balay 
598827bd09bSSatish Balay   /* hold #ints needed for processor masks */
599827bd09bSSatish Balay   gs->mask_sz=p_mask_size;
600827bd09bSSatish Balay 
601827bd09bSSatish Balay   /* init buffers */
602f1ed62a8SBarry Smith   ierr = ivec_zero(sh_proc_mask,p_mask_size);CHKERRQ(ierr);
603f1ed62a8SBarry Smith   ierr = ivec_zero(pw_sh_proc_mask,p_mask_size);CHKERRQ(ierr);
604f1ed62a8SBarry Smith   ierr = ivec_zero(ngh_buf,ngh_buf_size);CHKERRQ(ierr);
605827bd09bSSatish Balay 
606827bd09bSSatish Balay   /* HACK reset tree info */
607827bd09bSSatish Balay   tree_buf=NULL;
608827bd09bSSatish Balay   tree_buf_sz=ntree=0;
609827bd09bSSatish Balay 
610827bd09bSSatish Balay   /* ok do it */
611827bd09bSSatish Balay   for (ptr1=ngh_buf,ptr2=elms,end=gs->gl_min,or_ct=i=0; or_ct<num_loads; or_ct++)
612827bd09bSSatish Balay     {
613827bd09bSSatish Balay       /* identity for bitwise or is 000...000 */
614827bd09bSSatish Balay       ivec_zero(buf1,buf_size);
615827bd09bSSatish Balay 
616827bd09bSSatish Balay       /* load msg buffer */
617827bd09bSSatish Balay       for (start=end,end+=per_load,i_start=i; (offset=*ptr2)<end; i++, ptr2++)
618827bd09bSSatish Balay         {
619827bd09bSSatish Balay           offset = (offset-start)*p_mask_size;
620827bd09bSSatish Balay           ivec_copy(buf1+offset,p_mask,p_mask_size);
621827bd09bSSatish Balay         }
622827bd09bSSatish Balay 
623827bd09bSSatish Balay       /* GLOBAL: pass buffer */
624f1ed62a8SBarry Smith       ierr = giop(buf1,buf2,buf_size,&oper);CHKERRQ(ierr);
625827bd09bSSatish Balay 
626827bd09bSSatish Balay 
627827bd09bSSatish Balay       /* unload buffer into ngh_buf */
628827bd09bSSatish Balay       ptr2=(elms+i_start);
629827bd09bSSatish Balay       for(ptr3=buf1,j=start; j<end; ptr3+=p_mask_size,j++)
630827bd09bSSatish Balay         {
631827bd09bSSatish Balay           /* I own it ... may have to pairwise it */
632827bd09bSSatish Balay           if (j==*ptr2)
633827bd09bSSatish Balay             {
634827bd09bSSatish Balay               /* do i share it w/anyone? */
635a501084fSBarry Smith               ct1 = ct_bits((char *)ptr3,p_mask_size*sizeof(PetscInt));
636827bd09bSSatish Balay               /* guess not */
637827bd09bSSatish Balay               if (ct1<2)
638827bd09bSSatish Balay                 {ptr2++; ptr1+=p_mask_size; continue;}
639827bd09bSSatish Balay 
640827bd09bSSatish Balay               /* i do ... so keep info and turn off my bit */
641827bd09bSSatish Balay               ivec_copy(ptr1,ptr3,p_mask_size);
642f1ed62a8SBarry Smith               ierr = ivec_xor(ptr1,p_mask,p_mask_size);CHKERRQ(ierr);
643f1ed62a8SBarry Smith               ierr = ivec_or(sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr);
644827bd09bSSatish Balay 
645827bd09bSSatish Balay               /* is it to be done pairwise? */
646827bd09bSSatish Balay               if (--ct1<=level)
647827bd09bSSatish Balay                 {
648827bd09bSSatish Balay                   npw++;
649827bd09bSSatish Balay 
650827bd09bSSatish Balay                   /* turn on high bit to indicate pw need to process */
651827bd09bSSatish Balay                   *ptr2++ |= TOP_BIT;
652f1ed62a8SBarry Smith                   ierr = ivec_or(pw_sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr);
653827bd09bSSatish Balay                   ptr1+=p_mask_size;
654827bd09bSSatish Balay                   continue;
655827bd09bSSatish Balay                 }
656827bd09bSSatish Balay 
657827bd09bSSatish Balay               /* get set for next and note that I have a tree contribution */
658827bd09bSSatish Balay               /* could save exact elm index for tree here -> save a search */
659827bd09bSSatish Balay               ptr2++; ptr1+=p_mask_size; ntree_map++;
660827bd09bSSatish Balay             }
661827bd09bSSatish Balay           /* i don't but still might be involved in tree */
662827bd09bSSatish Balay           else
663827bd09bSSatish Balay             {
664827bd09bSSatish Balay 
665827bd09bSSatish Balay               /* shared by how many? */
666a501084fSBarry Smith               ct1 = ct_bits((char *)ptr3,p_mask_size*sizeof(PetscInt));
667827bd09bSSatish Balay 
668827bd09bSSatish Balay               /* none! */
669f1ed62a8SBarry Smith               if (ct1<2) continue;
670827bd09bSSatish Balay 
671827bd09bSSatish Balay               /* is it going to be done pairwise? but not by me of course!*/
672f1ed62a8SBarry Smith               if (--ct1<=level) continue;
673827bd09bSSatish Balay             }
674827bd09bSSatish Balay           /* LATER we're going to have to process it NOW */
675827bd09bSSatish Balay           /* nope ... tree it */
676f1ed62a8SBarry Smith           ierr = place_in_tree(j);CHKERRQ(ierr);
677827bd09bSSatish Balay         }
678827bd09bSSatish Balay     }
679827bd09bSSatish Balay 
680a501084fSBarry Smith   free((void*)t_mask);
681a501084fSBarry Smith   free((void*)buf1);
682827bd09bSSatish Balay 
683827bd09bSSatish Balay   gs->len_pw_list=npw;
684a501084fSBarry Smith   gs->num_nghs = ct_bits((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt));
685827bd09bSSatish Balay 
686827bd09bSSatish Balay   /* expand from bit mask list to int list and save ngh list */
68752f87cdaSBarry Smith   gs->nghs = (PetscInt*) malloc(gs->num_nghs * sizeof(PetscInt));
688a501084fSBarry Smith   bm_to_proc((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt),gs->nghs);
689827bd09bSSatish Balay 
690a501084fSBarry Smith   gs->num_pw_nghs = ct_bits((char *)pw_sh_proc_mask,p_mask_size*sizeof(PetscInt));
691827bd09bSSatish Balay 
692827bd09bSSatish Balay   oper = GL_MAX;
693827bd09bSSatish Balay   ct1 = gs->num_nghs;
694f1ed62a8SBarry Smith   ierr = giop(&ct1,&ct2,1,&oper);CHKERRQ(ierr);
695827bd09bSSatish Balay   gs->max_nghs = ct1;
696827bd09bSSatish Balay 
697827bd09bSSatish Balay   gs->tree_map_sz  = ntree_map;
698827bd09bSSatish Balay   gs->max_left_over=ntree;
699827bd09bSSatish Balay 
700a501084fSBarry Smith   free((void*)p_mask);
701a501084fSBarry Smith   free((void*)sh_proc_mask);
7023fdc5746SBarry Smith   PetscFunctionReturn(0);
703827bd09bSSatish Balay }
704827bd09bSSatish Balay 
705f1ed62a8SBarry Smith /******************************************************************************/
7060924e98cSBarry Smith static PetscErrorCode set_pairwise(gs_id *gs)
707827bd09bSSatish Balay {
70852f87cdaSBarry Smith    PetscInt i, j;
70952f87cdaSBarry Smith   PetscInt p_mask_size;
71052f87cdaSBarry Smith   PetscInt *p_mask, *sh_proc_mask, *tmp_proc_mask;
71152f87cdaSBarry Smith   PetscInt *ngh_buf, *buf2;
71252f87cdaSBarry Smith   PetscInt offset;
71352f87cdaSBarry Smith   PetscInt *msg_list, *msg_size, **msg_nodes, nprs;
71452f87cdaSBarry Smith   PetscInt *pairwise_elm_list, len_pair_list=0;
71552f87cdaSBarry Smith   PetscInt *iptr, t1, i_start, nel, *elms;
71652f87cdaSBarry Smith   PetscInt ct;
717f1ed62a8SBarry Smith   PetscErrorCode ierr;
718827bd09bSSatish Balay 
7193fdc5746SBarry Smith   PetscFunctionBegin;
720827bd09bSSatish Balay   /* to make life easier */
721827bd09bSSatish Balay   nel  = gs->nel;
722827bd09bSSatish Balay   elms = gs->elms;
723827bd09bSSatish Balay   ngh_buf = gs->ngh_buf;
724827bd09bSSatish Balay   sh_proc_mask  = gs->pw_nghs;
725827bd09bSSatish Balay 
726827bd09bSSatish Balay   /* need a few temp masks */
727827bd09bSSatish Balay   p_mask_size   = len_bit_mask(num_nodes);
72852f87cdaSBarry Smith   p_mask        = (PetscInt*) malloc(p_mask_size);
72952f87cdaSBarry Smith   tmp_proc_mask = (PetscInt*) malloc(p_mask_size);
730827bd09bSSatish Balay 
731827bd09bSSatish Balay   /* set mask to my my_id's bit mask */
732f1ed62a8SBarry Smith   ierr = set_bit_mask(p_mask,p_mask_size,my_id);CHKERRQ(ierr);
733827bd09bSSatish Balay 
734a501084fSBarry Smith   p_mask_size /= sizeof(PetscInt);
735827bd09bSSatish Balay 
736827bd09bSSatish Balay   len_pair_list=gs->len_pw_list;
73752f87cdaSBarry Smith   gs->pw_elm_list=pairwise_elm_list=(PetscInt*)malloc((len_pair_list+1)*sizeof(PetscInt));
738827bd09bSSatish Balay 
739827bd09bSSatish Balay   /* how many processors (nghs) do we have to exchange with? */
740a501084fSBarry Smith   nprs=gs->num_pairs=ct_bits((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt));
741827bd09bSSatish Balay 
742827bd09bSSatish Balay 
743827bd09bSSatish Balay   /* allocate space for gs_gop() info */
74452f87cdaSBarry Smith   gs->pair_list = msg_list = (PetscInt *)  malloc(sizeof(PetscInt)*nprs);
74552f87cdaSBarry Smith   gs->msg_sizes = msg_size  = (PetscInt *)  malloc(sizeof(PetscInt)*nprs);
74652f87cdaSBarry Smith   gs->node_list = msg_nodes = (PetscInt **) malloc(sizeof(PetscInt*)*(nprs+1));
747827bd09bSSatish Balay 
748827bd09bSSatish Balay   /* init msg_size list */
749f1ed62a8SBarry Smith   ierr = ivec_zero(msg_size,nprs);CHKERRQ(ierr);
750827bd09bSSatish Balay 
751827bd09bSSatish Balay   /* expand from bit mask list to int list */
752f1ed62a8SBarry Smith   ierr = bm_to_proc((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt),msg_list);CHKERRQ(ierr);
753827bd09bSSatish Balay 
754827bd09bSSatish Balay   /* keep list of elements being handled pairwise */
755827bd09bSSatish Balay   for (i=j=0;i<nel;i++)
756827bd09bSSatish Balay     {
757827bd09bSSatish Balay       if (elms[i] & TOP_BIT)
758827bd09bSSatish Balay         {elms[i] ^= TOP_BIT; pairwise_elm_list[j++] = i;}
759827bd09bSSatish Balay     }
760827bd09bSSatish Balay   pairwise_elm_list[j] = -1;
761827bd09bSSatish Balay 
762a501084fSBarry Smith   gs->msg_ids_out = (MPI_Request *)  malloc(sizeof(MPI_Request)*(nprs+1));
763827bd09bSSatish Balay   gs->msg_ids_out[nprs] = MPI_REQUEST_NULL;
764a501084fSBarry Smith   gs->msg_ids_in = (MPI_Request *)  malloc(sizeof(MPI_Request)*(nprs+1));
765827bd09bSSatish Balay   gs->msg_ids_in[nprs] = MPI_REQUEST_NULL;
766a501084fSBarry Smith   gs->pw_vals = (PetscScalar *) malloc(sizeof(PetscScalar)*len_pair_list*vec_sz);
767827bd09bSSatish Balay 
768827bd09bSSatish Balay   /* find who goes to each processor */
769827bd09bSSatish Balay   for (i_start=i=0;i<nprs;i++)
770827bd09bSSatish Balay     {
771827bd09bSSatish Balay       /* processor i's mask */
772f1ed62a8SBarry Smith       ierr = set_bit_mask(p_mask,p_mask_size*sizeof(PetscInt),msg_list[i]);CHKERRQ(ierr);
773827bd09bSSatish Balay 
774827bd09bSSatish Balay       /* det # going to processor i */
775827bd09bSSatish Balay       for (ct=j=0;j<len_pair_list;j++)
776827bd09bSSatish Balay         {
777827bd09bSSatish Balay           buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size);
778f1ed62a8SBarry Smith           ierr = ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr);
779a501084fSBarry Smith           if (ct_bits((char *)tmp_proc_mask,p_mask_size*sizeof(PetscInt)))
780827bd09bSSatish Balay             {ct++;}
781827bd09bSSatish Balay         }
782827bd09bSSatish Balay       msg_size[i] = ct;
78339945688SSatish Balay       i_start = PetscMax(i_start,ct);
784827bd09bSSatish Balay 
785827bd09bSSatish Balay       /*space to hold nodes in message to first neighbor */
78652f87cdaSBarry Smith       msg_nodes[i] = iptr = (PetscInt*) malloc(sizeof(PetscInt)*(ct+1));
787827bd09bSSatish Balay 
788827bd09bSSatish Balay       for (j=0;j<len_pair_list;j++)
789827bd09bSSatish Balay         {
790827bd09bSSatish Balay           buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size);
791f1ed62a8SBarry Smith           ierr = ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr);
792a501084fSBarry Smith           if (ct_bits((char *)tmp_proc_mask,p_mask_size*sizeof(PetscInt)))
793827bd09bSSatish Balay             {*iptr++ = j;}
794827bd09bSSatish Balay         }
795827bd09bSSatish Balay       *iptr = -1;
796827bd09bSSatish Balay     }
797827bd09bSSatish Balay   msg_nodes[nprs] = NULL;
798827bd09bSSatish Balay 
799827bd09bSSatish Balay   j=gs->loc_node_pairs=i_start;
800827bd09bSSatish Balay   t1 = GL_MAX;
801f1ed62a8SBarry Smith   ierr = giop(&i_start,&offset,1,&t1);CHKERRQ(ierr);
802827bd09bSSatish Balay   gs->max_node_pairs = i_start;
803827bd09bSSatish Balay 
804827bd09bSSatish Balay   i_start=j;
805827bd09bSSatish Balay   t1 = GL_MIN;
806f1ed62a8SBarry Smith   ierr = giop(&i_start,&offset,1,&t1);CHKERRQ(ierr);
807827bd09bSSatish Balay   gs->min_node_pairs = i_start;
808827bd09bSSatish Balay 
809827bd09bSSatish Balay   i_start=j;
810827bd09bSSatish Balay   t1 = GL_ADD;
811f1ed62a8SBarry Smith   ierr = giop(&i_start,&offset,1,&t1);CHKERRQ(ierr);
812827bd09bSSatish Balay   gs->avg_node_pairs = i_start/num_nodes + 1;
813827bd09bSSatish Balay 
814827bd09bSSatish Balay   i_start=nprs;
815827bd09bSSatish Balay   t1 = GL_MAX;
816827bd09bSSatish Balay   giop(&i_start,&offset,1,&t1);
817827bd09bSSatish Balay   gs->max_pairs = i_start;
818827bd09bSSatish Balay 
819827bd09bSSatish Balay 
820827bd09bSSatish Balay   /* remap pairwise in tail of gsi_via_bit_mask() */
821827bd09bSSatish Balay   gs->msg_total = ivec_sum(gs->msg_sizes,nprs);
822a501084fSBarry Smith   gs->out = (PetscScalar *) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz);
823a501084fSBarry Smith   gs->in  = (PetscScalar *) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz);
824827bd09bSSatish Balay 
825827bd09bSSatish Balay   /* reset malloc pool */
826a501084fSBarry Smith   free((void*)p_mask);
827a501084fSBarry Smith   free((void*)tmp_proc_mask);
8283fdc5746SBarry Smith   PetscFunctionReturn(0);
829827bd09bSSatish Balay }
830827bd09bSSatish Balay 
831f1ed62a8SBarry Smith /* to do pruned tree just save ngh buf copy for each one and decode here!
832827bd09bSSatish Balay ******************************************************************************/
8330924e98cSBarry Smith static PetscErrorCode set_tree(gs_id *gs)
834827bd09bSSatish Balay {
83552f87cdaSBarry Smith   PetscInt i, j, n, nel;
83652f87cdaSBarry Smith   PetscInt *iptr_in, *iptr_out, *tree_elms, *elms;
837827bd09bSSatish Balay 
8383fdc5746SBarry Smith   PetscFunctionBegin;
839827bd09bSSatish Balay   /* local work ptrs */
840827bd09bSSatish Balay   elms = gs->elms;
841827bd09bSSatish Balay   nel     = gs->nel;
842827bd09bSSatish Balay 
843827bd09bSSatish Balay   /* how many via tree */
844827bd09bSSatish Balay   gs->tree_nel  = n = ntree;
845827bd09bSSatish Balay   gs->tree_elms = tree_elms = iptr_in = tree_buf;
846a501084fSBarry Smith   gs->tree_buf  = (PetscScalar *) malloc(sizeof(PetscScalar)*n*vec_sz);
847a501084fSBarry Smith   gs->tree_work = (PetscScalar *) malloc(sizeof(PetscScalar)*n*vec_sz);
848827bd09bSSatish Balay   j=gs->tree_map_sz;
84952f87cdaSBarry Smith   gs->tree_map_in = iptr_in  = (PetscInt*) malloc(sizeof(PetscInt)*(j+1));
85052f87cdaSBarry Smith   gs->tree_map_out = iptr_out = (PetscInt*) malloc(sizeof(PetscInt)*(j+1));
851827bd09bSSatish Balay 
852827bd09bSSatish Balay   /* search the longer of the two lists */
853827bd09bSSatish Balay   /* note ... could save this info in get_ngh_buf and save searches */
854827bd09bSSatish Balay   if (n<=nel)
855827bd09bSSatish Balay     {
856827bd09bSSatish Balay       /* bijective fct w/remap - search elm list */
857827bd09bSSatish Balay       for (i=0; i<n; i++)
858827bd09bSSatish Balay         {
859827bd09bSSatish Balay           if ((j=ivec_binary_search(*tree_elms++,elms,nel))>=0)
860827bd09bSSatish Balay             {*iptr_in++ = j; *iptr_out++ = i;}
861827bd09bSSatish Balay         }
862827bd09bSSatish Balay     }
863827bd09bSSatish Balay   else
864827bd09bSSatish Balay     {
865827bd09bSSatish Balay       for (i=0; i<nel; i++)
866827bd09bSSatish Balay         {
867827bd09bSSatish Balay           if ((j=ivec_binary_search(*elms++,tree_elms,n))>=0)
868827bd09bSSatish Balay             {*iptr_in++ = i; *iptr_out++ = j;}
869827bd09bSSatish Balay         }
870827bd09bSSatish Balay     }
871827bd09bSSatish Balay 
872827bd09bSSatish Balay   /* sentinel */
873827bd09bSSatish Balay   *iptr_in = *iptr_out = -1;
8743fdc5746SBarry Smith   PetscFunctionReturn(0);
875827bd09bSSatish Balay }
876827bd09bSSatish Balay 
877f1ed62a8SBarry Smith /******************************************************************************/
8780924e98cSBarry Smith static PetscErrorCode gs_gop_local_out( gs_id *gs,  PetscScalar *vals)
879827bd09bSSatish Balay {
88052f87cdaSBarry Smith   PetscInt *num, *map, **reduce;
881a501084fSBarry Smith   PetscScalar tmp;
882827bd09bSSatish Balay 
8833fdc5746SBarry Smith   PetscFunctionBegin;
884827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
885827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
886827bd09bSSatish Balay   while ((map = *reduce++))
887827bd09bSSatish Balay     {
888827bd09bSSatish Balay       /* wall */
889827bd09bSSatish Balay       if (*num == 2)
890827bd09bSSatish Balay         {
891827bd09bSSatish Balay           num ++;
892827bd09bSSatish Balay           vals[map[1]] = vals[map[0]];
893827bd09bSSatish Balay         }
894827bd09bSSatish Balay       /* corner shared by three elements */
895827bd09bSSatish Balay       else if (*num == 3)
896827bd09bSSatish Balay         {
897827bd09bSSatish Balay           num ++;
898827bd09bSSatish Balay           vals[map[2]] = vals[map[1]] = vals[map[0]];
899827bd09bSSatish Balay         }
900827bd09bSSatish Balay       /* corner shared by four elements */
901827bd09bSSatish Balay       else if (*num == 4)
902827bd09bSSatish Balay         {
903827bd09bSSatish Balay           num ++;
904827bd09bSSatish Balay           vals[map[3]] = vals[map[2]] = vals[map[1]] = vals[map[0]];
905827bd09bSSatish Balay         }
906827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
907827bd09bSSatish Balay       else
908827bd09bSSatish Balay         {
909827bd09bSSatish Balay           num++;
910827bd09bSSatish Balay           tmp = *(vals + *map++);
911827bd09bSSatish Balay           while (*map >= 0)
912827bd09bSSatish Balay             {*(vals + *map++) = tmp;}
913827bd09bSSatish Balay         }
914827bd09bSSatish Balay     }
9153fdc5746SBarry Smith   PetscFunctionReturn(0);
916827bd09bSSatish Balay }
917827bd09bSSatish Balay 
9187b1ae94cSBarry Smith /******************************************************************************/
9190924e98cSBarry Smith static PetscErrorCode gs_gop_local_plus( gs_id *gs,  PetscScalar *vals)
920827bd09bSSatish Balay {
92152f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
922a501084fSBarry Smith    PetscScalar tmp;
923827bd09bSSatish Balay 
9243fdc5746SBarry Smith   PetscFunctionBegin;
925827bd09bSSatish Balay   num    = gs->num_local_reduce;
926827bd09bSSatish Balay   reduce = gs->local_reduce;
927827bd09bSSatish Balay   while ((map = *reduce))
928827bd09bSSatish Balay     {
929827bd09bSSatish Balay       /* wall */
930827bd09bSSatish Balay       if (*num == 2)
931827bd09bSSatish Balay         {
932827bd09bSSatish Balay           num ++; reduce++;
933827bd09bSSatish Balay           vals[map[1]] = vals[map[0]] += vals[map[1]];
934827bd09bSSatish Balay         }
935827bd09bSSatish Balay       /* corner shared by three elements */
936827bd09bSSatish Balay       else if (*num == 3)
937827bd09bSSatish Balay         {
938827bd09bSSatish Balay           num ++; reduce++;
939827bd09bSSatish Balay           vals[map[2]]=vals[map[1]]=vals[map[0]]+=(vals[map[1]]+vals[map[2]]);
940827bd09bSSatish Balay         }
941827bd09bSSatish Balay       /* corner shared by four elements */
942827bd09bSSatish Balay       else if (*num == 4)
943827bd09bSSatish Balay         {
944827bd09bSSatish Balay           num ++; reduce++;
945827bd09bSSatish Balay           vals[map[1]]=vals[map[2]]=vals[map[3]]=vals[map[0]] +=
946827bd09bSSatish Balay                                  (vals[map[1]] + vals[map[2]] + vals[map[3]]);
947827bd09bSSatish Balay         }
948827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
949827bd09bSSatish Balay       else
950827bd09bSSatish Balay         {
951827bd09bSSatish Balay           num ++;
952827bd09bSSatish Balay           tmp = 0.0;
953827bd09bSSatish Balay           while (*map >= 0)
954827bd09bSSatish Balay             {tmp += *(vals + *map++);}
955827bd09bSSatish Balay 
956827bd09bSSatish Balay           map = *reduce++;
957827bd09bSSatish Balay           while (*map >= 0)
958827bd09bSSatish Balay             {*(vals + *map++) = tmp;}
959827bd09bSSatish Balay         }
960827bd09bSSatish Balay     }
9613fdc5746SBarry Smith   PetscFunctionReturn(0);
962827bd09bSSatish Balay }
963827bd09bSSatish Balay 
9647b1ae94cSBarry Smith /******************************************************************************/
9650924e98cSBarry Smith static PetscErrorCode gs_gop_local_in_plus( gs_id *gs,  PetscScalar *vals)
966827bd09bSSatish Balay {
96752f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
968a501084fSBarry Smith    PetscScalar *base;
969827bd09bSSatish Balay 
9703fdc5746SBarry Smith   PetscFunctionBegin;
971827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
972827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
973827bd09bSSatish Balay   while ((map = *reduce++))
974827bd09bSSatish Balay     {
975827bd09bSSatish Balay       /* wall */
976827bd09bSSatish Balay       if (*num == 2)
977827bd09bSSatish Balay         {
978827bd09bSSatish Balay           num ++;
979827bd09bSSatish Balay           vals[map[0]] += vals[map[1]];
980827bd09bSSatish Balay         }
981827bd09bSSatish Balay       /* corner shared by three elements */
982827bd09bSSatish Balay       else if (*num == 3)
983827bd09bSSatish Balay         {
984827bd09bSSatish Balay           num ++;
985827bd09bSSatish Balay           vals[map[0]] += (vals[map[1]] + vals[map[2]]);
986827bd09bSSatish Balay         }
987827bd09bSSatish Balay       /* corner shared by four elements */
988827bd09bSSatish Balay       else if (*num == 4)
989827bd09bSSatish Balay         {
990827bd09bSSatish Balay           num ++;
991827bd09bSSatish Balay           vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]);
992827bd09bSSatish Balay         }
993827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
994827bd09bSSatish Balay       else
995827bd09bSSatish Balay         {
996827bd09bSSatish Balay           num++;
997827bd09bSSatish Balay           base = vals + *map++;
998827bd09bSSatish Balay           while (*map >= 0)
999827bd09bSSatish Balay             {*base += *(vals + *map++);}
1000827bd09bSSatish Balay         }
1001827bd09bSSatish Balay     }
10023fdc5746SBarry Smith   PetscFunctionReturn(0);
1003827bd09bSSatish Balay }
1004827bd09bSSatish Balay 
10057b1ae94cSBarry Smith /******************************************************************************/
10060924e98cSBarry Smith PetscErrorCode gs_free( gs_id *gs)
1007827bd09bSSatish Balay {
100852f87cdaSBarry Smith    PetscInt i;
1009827bd09bSSatish Balay 
10103fdc5746SBarry Smith   PetscFunctionBegin;
1011a501084fSBarry Smith   if (gs->nghs) {free((void*) gs->nghs);}
1012a501084fSBarry Smith   if (gs->pw_nghs) {free((void*) gs->pw_nghs);}
1013827bd09bSSatish Balay 
1014827bd09bSSatish Balay   /* tree */
1015827bd09bSSatish Balay   if (gs->max_left_over)
1016827bd09bSSatish Balay     {
1017a501084fSBarry Smith       if (gs->tree_elms) {free((void*) gs->tree_elms);}
1018a501084fSBarry Smith       if (gs->tree_buf) {free((void*) gs->tree_buf);}
1019a501084fSBarry Smith       if (gs->tree_work) {free((void*) gs->tree_work);}
1020a501084fSBarry Smith       if (gs->tree_map_in) {free((void*) gs->tree_map_in);}
1021a501084fSBarry Smith       if (gs->tree_map_out) {free((void*) gs->tree_map_out);}
1022827bd09bSSatish Balay     }
1023827bd09bSSatish Balay 
1024827bd09bSSatish Balay   /* pairwise info */
1025827bd09bSSatish Balay   if (gs->num_pairs)
1026827bd09bSSatish Balay     {
1027827bd09bSSatish Balay       /* should be NULL already */
1028a501084fSBarry Smith       if (gs->ngh_buf) {free((void*) gs->ngh_buf);}
1029a501084fSBarry Smith       if (gs->elms) {free((void*) gs->elms);}
1030a501084fSBarry Smith       if (gs->local_elms) {free((void*) gs->local_elms);}
1031a501084fSBarry Smith       if (gs->companion) {free((void*) gs->companion);}
1032827bd09bSSatish Balay 
1033827bd09bSSatish Balay       /* only set if pairwise */
1034a501084fSBarry Smith       if (gs->vals) {free((void*) gs->vals);}
1035a501084fSBarry Smith       if (gs->in) {free((void*) gs->in);}
1036a501084fSBarry Smith       if (gs->out) {free((void*) gs->out);}
1037a501084fSBarry Smith       if (gs->msg_ids_in) {free((void*) gs->msg_ids_in);}
1038a501084fSBarry Smith       if (gs->msg_ids_out) {free((void*) gs->msg_ids_out);}
1039a501084fSBarry Smith       if (gs->pw_vals) {free((void*) gs->pw_vals);}
1040a501084fSBarry Smith       if (gs->pw_elm_list) {free((void*) gs->pw_elm_list);}
1041827bd09bSSatish Balay       if (gs->node_list)
1042827bd09bSSatish Balay         {
1043827bd09bSSatish Balay           for (i=0;i<gs->num_pairs;i++)
1044a501084fSBarry Smith             {if (gs->node_list[i]) {free((void*) gs->node_list[i]);}}
1045a501084fSBarry Smith           free((void*) gs->node_list);
1046827bd09bSSatish Balay         }
1047a501084fSBarry Smith       if (gs->msg_sizes) {free((void*) gs->msg_sizes);}
1048a501084fSBarry Smith       if (gs->pair_list) {free((void*) gs->pair_list);}
1049827bd09bSSatish Balay     }
1050827bd09bSSatish Balay 
1051827bd09bSSatish Balay   /* local info */
1052827bd09bSSatish Balay   if (gs->num_local_total>=0)
1053827bd09bSSatish Balay     {
1054827bd09bSSatish Balay       for (i=0;i<gs->num_local_total+1;i++)
1055827bd09bSSatish Balay         /*      for (i=0;i<gs->num_local_total;i++) */
1056827bd09bSSatish Balay         {
1057827bd09bSSatish Balay           if (gs->num_gop_local_reduce[i])
1058a501084fSBarry Smith             {free((void*) gs->gop_local_reduce[i]);}
1059827bd09bSSatish Balay         }
1060827bd09bSSatish Balay     }
1061827bd09bSSatish Balay 
1062827bd09bSSatish Balay   /* if intersection tree/pairwise and local isn't empty */
1063a501084fSBarry Smith   if (gs->gop_local_reduce) {free((void*) gs->gop_local_reduce);}
1064a501084fSBarry Smith   if (gs->num_gop_local_reduce) {free((void*) gs->num_gop_local_reduce);}
1065827bd09bSSatish Balay 
1066a501084fSBarry Smith   free((void*) gs);
10673fdc5746SBarry Smith   PetscFunctionReturn(0);
1068827bd09bSSatish Balay }
1069827bd09bSSatish Balay 
10707b1ae94cSBarry Smith /******************************************************************************/
107152f87cdaSBarry Smith PetscErrorCode gs_gop_vec( gs_id *gs,  PetscScalar *vals,  const char *op,  PetscInt step)
1072827bd09bSSatish Balay {
1073d1528f56SBarry Smith   PetscErrorCode ierr;
1074d1528f56SBarry Smith 
10753fdc5746SBarry Smith   PetscFunctionBegin;
1076827bd09bSSatish Balay   switch (*op) {
1077827bd09bSSatish Balay   case '+':
1078827bd09bSSatish Balay     gs_gop_vec_plus(gs,vals,step);
1079827bd09bSSatish Balay     break;
1080827bd09bSSatish Balay   default:
1081f1ed62a8SBarry Smith     ierr = PetscInfo1(0,"gs_gop_vec() :: %c is not a valid op",op[0]);CHKERRQ(ierr);
1082f1ed62a8SBarry Smith     ierr = PetscInfo(0,"gs_gop_vec() :: default :: plus");CHKERRQ(ierr);
1083827bd09bSSatish Balay     gs_gop_vec_plus(gs,vals,step);
1084827bd09bSSatish Balay     break;
1085827bd09bSSatish Balay   }
10863fdc5746SBarry Smith   PetscFunctionReturn(0);
1087827bd09bSSatish Balay }
1088827bd09bSSatish Balay 
10897b1ae94cSBarry Smith /******************************************************************************/
109052f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_plus( gs_id *gs,  PetscScalar *vals,  PetscInt step)
1091827bd09bSSatish Balay {
10923fdc5746SBarry Smith   PetscFunctionBegin;
1093388eb383SBarry Smith   if (!gs) {SETERRQ(PETSC_ERR_PLIB,"gs_gop_vec() passed NULL gs handle!!!");}
1094827bd09bSSatish Balay 
1095827bd09bSSatish Balay   /* local only operations!!! */
1096827bd09bSSatish Balay   if (gs->num_local)
1097827bd09bSSatish Balay     {gs_gop_vec_local_plus(gs,vals,step);}
1098827bd09bSSatish Balay 
1099827bd09bSSatish Balay   /* if intersection tree/pairwise and local isn't empty */
1100827bd09bSSatish Balay   if (gs->num_local_gop)
1101827bd09bSSatish Balay     {
1102827bd09bSSatish Balay       gs_gop_vec_local_in_plus(gs,vals,step);
1103827bd09bSSatish Balay 
1104827bd09bSSatish Balay       /* pairwise */
1105827bd09bSSatish Balay       if (gs->num_pairs)
1106827bd09bSSatish Balay         {gs_gop_vec_pairwise_plus(gs,vals,step);}
1107827bd09bSSatish Balay 
1108827bd09bSSatish Balay       /* tree */
1109827bd09bSSatish Balay       else if (gs->max_left_over)
1110827bd09bSSatish Balay         {gs_gop_vec_tree_plus(gs,vals,step);}
1111827bd09bSSatish Balay 
1112827bd09bSSatish Balay       gs_gop_vec_local_out(gs,vals,step);
1113827bd09bSSatish Balay     }
1114827bd09bSSatish Balay   /* if intersection tree/pairwise and local is empty */
1115827bd09bSSatish Balay   else
1116827bd09bSSatish Balay     {
1117827bd09bSSatish Balay       /* pairwise */
1118827bd09bSSatish Balay       if (gs->num_pairs)
1119827bd09bSSatish Balay         {gs_gop_vec_pairwise_plus(gs,vals,step);}
1120827bd09bSSatish Balay 
1121827bd09bSSatish Balay       /* tree */
1122827bd09bSSatish Balay       else if (gs->max_left_over)
1123827bd09bSSatish Balay         {gs_gop_vec_tree_plus(gs,vals,step);}
1124827bd09bSSatish Balay     }
11253fdc5746SBarry Smith   PetscFunctionReturn(0);
1126827bd09bSSatish Balay }
1127827bd09bSSatish Balay 
11287b1ae94cSBarry Smith /******************************************************************************/
112952f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_local_plus( gs_id *gs,  PetscScalar *vals, PetscInt step)
1130827bd09bSSatish Balay {
113152f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
1132a501084fSBarry Smith    PetscScalar *base;
1133827bd09bSSatish Balay 
11343fdc5746SBarry Smith   PetscFunctionBegin;
1135827bd09bSSatish Balay   num    = gs->num_local_reduce;
1136827bd09bSSatish Balay   reduce = gs->local_reduce;
1137827bd09bSSatish Balay   while ((map = *reduce))
1138827bd09bSSatish Balay     {
1139827bd09bSSatish Balay       base = vals + map[0] * step;
1140827bd09bSSatish Balay 
1141827bd09bSSatish Balay       /* wall */
1142827bd09bSSatish Balay       if (*num == 2)
1143827bd09bSSatish Balay         {
1144827bd09bSSatish Balay           num++; reduce++;
1145827bd09bSSatish Balay           rvec_add (base,vals+map[1]*step,step);
1146827bd09bSSatish Balay           rvec_copy(vals+map[1]*step,base,step);
1147827bd09bSSatish Balay         }
1148827bd09bSSatish Balay       /* corner shared by three elements */
1149827bd09bSSatish Balay       else if (*num == 3)
1150827bd09bSSatish Balay         {
1151827bd09bSSatish Balay           num++; reduce++;
1152827bd09bSSatish Balay           rvec_add (base,vals+map[1]*step,step);
1153827bd09bSSatish Balay           rvec_add (base,vals+map[2]*step,step);
1154827bd09bSSatish Balay           rvec_copy(vals+map[2]*step,base,step);
1155827bd09bSSatish Balay           rvec_copy(vals+map[1]*step,base,step);
1156827bd09bSSatish Balay         }
1157827bd09bSSatish Balay       /* corner shared by four elements */
1158827bd09bSSatish Balay       else if (*num == 4)
1159827bd09bSSatish Balay         {
1160827bd09bSSatish Balay           num++; reduce++;
1161827bd09bSSatish Balay           rvec_add (base,vals+map[1]*step,step);
1162827bd09bSSatish Balay           rvec_add (base,vals+map[2]*step,step);
1163827bd09bSSatish Balay           rvec_add (base,vals+map[3]*step,step);
1164827bd09bSSatish Balay           rvec_copy(vals+map[3]*step,base,step);
1165827bd09bSSatish Balay           rvec_copy(vals+map[2]*step,base,step);
1166827bd09bSSatish Balay           rvec_copy(vals+map[1]*step,base,step);
1167827bd09bSSatish Balay         }
1168827bd09bSSatish Balay       /* general case ... odd geoms ... 3D */
1169827bd09bSSatish Balay       else
1170827bd09bSSatish Balay         {
1171827bd09bSSatish Balay           num++;
1172827bd09bSSatish Balay           while (*++map >= 0)
1173827bd09bSSatish Balay             {rvec_add (base,vals+*map*step,step);}
1174827bd09bSSatish Balay 
1175827bd09bSSatish Balay           map = *reduce;
1176827bd09bSSatish Balay           while (*++map >= 0)
1177827bd09bSSatish Balay             {rvec_copy(vals+*map*step,base,step);}
1178827bd09bSSatish Balay 
1179827bd09bSSatish Balay           reduce++;
1180827bd09bSSatish Balay         }
1181827bd09bSSatish Balay     }
11823fdc5746SBarry Smith   PetscFunctionReturn(0);
1183827bd09bSSatish Balay }
1184827bd09bSSatish Balay 
11857b1ae94cSBarry Smith /******************************************************************************/
118652f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_local_in_plus( gs_id *gs,  PetscScalar *vals, PetscInt step)
1187827bd09bSSatish Balay {
118852f87cdaSBarry Smith    PetscInt  *num, *map, **reduce;
1189a501084fSBarry Smith    PetscScalar *base;
11903fdc5746SBarry Smith   PetscFunctionBegin;
1191827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
1192827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
1193827bd09bSSatish Balay   while ((map = *reduce++))
1194827bd09bSSatish Balay     {
1195827bd09bSSatish Balay       base = vals + map[0] * step;
1196827bd09bSSatish Balay 
1197827bd09bSSatish Balay       /* wall */
1198827bd09bSSatish Balay       if (*num == 2)
1199827bd09bSSatish Balay         {
1200827bd09bSSatish Balay           num ++;
1201827bd09bSSatish Balay           rvec_add(base,vals+map[1]*step,step);
1202827bd09bSSatish Balay         }
1203827bd09bSSatish Balay       /* corner shared by three elements */
1204827bd09bSSatish Balay       else if (*num == 3)
1205827bd09bSSatish Balay         {
1206827bd09bSSatish Balay           num ++;
1207827bd09bSSatish Balay           rvec_add(base,vals+map[1]*step,step);
1208827bd09bSSatish Balay           rvec_add(base,vals+map[2]*step,step);
1209827bd09bSSatish Balay         }
1210827bd09bSSatish Balay       /* corner shared by four elements */
1211827bd09bSSatish Balay       else if (*num == 4)
1212827bd09bSSatish Balay         {
1213827bd09bSSatish Balay           num ++;
1214827bd09bSSatish Balay           rvec_add(base,vals+map[1]*step,step);
1215827bd09bSSatish Balay           rvec_add(base,vals+map[2]*step,step);
1216827bd09bSSatish Balay           rvec_add(base,vals+map[3]*step,step);
1217827bd09bSSatish Balay         }
1218827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
1219827bd09bSSatish Balay       else
1220827bd09bSSatish Balay         {
1221827bd09bSSatish Balay           num++;
1222827bd09bSSatish Balay           while (*++map >= 0)
1223827bd09bSSatish Balay             {rvec_add(base,vals+*map*step,step);}
1224827bd09bSSatish Balay         }
1225827bd09bSSatish Balay     }
12263fdc5746SBarry Smith   PetscFunctionReturn(0);
1227827bd09bSSatish Balay }
1228827bd09bSSatish Balay 
12297b1ae94cSBarry Smith /******************************************************************************/
123052f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_local_out( gs_id *gs,  PetscScalar *vals, PetscInt step)
1231827bd09bSSatish Balay {
123252f87cdaSBarry Smith    PetscInt *num, *map, **reduce;
1233a501084fSBarry Smith    PetscScalar *base;
1234827bd09bSSatish Balay 
12353fdc5746SBarry Smith   PetscFunctionBegin;
1236827bd09bSSatish Balay   num    = gs->num_gop_local_reduce;
1237827bd09bSSatish Balay   reduce = gs->gop_local_reduce;
1238827bd09bSSatish Balay   while ((map = *reduce++))
1239827bd09bSSatish Balay     {
1240827bd09bSSatish Balay       base = vals + map[0] * step;
1241827bd09bSSatish Balay 
1242827bd09bSSatish Balay       /* wall */
1243827bd09bSSatish Balay       if (*num == 2)
1244827bd09bSSatish Balay         {
1245827bd09bSSatish Balay           num ++;
1246827bd09bSSatish Balay           rvec_copy(vals+map[1]*step,base,step);
1247827bd09bSSatish Balay         }
1248827bd09bSSatish Balay       /* corner shared by three elements */
1249827bd09bSSatish Balay       else if (*num == 3)
1250827bd09bSSatish Balay         {
1251827bd09bSSatish Balay           num ++;
1252827bd09bSSatish Balay           rvec_copy(vals+map[1]*step,base,step);
1253827bd09bSSatish Balay           rvec_copy(vals+map[2]*step,base,step);
1254827bd09bSSatish Balay         }
1255827bd09bSSatish Balay       /* corner shared by four elements */
1256827bd09bSSatish Balay       else if (*num == 4)
1257827bd09bSSatish Balay         {
1258827bd09bSSatish Balay           num ++;
1259827bd09bSSatish Balay           rvec_copy(vals+map[1]*step,base,step);
1260827bd09bSSatish Balay           rvec_copy(vals+map[2]*step,base,step);
1261827bd09bSSatish Balay           rvec_copy(vals+map[3]*step,base,step);
1262827bd09bSSatish Balay         }
1263827bd09bSSatish Balay       /* general case ... odd geoms ... 3D*/
1264827bd09bSSatish Balay       else
1265827bd09bSSatish Balay         {
1266827bd09bSSatish Balay           num++;
1267827bd09bSSatish Balay           while (*++map >= 0)
1268827bd09bSSatish Balay             {rvec_copy(vals+*map*step,base,step);}
1269827bd09bSSatish Balay         }
1270827bd09bSSatish Balay     }
12713fdc5746SBarry Smith   PetscFunctionReturn(0);
1272827bd09bSSatish Balay }
1273827bd09bSSatish Balay 
12747b1ae94cSBarry Smith /******************************************************************************/
127552f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_pairwise_plus( gs_id *gs,  PetscScalar *in_vals, PetscInt step)
1276827bd09bSSatish Balay {
1277a501084fSBarry Smith   PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2;
127852f87cdaSBarry Smith   PetscInt *iptr, *msg_list, *msg_size, **msg_nodes;
127952f87cdaSBarry Smith   PetscInt *pw, *list, *size, **nodes;
1280827bd09bSSatish Balay   MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out;
1281827bd09bSSatish Balay   MPI_Status status;
12820805154bSBarry Smith   PetscBLASInt i1 = 1,dstep;
12833fdc5746SBarry Smith   PetscErrorCode ierr;
1284827bd09bSSatish Balay 
12853fdc5746SBarry Smith   PetscFunctionBegin;
1286a501084fSBarry Smith   /* strip and load s */
1287827bd09bSSatish Balay   msg_list =list         = gs->pair_list;
1288827bd09bSSatish Balay   msg_size =size         = gs->msg_sizes;
1289827bd09bSSatish Balay   msg_nodes=nodes        = gs->node_list;
1290827bd09bSSatish Balay   iptr=pw                = gs->pw_elm_list;
1291827bd09bSSatish Balay   dptr1=dptr3            = gs->pw_vals;
1292827bd09bSSatish Balay   msg_ids_in  = ids_in   = gs->msg_ids_in;
1293827bd09bSSatish Balay   msg_ids_out = ids_out  = gs->msg_ids_out;
1294827bd09bSSatish Balay   dptr2                  = gs->out;
1295827bd09bSSatish Balay   in1=in2                = gs->in;
1296827bd09bSSatish Balay 
1297827bd09bSSatish Balay   /* post the receives */
1298827bd09bSSatish Balay   /*  msg_nodes=nodes; */
1299827bd09bSSatish Balay   do
1300827bd09bSSatish Balay     {
1301827bd09bSSatish Balay       /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the
1302827bd09bSSatish Balay          second one *list and do list++ afterwards */
13033fdc5746SBarry Smith       ierr = MPI_Irecv(in1, *size *step, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list++, gs->gs_comm, msg_ids_in++);CHKERRQ(ierr);
1304827bd09bSSatish Balay       in1 += *size++ *step;
1305827bd09bSSatish Balay     }
1306827bd09bSSatish Balay   while (*++msg_nodes);
1307827bd09bSSatish Balay   msg_nodes=nodes;
1308827bd09bSSatish Balay 
1309827bd09bSSatish Balay   /* load gs values into in out gs buffers */
1310827bd09bSSatish Balay   while (*iptr >= 0)
1311827bd09bSSatish Balay     {
1312827bd09bSSatish Balay       rvec_copy(dptr3,in_vals + *iptr*step,step);
1313827bd09bSSatish Balay       dptr3+=step;
1314827bd09bSSatish Balay       iptr++;
1315827bd09bSSatish Balay     }
1316827bd09bSSatish Balay 
1317827bd09bSSatish Balay   /* load out buffers and post the sends */
1318827bd09bSSatish Balay   while ((iptr = *msg_nodes++))
1319827bd09bSSatish Balay     {
1320827bd09bSSatish Balay       dptr3 = dptr2;
1321827bd09bSSatish Balay       while (*iptr >= 0)
1322827bd09bSSatish Balay         {
1323827bd09bSSatish Balay           rvec_copy(dptr2,dptr1 + *iptr*step,step);
1324827bd09bSSatish Balay           dptr2+=step;
1325827bd09bSSatish Balay           iptr++;
1326827bd09bSSatish Balay         }
13273fdc5746SBarry Smith       ierr = MPI_Isend(dptr3, *msg_size++ *step, MPIU_SCALAR, *msg_list++, MSGTAG1+my_id, gs->gs_comm, msg_ids_out++);CHKERRQ(ierr);
1328827bd09bSSatish Balay     }
1329827bd09bSSatish Balay 
1330827bd09bSSatish Balay   /* tree */
1331827bd09bSSatish Balay   if (gs->max_left_over)
1332827bd09bSSatish Balay     {gs_gop_vec_tree_plus(gs,in_vals,step);}
1333827bd09bSSatish Balay 
1334827bd09bSSatish Balay   /* process the received data */
1335827bd09bSSatish Balay   msg_nodes=nodes;
1336a501084fSBarry Smith   while ((iptr = *nodes++)){
1337a501084fSBarry Smith     PetscScalar d1 = 1.0;
1338827bd09bSSatish Balay       /* Should I check the return value of MPI_Wait() or status? */
1339827bd09bSSatish Balay       /* Can this loop be replaced by a call to MPI_Waitall()? */
13403fdc5746SBarry Smith       ierr = MPI_Wait(ids_in++, &status);CHKERRQ(ierr);
1341a501084fSBarry Smith       while (*iptr >= 0) {
13420805154bSBarry Smith 	dstep = PetscBLASIntCast(step);
13434a0de3f6SBarry Smith         BLASaxpy_(&dstep,&d1,in2,&i1,dptr1 + *iptr*step,&i1);
1344827bd09bSSatish Balay 	in2+=step;
1345827bd09bSSatish Balay 	iptr++;
1346827bd09bSSatish Balay       }
1347827bd09bSSatish Balay   }
1348827bd09bSSatish Balay 
1349827bd09bSSatish Balay   /* replace vals */
1350827bd09bSSatish Balay   while (*pw >= 0)
1351827bd09bSSatish Balay     {
1352827bd09bSSatish Balay       rvec_copy(in_vals + *pw*step,dptr1,step);
1353827bd09bSSatish Balay       dptr1+=step;
1354827bd09bSSatish Balay       pw++;
1355827bd09bSSatish Balay     }
1356827bd09bSSatish Balay 
1357827bd09bSSatish Balay   /* clear isend message handles */
1358827bd09bSSatish Balay   /* This changed for clarity though it could be the same */
1359827bd09bSSatish Balay   while (*msg_nodes++)
1360827bd09bSSatish Balay     /* Should I check the return value of MPI_Wait() or status? */
1361827bd09bSSatish Balay     /* Can this loop be replaced by a call to MPI_Waitall()? */
13623fdc5746SBarry Smith     {ierr = MPI_Wait(ids_out++, &status);CHKERRQ(ierr);}
1363827bd09bSSatish Balay 
13643fdc5746SBarry Smith   PetscFunctionReturn(0);
1365827bd09bSSatish Balay }
1366827bd09bSSatish Balay 
13677b1ae94cSBarry Smith /******************************************************************************/
136852f87cdaSBarry Smith static PetscErrorCode gs_gop_vec_tree_plus( gs_id *gs,  PetscScalar *vals,  PetscInt step)
1369827bd09bSSatish Balay {
137052f87cdaSBarry Smith   PetscInt size, *in, *out;
1371a501084fSBarry Smith   PetscScalar *buf, *work;
137252f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
1373a501084fSBarry Smith   PetscBLASInt i1 = 1;
1374827bd09bSSatish Balay 
13753fdc5746SBarry Smith   PetscFunctionBegin;
1376827bd09bSSatish Balay   /* copy over to local variables */
1377827bd09bSSatish Balay   in   = gs->tree_map_in;
1378827bd09bSSatish Balay   out  = gs->tree_map_out;
1379827bd09bSSatish Balay   buf  = gs->tree_buf;
1380827bd09bSSatish Balay   work = gs->tree_work;
1381827bd09bSSatish Balay   size = gs->tree_nel*step;
1382827bd09bSSatish Balay 
1383827bd09bSSatish Balay   /* zero out collection buffer */
1384827bd09bSSatish Balay   rvec_zero(buf,size);
1385827bd09bSSatish Balay 
1386827bd09bSSatish Balay 
1387827bd09bSSatish Balay   /* copy over my contributions */
1388827bd09bSSatish Balay   while (*in >= 0)
1389827bd09bSSatish Balay     {
13900805154bSBarry Smith       PetscBLASInt dstep = PetscBLASIntCast(step);
13916e4f4d19SBarry Smith       BLAScopy_(&dstep,vals + *in++*step,&i1,buf + *out++*step,&i1);
1392827bd09bSSatish Balay     }
1393827bd09bSSatish Balay 
1394827bd09bSSatish Balay   /* perform fan in/out on full buffer */
1395827bd09bSSatish Balay   /* must change grop to handle the blas */
1396827bd09bSSatish Balay   grop(buf,work,size,op);
1397827bd09bSSatish Balay 
1398827bd09bSSatish Balay   /* reset */
1399827bd09bSSatish Balay   in   = gs->tree_map_in;
1400827bd09bSSatish Balay   out  = gs->tree_map_out;
1401827bd09bSSatish Balay 
1402827bd09bSSatish Balay   /* get the portion of the results I need */
1403827bd09bSSatish Balay   while (*in >= 0)
1404827bd09bSSatish Balay     {
14050805154bSBarry Smith       PetscBLASInt dstep = PetscBLASIntCast(step);
14066e4f4d19SBarry Smith       BLAScopy_(&dstep,buf + *out++*step,&i1,vals + *in++*step,&i1);
1407827bd09bSSatish Balay     }
14083fdc5746SBarry Smith   PetscFunctionReturn(0);
1409827bd09bSSatish Balay }
1410827bd09bSSatish Balay 
14117b1ae94cSBarry Smith /******************************************************************************/
141252f87cdaSBarry Smith PetscErrorCode gs_gop_hc( gs_id *gs,  PetscScalar *vals,  const char *op,  PetscInt dim)
1413827bd09bSSatish Balay {
1414d1528f56SBarry Smith   PetscErrorCode ierr;
1415d1528f56SBarry Smith 
14163fdc5746SBarry Smith   PetscFunctionBegin;
1417827bd09bSSatish Balay   switch (*op) {
1418827bd09bSSatish Balay   case '+':
1419827bd09bSSatish Balay     gs_gop_plus_hc(gs,vals,dim);
1420827bd09bSSatish Balay     break;
1421827bd09bSSatish Balay   default:
1422f1ed62a8SBarry Smith     ierr = PetscInfo1(0,"gs_gop_hc() :: %c is not a valid op",op[0]);CHKERRQ(ierr);
1423f1ed62a8SBarry Smith     ierr = PetscInfo(0,"gs_gop_hc() :: default :: plus\n");CHKERRQ(ierr);
1424827bd09bSSatish Balay     gs_gop_plus_hc(gs,vals,dim);
1425827bd09bSSatish Balay     break;
1426827bd09bSSatish Balay   }
14273fdc5746SBarry Smith   PetscFunctionReturn(0);
1428827bd09bSSatish Balay }
1429827bd09bSSatish Balay 
14307b1ae94cSBarry Smith /******************************************************************************/
143152f87cdaSBarry Smith static PetscErrorCode gs_gop_plus_hc( gs_id *gs,  PetscScalar *vals, PetscInt dim)
1432827bd09bSSatish Balay {
14333fdc5746SBarry Smith   PetscFunctionBegin;
1434827bd09bSSatish Balay   /* if there's nothing to do return */
1435827bd09bSSatish Balay   if (dim<=0)
14363fdc5746SBarry Smith     {  PetscFunctionReturn(0);}
1437827bd09bSSatish Balay 
1438827bd09bSSatish Balay   /* can't do more dimensions then exist */
143939945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
1440827bd09bSSatish Balay 
1441827bd09bSSatish Balay   /* local only operations!!! */
1442827bd09bSSatish Balay   if (gs->num_local)
1443827bd09bSSatish Balay     {gs_gop_local_plus(gs,vals);}
1444827bd09bSSatish Balay 
1445827bd09bSSatish Balay   /* if intersection tree/pairwise and local isn't empty */
1446827bd09bSSatish Balay   if (gs->num_local_gop)
1447827bd09bSSatish Balay     {
1448827bd09bSSatish Balay       gs_gop_local_in_plus(gs,vals);
1449827bd09bSSatish Balay 
1450827bd09bSSatish Balay       /* pairwise will do tree inside ... */
1451827bd09bSSatish Balay       if (gs->num_pairs)
1452827bd09bSSatish Balay         {gs_gop_pairwise_plus_hc(gs,vals,dim);}
1453827bd09bSSatish Balay 
1454827bd09bSSatish Balay       /* tree only */
1455827bd09bSSatish Balay       else if (gs->max_left_over)
1456827bd09bSSatish Balay         {gs_gop_tree_plus_hc(gs,vals,dim);}
1457827bd09bSSatish Balay 
1458827bd09bSSatish Balay       gs_gop_local_out(gs,vals);
1459827bd09bSSatish Balay     }
1460827bd09bSSatish Balay   /* if intersection tree/pairwise and local is empty */
1461827bd09bSSatish Balay   else
1462827bd09bSSatish Balay     {
1463827bd09bSSatish Balay       /* pairwise will do tree inside */
1464827bd09bSSatish Balay       if (gs->num_pairs)
1465827bd09bSSatish Balay         {gs_gop_pairwise_plus_hc(gs,vals,dim);}
1466827bd09bSSatish Balay 
1467827bd09bSSatish Balay       /* tree */
1468827bd09bSSatish Balay       else if (gs->max_left_over)
1469827bd09bSSatish Balay         {gs_gop_tree_plus_hc(gs,vals,dim);}
1470827bd09bSSatish Balay     }
14713fdc5746SBarry Smith   PetscFunctionReturn(0);
1472827bd09bSSatish Balay }
1473827bd09bSSatish Balay 
14747b1ae94cSBarry Smith /******************************************************************************/
147552f87cdaSBarry Smith static PetscErrorCode gs_gop_pairwise_plus_hc( gs_id *gs,  PetscScalar *in_vals, PetscInt dim)
1476827bd09bSSatish Balay {
1477a501084fSBarry Smith    PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2;
147852f87cdaSBarry Smith    PetscInt *iptr, *msg_list, *msg_size, **msg_nodes;
147952f87cdaSBarry Smith    PetscInt *pw, *list, *size, **nodes;
1480827bd09bSSatish Balay   MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out;
1481827bd09bSSatish Balay   MPI_Status status;
148252f87cdaSBarry Smith   PetscInt i, mask=1;
14833fdc5746SBarry Smith   PetscErrorCode ierr;
1484827bd09bSSatish Balay 
14853fdc5746SBarry Smith   PetscFunctionBegin;
1486827bd09bSSatish Balay   for (i=1; i<dim; i++)
1487827bd09bSSatish Balay     {mask<<=1; mask++;}
1488827bd09bSSatish Balay 
1489827bd09bSSatish Balay 
1490a501084fSBarry Smith   /* strip and load s */
1491827bd09bSSatish Balay   msg_list =list         = gs->pair_list;
1492827bd09bSSatish Balay   msg_size =size         = gs->msg_sizes;
1493827bd09bSSatish Balay   msg_nodes=nodes        = gs->node_list;
1494827bd09bSSatish Balay   iptr=pw                = gs->pw_elm_list;
1495827bd09bSSatish Balay   dptr1=dptr3            = gs->pw_vals;
1496827bd09bSSatish Balay   msg_ids_in  = ids_in   = gs->msg_ids_in;
1497827bd09bSSatish Balay   msg_ids_out = ids_out  = gs->msg_ids_out;
1498827bd09bSSatish Balay   dptr2                  = gs->out;
1499827bd09bSSatish Balay   in1=in2                = gs->in;
1500827bd09bSSatish Balay 
1501827bd09bSSatish Balay   /* post the receives */
1502827bd09bSSatish Balay   /*  msg_nodes=nodes; */
1503827bd09bSSatish Balay   do
1504827bd09bSSatish Balay     {
1505827bd09bSSatish Balay       /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the
1506827bd09bSSatish Balay          second one *list and do list++ afterwards */
1507827bd09bSSatish Balay       if ((my_id|mask)==(*list|mask))
1508827bd09bSSatish Balay         {
15093fdc5746SBarry Smith           ierr = MPI_Irecv(in1, *size, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list++, gs->gs_comm, msg_ids_in++);CHKERRQ(ierr);
1510827bd09bSSatish Balay           in1 += *size++;
1511827bd09bSSatish Balay         }
1512827bd09bSSatish Balay       else
1513827bd09bSSatish Balay         {list++; size++;}
1514827bd09bSSatish Balay     }
1515827bd09bSSatish Balay   while (*++msg_nodes);
1516827bd09bSSatish Balay 
1517827bd09bSSatish Balay   /* load gs values into in out gs buffers */
1518827bd09bSSatish Balay   while (*iptr >= 0)
1519827bd09bSSatish Balay     {*dptr3++ = *(in_vals + *iptr++);}
1520827bd09bSSatish Balay 
1521827bd09bSSatish Balay   /* load out buffers and post the sends */
1522827bd09bSSatish Balay   msg_nodes=nodes;
1523827bd09bSSatish Balay   list = msg_list;
1524827bd09bSSatish Balay   while ((iptr = *msg_nodes++))
1525827bd09bSSatish Balay     {
1526827bd09bSSatish Balay       if ((my_id|mask)==(*list|mask))
1527827bd09bSSatish Balay         {
1528827bd09bSSatish Balay           dptr3 = dptr2;
1529827bd09bSSatish Balay           while (*iptr >= 0)
1530827bd09bSSatish Balay             {*dptr2++ = *(dptr1 + *iptr++);}
1531827bd09bSSatish Balay           /* CHECK PERSISTENT COMMS MODE FOR ALL THIS STUFF */
1532827bd09bSSatish Balay           /* is msg_ids_out++ correct? */
15333fdc5746SBarry Smith           ierr = MPI_Isend(dptr3, *msg_size++, MPIU_SCALAR, *list++, MSGTAG1+my_id, gs->gs_comm, msg_ids_out++);CHKERRQ(ierr);
1534827bd09bSSatish Balay         }
1535827bd09bSSatish Balay       else
1536827bd09bSSatish Balay         {list++; msg_size++;}
1537827bd09bSSatish Balay     }
1538827bd09bSSatish Balay 
1539827bd09bSSatish Balay   /* do the tree while we're waiting */
1540827bd09bSSatish Balay   if (gs->max_left_over)
1541827bd09bSSatish Balay     {gs_gop_tree_plus_hc(gs,in_vals,dim);}
1542827bd09bSSatish Balay 
1543827bd09bSSatish Balay   /* process the received data */
1544827bd09bSSatish Balay   msg_nodes=nodes;
1545827bd09bSSatish Balay   list = msg_list;
1546827bd09bSSatish Balay   while ((iptr = *nodes++))
1547827bd09bSSatish Balay     {
1548827bd09bSSatish Balay       if ((my_id|mask)==(*list|mask))
1549827bd09bSSatish Balay         {
1550827bd09bSSatish Balay           /* Should I check the return value of MPI_Wait() or status? */
1551827bd09bSSatish Balay           /* Can this loop be replaced by a call to MPI_Waitall()? */
15523fdc5746SBarry Smith           ierr = MPI_Wait(ids_in++, &status);CHKERRQ(ierr);
1553827bd09bSSatish Balay           while (*iptr >= 0)
1554827bd09bSSatish Balay             {*(dptr1 + *iptr++) += *in2++;}
1555827bd09bSSatish Balay         }
1556827bd09bSSatish Balay       list++;
1557827bd09bSSatish Balay     }
1558827bd09bSSatish Balay 
1559827bd09bSSatish Balay   /* replace vals */
1560827bd09bSSatish Balay   while (*pw >= 0)
1561827bd09bSSatish Balay     {*(in_vals + *pw++) = *dptr1++;}
1562827bd09bSSatish Balay 
1563827bd09bSSatish Balay   /* clear isend message handles */
1564827bd09bSSatish Balay   /* This changed for clarity though it could be the same */
1565827bd09bSSatish Balay   while (*msg_nodes++)
1566827bd09bSSatish Balay     {
1567827bd09bSSatish Balay       if ((my_id|mask)==(*msg_list|mask))
1568827bd09bSSatish Balay         {
1569827bd09bSSatish Balay           /* Should I check the return value of MPI_Wait() or status? */
1570827bd09bSSatish Balay           /* Can this loop be replaced by a call to MPI_Waitall()? */
15713fdc5746SBarry Smith           ierr = MPI_Wait(ids_out++, &status);CHKERRQ(ierr);
1572827bd09bSSatish Balay         }
1573827bd09bSSatish Balay       msg_list++;
1574827bd09bSSatish Balay     }
1575827bd09bSSatish Balay 
15763fdc5746SBarry Smith   PetscFunctionReturn(0);
1577827bd09bSSatish Balay }
1578827bd09bSSatish Balay 
15797b1ae94cSBarry Smith /******************************************************************************/
158052f87cdaSBarry Smith static PetscErrorCode gs_gop_tree_plus_hc(gs_id *gs, PetscScalar *vals, PetscInt dim)
1581827bd09bSSatish Balay {
158252f87cdaSBarry Smith   PetscInt size;
158352f87cdaSBarry Smith   PetscInt *in, *out;
1584a501084fSBarry Smith   PetscScalar *buf, *work;
158552f87cdaSBarry Smith   PetscInt op[] = {GL_ADD,0};
1586827bd09bSSatish Balay 
15873fdc5746SBarry Smith   PetscFunctionBegin;
1588827bd09bSSatish Balay   in   = gs->tree_map_in;
1589827bd09bSSatish Balay   out  = gs->tree_map_out;
1590827bd09bSSatish Balay   buf  = gs->tree_buf;
1591827bd09bSSatish Balay   work = gs->tree_work;
1592827bd09bSSatish Balay   size = gs->tree_nel;
1593827bd09bSSatish Balay 
1594827bd09bSSatish Balay   rvec_zero(buf,size);
1595827bd09bSSatish Balay 
1596827bd09bSSatish Balay   while (*in >= 0)
1597827bd09bSSatish Balay     {*(buf + *out++) = *(vals + *in++);}
1598827bd09bSSatish Balay 
1599827bd09bSSatish Balay   in   = gs->tree_map_in;
1600827bd09bSSatish Balay   out  = gs->tree_map_out;
1601827bd09bSSatish Balay 
1602827bd09bSSatish Balay   grop_hc(buf,work,size,op,dim);
1603827bd09bSSatish Balay 
1604827bd09bSSatish Balay   while (*in >= 0)
1605827bd09bSSatish Balay     {*(vals + *in++) = *(buf + *out++);}
16063fdc5746SBarry Smith   PetscFunctionReturn(0);
1607827bd09bSSatish Balay }
1608827bd09bSSatish Balay 
1609827bd09bSSatish Balay 
1610827bd09bSSatish Balay 
1611