xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision dba47a550923b04c7c4ebbb735eb62a1b3e4e9ae)
1*dba47a55SKris Buschelman #define PETSCKSP_DLL
2827bd09bSSatish Balay 
3827bd09bSSatish Balay /***********************************comm.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 11.21.97
16827bd09bSSatish Balay ***********************************comm.c*************************************/
17827bd09bSSatish Balay 
18827bd09bSSatish Balay /***********************************comm.c*************************************
19827bd09bSSatish Balay File Description:
20827bd09bSSatish Balay -----------------
21827bd09bSSatish Balay 
22827bd09bSSatish Balay ***********************************comm.c*************************************/
237758a8cdSBarry Smith #include "src/ksp/pc/impls/tfs/tfs.h"
24827bd09bSSatish Balay 
25827bd09bSSatish Balay 
26827bd09bSSatish Balay /* global program control variables - explicitly exported */
27827bd09bSSatish Balay int my_id            = 0;
28827bd09bSSatish Balay int num_nodes        = 1;
29827bd09bSSatish Balay int floor_num_nodes  = 0;
30827bd09bSSatish Balay int i_log2_num_nodes = 0;
31827bd09bSSatish Balay 
32827bd09bSSatish Balay /* global program control variables */
33827bd09bSSatish Balay static int p_init = 0;
34827bd09bSSatish Balay static int modfl_num_nodes;
35827bd09bSSatish Balay static int edge_not_pow_2;
36827bd09bSSatish Balay 
37a501084fSBarry Smith static unsigned int edge_node[sizeof(PetscInt)*32];
38827bd09bSSatish Balay 
39827bd09bSSatish Balay /***********************************comm.c*************************************
40827bd09bSSatish Balay Function: giop()
41827bd09bSSatish Balay 
42827bd09bSSatish Balay Input :
43827bd09bSSatish Balay Output:
44827bd09bSSatish Balay Return:
45827bd09bSSatish Balay Description:
46827bd09bSSatish Balay ***********************************comm.c*************************************/
47827bd09bSSatish Balay void
48827bd09bSSatish Balay comm_init (void)
49827bd09bSSatish Balay {
50827bd09bSSatish Balay 
51827bd09bSSatish Balay   if (p_init++) return;
52827bd09bSSatish Balay 
53827bd09bSSatish Balay #if PETSC_SIZEOF_INT != 4
54827bd09bSSatish Balay   error_msg_warning("Int != 4 Bytes!");
55827bd09bSSatish Balay #endif
56827bd09bSSatish Balay 
57827bd09bSSatish Balay 
58827bd09bSSatish Balay   MPI_Comm_size(MPI_COMM_WORLD,&num_nodes);
59827bd09bSSatish Balay   MPI_Comm_rank(MPI_COMM_WORLD,&my_id);
60827bd09bSSatish Balay 
61827bd09bSSatish Balay   if (num_nodes> (INT_MAX >> 1))
62827bd09bSSatish Balay   {error_msg_fatal("Can't have more then MAX_INT/2 nodes!!!");}
63827bd09bSSatish Balay 
64a501084fSBarry Smith   ivec_zero((int*)edge_node,sizeof(PetscInt)*32);
65827bd09bSSatish Balay 
66827bd09bSSatish Balay   floor_num_nodes = 1;
67827bd09bSSatish Balay   i_log2_num_nodes = modfl_num_nodes = 0;
68827bd09bSSatish Balay   while (floor_num_nodes <= num_nodes)
69827bd09bSSatish Balay     {
70827bd09bSSatish Balay       edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes;
71827bd09bSSatish Balay       floor_num_nodes <<= 1;
72827bd09bSSatish Balay       i_log2_num_nodes++;
73827bd09bSSatish Balay     }
74827bd09bSSatish Balay 
75827bd09bSSatish Balay   i_log2_num_nodes--;
76827bd09bSSatish Balay   floor_num_nodes >>= 1;
77827bd09bSSatish Balay   modfl_num_nodes = (num_nodes - floor_num_nodes);
78827bd09bSSatish Balay 
79827bd09bSSatish Balay   if ((my_id > 0) && (my_id <= modfl_num_nodes))
80827bd09bSSatish Balay     {edge_not_pow_2=((my_id|floor_num_nodes)-1);}
81827bd09bSSatish Balay   else if (my_id >= floor_num_nodes)
82827bd09bSSatish Balay     {edge_not_pow_2=((my_id^floor_num_nodes)+1);
83827bd09bSSatish Balay     }
84827bd09bSSatish Balay   else
85827bd09bSSatish Balay     {edge_not_pow_2 = 0;}
86827bd09bSSatish Balay 
87827bd09bSSatish Balay }
88827bd09bSSatish Balay 
89827bd09bSSatish Balay 
90827bd09bSSatish Balay 
91827bd09bSSatish Balay /***********************************comm.c*************************************
92827bd09bSSatish Balay Function: giop()
93827bd09bSSatish Balay 
94827bd09bSSatish Balay Input :
95827bd09bSSatish Balay Output:
96827bd09bSSatish Balay Return:
97827bd09bSSatish Balay Description: fan-in/out version
98827bd09bSSatish Balay ***********************************comm.c*************************************/
99827bd09bSSatish Balay void
100827bd09bSSatish Balay giop(int *vals, int *work, int n, int *oprs)
101827bd09bSSatish Balay {
102a501084fSBarry Smith    int mask, edge;
103827bd09bSSatish Balay   int type, dest;
104827bd09bSSatish Balay   vfp fp;
105827bd09bSSatish Balay   MPI_Status  status;
106827bd09bSSatish Balay 
107827bd09bSSatish Balay 
108827bd09bSSatish Balay #ifdef SAFE
109827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
110827bd09bSSatish Balay   if (!vals||!work||!oprs)
11177431f27SBarry Smith     {error_msg_fatal("giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
112827bd09bSSatish Balay 
113827bd09bSSatish Balay   /* non-uniform should have at least two entries */
114827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
115827bd09bSSatish Balay     {error_msg_fatal("giop() :: non_uniform and n=0,1?");}
116827bd09bSSatish Balay 
117827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
118827bd09bSSatish Balay   if (!p_init)
119827bd09bSSatish Balay     {comm_init();}
120827bd09bSSatish Balay #endif
121827bd09bSSatish Balay 
122827bd09bSSatish Balay   /* if there's nothing to do return */
123827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
124827bd09bSSatish Balay     {
125827bd09bSSatish Balay       return;
126827bd09bSSatish Balay     }
127827bd09bSSatish Balay 
128827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
129827bd09bSSatish Balay   if (n<0)
13077431f27SBarry Smith     {error_msg_fatal("giop() :: n=%D<0?",n);}
131827bd09bSSatish Balay 
132827bd09bSSatish Balay   /* advance to list of n operations for custom */
133827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
134827bd09bSSatish Balay     {oprs++;}
135827bd09bSSatish Balay 
136827bd09bSSatish Balay   /* major league hack */
137d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))) {
138827bd09bSSatish Balay     error_msg_warning("giop() :: hope you passed in a rbfp!\n");
139827bd09bSSatish Balay     fp = (vfp) oprs;
140827bd09bSSatish Balay   }
141827bd09bSSatish Balay 
142827bd09bSSatish Balay   /* all msgs will be of the same length */
143827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
144827bd09bSSatish Balay   if (edge_not_pow_2)
145827bd09bSSatish Balay     {
146827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
147827bd09bSSatish Balay 	{MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);}
148827bd09bSSatish Balay       else
149827bd09bSSatish Balay 	{
150827bd09bSSatish Balay 	  MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
151827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
152827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
153827bd09bSSatish Balay 	}
154827bd09bSSatish Balay     }
155827bd09bSSatish Balay 
156827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
157827bd09bSSatish Balay   if (my_id<floor_num_nodes)
158827bd09bSSatish Balay     {
159827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
160827bd09bSSatish Balay 	{
161827bd09bSSatish Balay 	  dest = my_id^mask;
162827bd09bSSatish Balay 	  if (my_id > dest)
163827bd09bSSatish Balay 	    {MPI_Send(vals,n,MPI_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
164827bd09bSSatish Balay 	  else
165827bd09bSSatish Balay 	    {
166827bd09bSSatish Balay 	      MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG2+dest,
167827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
168827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
169827bd09bSSatish Balay 	    }
170827bd09bSSatish Balay 	}
171827bd09bSSatish Balay 
172827bd09bSSatish Balay       mask=floor_num_nodes>>1;
173827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
174827bd09bSSatish Balay 	{
175827bd09bSSatish Balay 	  if (my_id%mask)
176827bd09bSSatish Balay 	    {continue;}
177827bd09bSSatish Balay 
178827bd09bSSatish Balay 	  dest = my_id^mask;
179827bd09bSSatish Balay 	  if (my_id < dest)
180827bd09bSSatish Balay 	    {MPI_Send(vals,n,MPI_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
181827bd09bSSatish Balay 	  else
182827bd09bSSatish Balay 	    {
183827bd09bSSatish Balay 	      MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG4+dest,
184827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
185827bd09bSSatish Balay 	    }
186827bd09bSSatish Balay 	}
187827bd09bSSatish Balay     }
188827bd09bSSatish Balay 
189827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
190827bd09bSSatish Balay   if (edge_not_pow_2)
191827bd09bSSatish Balay     {
192827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
193827bd09bSSatish Balay 	{
194827bd09bSSatish Balay 	  MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
195827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
196827bd09bSSatish Balay 	}
197827bd09bSSatish Balay       else
198827bd09bSSatish Balay 	{MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);}
199827bd09bSSatish Balay     }
200827bd09bSSatish Balay }
201827bd09bSSatish Balay 
202827bd09bSSatish Balay /***********************************comm.c*************************************
203827bd09bSSatish Balay Function: grop()
204827bd09bSSatish Balay 
205827bd09bSSatish Balay Input :
206827bd09bSSatish Balay Output:
207827bd09bSSatish Balay Return:
208827bd09bSSatish Balay Description: fan-in/out version
209827bd09bSSatish Balay ***********************************comm.c*************************************/
210827bd09bSSatish Balay void
211a501084fSBarry Smith grop(PetscScalar *vals, PetscScalar *work, int n, int *oprs)
212827bd09bSSatish Balay {
213a501084fSBarry Smith    int mask, edge;
214827bd09bSSatish Balay   int type, dest;
215827bd09bSSatish Balay   vfp fp;
216827bd09bSSatish Balay   MPI_Status  status;
217827bd09bSSatish Balay 
218827bd09bSSatish Balay 
219827bd09bSSatish Balay #ifdef SAFE
220827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
221827bd09bSSatish Balay   if (!vals||!work||!oprs)
22277431f27SBarry Smith     {error_msg_fatal("grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
223827bd09bSSatish Balay 
224827bd09bSSatish Balay   /* non-uniform should have at least two entries */
225827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
226827bd09bSSatish Balay     {error_msg_fatal("grop() :: non_uniform and n=0,1?");}
227827bd09bSSatish Balay 
228827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
229827bd09bSSatish Balay   if (!p_init)
230827bd09bSSatish Balay     {comm_init();}
231827bd09bSSatish Balay #endif
232827bd09bSSatish Balay 
233827bd09bSSatish Balay   /* if there's nothing to do return */
234827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
235827bd09bSSatish Balay     {return;}
236827bd09bSSatish Balay 
237827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
238827bd09bSSatish Balay   if (n<0)
23977431f27SBarry Smith     {error_msg_fatal("gdop() :: n=%D<0?",n);}
240827bd09bSSatish Balay 
241827bd09bSSatish Balay   /* advance to list of n operations for custom */
242827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
243827bd09bSSatish Balay     {oprs++;}
244827bd09bSSatish Balay 
245d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
246827bd09bSSatish Balay     error_msg_warning("grop() :: hope you passed in a rbfp!\n");
247827bd09bSSatish Balay     fp = (vfp) oprs;
248827bd09bSSatish Balay   }
249827bd09bSSatish Balay 
250827bd09bSSatish Balay   /* all msgs will be of the same length */
251827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
252827bd09bSSatish Balay   if (edge_not_pow_2)
253827bd09bSSatish Balay     {
254827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
255a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+my_id,
256827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
257827bd09bSSatish Balay       else
258827bd09bSSatish Balay 	{
259a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
260827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
261827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
262827bd09bSSatish Balay 	}
263827bd09bSSatish Balay     }
264827bd09bSSatish Balay 
265827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
266827bd09bSSatish Balay   if (my_id<floor_num_nodes)
267827bd09bSSatish Balay     {
268827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
269827bd09bSSatish Balay 	{
270827bd09bSSatish Balay 	  dest = my_id^mask;
271827bd09bSSatish Balay 	  if (my_id > dest)
272a501084fSBarry Smith 	    {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
273827bd09bSSatish Balay 	  else
274827bd09bSSatish Balay 	    {
275a501084fSBarry Smith 	      MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,
276827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
277827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
278827bd09bSSatish Balay 	    }
279827bd09bSSatish Balay 	}
280827bd09bSSatish Balay 
281827bd09bSSatish Balay       mask=floor_num_nodes>>1;
282827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
283827bd09bSSatish Balay 	{
284827bd09bSSatish Balay 	  if (my_id%mask)
285827bd09bSSatish Balay 	    {continue;}
286827bd09bSSatish Balay 
287827bd09bSSatish Balay 	  dest = my_id^mask;
288827bd09bSSatish Balay 	  if (my_id < dest)
289a501084fSBarry Smith 	    {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
290827bd09bSSatish Balay 	  else
291827bd09bSSatish Balay 	    {
292a501084fSBarry Smith 	      MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,
293827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
294827bd09bSSatish Balay 	    }
295827bd09bSSatish Balay 	}
296827bd09bSSatish Balay     }
297827bd09bSSatish Balay 
298827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
299827bd09bSSatish Balay   if (edge_not_pow_2)
300827bd09bSSatish Balay     {
301827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
302827bd09bSSatish Balay 	{
303a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
304827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
305827bd09bSSatish Balay 	}
306827bd09bSSatish Balay       else
307a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+my_id,
308827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
309827bd09bSSatish Balay     }
310827bd09bSSatish Balay }
311827bd09bSSatish Balay 
312827bd09bSSatish Balay 
313827bd09bSSatish Balay /***********************************comm.c*************************************
314827bd09bSSatish Balay Function: grop()
315827bd09bSSatish Balay 
316827bd09bSSatish Balay Input :
317827bd09bSSatish Balay Output:
318827bd09bSSatish Balay Return:
319827bd09bSSatish Balay Description: fan-in/out version
320827bd09bSSatish Balay 
321827bd09bSSatish Balay note good only for num_nodes=2^k!!!
322827bd09bSSatish Balay 
323827bd09bSSatish Balay ***********************************comm.c*************************************/
324827bd09bSSatish Balay void
325a501084fSBarry Smith grop_hc(PetscScalar *vals, PetscScalar *work, int n, int *oprs, int dim)
326827bd09bSSatish Balay {
327a501084fSBarry Smith    int mask, edge;
328827bd09bSSatish Balay   int type, dest;
329827bd09bSSatish Balay   vfp fp;
330827bd09bSSatish Balay   MPI_Status  status;
331827bd09bSSatish Balay 
332827bd09bSSatish Balay #ifdef SAFE
333827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
334827bd09bSSatish Balay   if (!vals||!work||!oprs)
33577431f27SBarry Smith     {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
336827bd09bSSatish Balay 
337827bd09bSSatish Balay   /* non-uniform should have at least two entries */
338827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
339827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: non_uniform and n=0,1?");}
340827bd09bSSatish Balay 
341827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
342827bd09bSSatish Balay   if (!p_init)
343827bd09bSSatish Balay     {comm_init();}
344827bd09bSSatish Balay #endif
345827bd09bSSatish Balay 
346827bd09bSSatish Balay   /* if there's nothing to do return */
347827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
348827bd09bSSatish Balay     {return;}
349827bd09bSSatish Balay 
350827bd09bSSatish Balay   /* the error msg says it all!!! */
351827bd09bSSatish Balay   if (modfl_num_nodes)
352827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");}
353827bd09bSSatish Balay 
354827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
355827bd09bSSatish Balay   if (n<0)
35677431f27SBarry Smith     {error_msg_fatal("grop_hc() :: n=%D<0?",n);}
357827bd09bSSatish Balay 
358827bd09bSSatish Balay   /* can't do more dimensions then exist */
35939945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
360827bd09bSSatish Balay 
361827bd09bSSatish Balay   /* advance to list of n operations for custom */
362827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
363827bd09bSSatish Balay     {oprs++;}
364827bd09bSSatish Balay 
365d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
366827bd09bSSatish Balay     error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n");
367827bd09bSSatish Balay     fp = (vfp) oprs;
368827bd09bSSatish Balay   }
369827bd09bSSatish Balay 
370827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
371827bd09bSSatish Balay     {
372827bd09bSSatish Balay       dest = my_id^mask;
373827bd09bSSatish Balay       if (my_id > dest)
374a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
375827bd09bSSatish Balay       else
376827bd09bSSatish Balay 	{
377a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,
378827bd09bSSatish Balay 		   &status);
379827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
380827bd09bSSatish Balay 	}
381827bd09bSSatish Balay     }
382827bd09bSSatish Balay 
383827bd09bSSatish Balay   if (edge==dim)
384827bd09bSSatish Balay     {mask>>=1;}
385827bd09bSSatish Balay   else
386827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
387827bd09bSSatish Balay 
388827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
389827bd09bSSatish Balay     {
390827bd09bSSatish Balay       if (my_id%mask)
391827bd09bSSatish Balay 	{continue;}
392827bd09bSSatish Balay 
393827bd09bSSatish Balay       dest = my_id^mask;
394827bd09bSSatish Balay       if (my_id < dest)
395a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
396827bd09bSSatish Balay       else
397827bd09bSSatish Balay 	{
398a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,
399827bd09bSSatish Balay 		   &status);
400827bd09bSSatish Balay 	}
401827bd09bSSatish Balay     }
402827bd09bSSatish Balay }
403827bd09bSSatish Balay 
404827bd09bSSatish Balay 
405827bd09bSSatish Balay /***********************************comm.c*************************************
406827bd09bSSatish Balay Function: gop()
407827bd09bSSatish Balay 
408827bd09bSSatish Balay Input :
409827bd09bSSatish Balay Output:
410827bd09bSSatish Balay Return:
411827bd09bSSatish Balay Description: fan-in/out version
412827bd09bSSatish Balay ***********************************comm.c*************************************/
413a501084fSBarry Smith void gfop(void *vals, void *work, int n, vbfp fp, MPI_Datatype dt, int comm_type)
414827bd09bSSatish Balay {
415a501084fSBarry Smith    int mask, edge;
416827bd09bSSatish Balay   int dest;
417827bd09bSSatish Balay   MPI_Status  status;
418827bd09bSSatish Balay   MPI_Op op;
419827bd09bSSatish Balay 
420827bd09bSSatish Balay #ifdef SAFE
421827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
422827bd09bSSatish Balay   if (!p_init)
423827bd09bSSatish Balay     {comm_init();}
424827bd09bSSatish Balay 
425827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
426827bd09bSSatish Balay   if (!vals||!work||!fp)
42777431f27SBarry Smith     {error_msg_fatal("gop() :: v=%D, w=%D, f=%D",vals,work,fp);}
428827bd09bSSatish Balay #endif
429827bd09bSSatish Balay 
430827bd09bSSatish Balay   /* if there's nothing to do return */
431827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
432827bd09bSSatish Balay     {return;}
433827bd09bSSatish Balay 
434827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
435827bd09bSSatish Balay   if (n<0)
43677431f27SBarry Smith     {error_msg_fatal("gop() :: n=%D<0?",n);}
437827bd09bSSatish Balay 
438827bd09bSSatish Balay   if (comm_type==MPI)
439827bd09bSSatish Balay     {
440827bd09bSSatish Balay       MPI_Op_create(fp,TRUE,&op);
441827bd09bSSatish Balay       MPI_Allreduce (vals, work, n, dt, op, MPI_COMM_WORLD);
442827bd09bSSatish Balay       MPI_Op_free(&op);
443827bd09bSSatish Balay       return;
444827bd09bSSatish Balay     }
445827bd09bSSatish Balay 
446827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
447827bd09bSSatish Balay   if (edge_not_pow_2)
448827bd09bSSatish Balay     {
449827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
450827bd09bSSatish Balay 	{MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG0+my_id,
451827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
452827bd09bSSatish Balay       else
453827bd09bSSatish Balay 	{
454827bd09bSSatish Balay 	  MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
455827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
456827bd09bSSatish Balay 	  (*fp)(vals,work,&n,&dt);
457827bd09bSSatish Balay 	}
458827bd09bSSatish Balay     }
459827bd09bSSatish Balay 
460827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
461827bd09bSSatish Balay   if (my_id<floor_num_nodes)
462827bd09bSSatish Balay     {
463827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
464827bd09bSSatish Balay 	{
465827bd09bSSatish Balay 	  dest = my_id^mask;
466827bd09bSSatish Balay 	  if (my_id > dest)
467827bd09bSSatish Balay 	    {MPI_Send(vals,n,dt,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
468827bd09bSSatish Balay 	  else
469827bd09bSSatish Balay 	    {
470827bd09bSSatish Balay 	      MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG2+dest,
471827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
472827bd09bSSatish Balay 	      (*fp)(vals, work, &n, &dt);
473827bd09bSSatish Balay 	    }
474827bd09bSSatish Balay 	}
475827bd09bSSatish Balay 
476827bd09bSSatish Balay       mask=floor_num_nodes>>1;
477827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
478827bd09bSSatish Balay 	{
479827bd09bSSatish Balay 	  if (my_id%mask)
480827bd09bSSatish Balay 	    {continue;}
481827bd09bSSatish Balay 
482827bd09bSSatish Balay 	  dest = my_id^mask;
483827bd09bSSatish Balay 	  if (my_id < dest)
484827bd09bSSatish Balay 	    {MPI_Send(vals,n,dt,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
485827bd09bSSatish Balay 	  else
486827bd09bSSatish Balay 	    {
487827bd09bSSatish Balay 	      MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG4+dest,
488827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
489827bd09bSSatish Balay 	    }
490827bd09bSSatish Balay 	}
491827bd09bSSatish Balay     }
492827bd09bSSatish Balay 
493827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
494827bd09bSSatish Balay   if (edge_not_pow_2)
495827bd09bSSatish Balay     {
496827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
497827bd09bSSatish Balay 	{
498827bd09bSSatish Balay 	  MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
499827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
500827bd09bSSatish Balay 	}
501827bd09bSSatish Balay       else
502827bd09bSSatish Balay 	{MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG5+my_id,
503827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
504827bd09bSSatish Balay     }
505827bd09bSSatish Balay }
506827bd09bSSatish Balay 
507827bd09bSSatish Balay 
508827bd09bSSatish Balay 
509827bd09bSSatish Balay 
510827bd09bSSatish Balay 
511827bd09bSSatish Balay 
512827bd09bSSatish Balay /******************************************************************************
513827bd09bSSatish Balay Function: giop()
514827bd09bSSatish Balay 
515827bd09bSSatish Balay Input :
516827bd09bSSatish Balay Output:
517827bd09bSSatish Balay Return:
518827bd09bSSatish Balay Description:
519827bd09bSSatish Balay 
520827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii
521827bd09bSSatish Balay 
522827bd09bSSatish Balay ******************************************************************************/
523827bd09bSSatish Balay void
524a501084fSBarry Smith ssgl_radd( PetscScalar *vals,  PetscScalar *work,  int level,
525a501084fSBarry Smith 	   int *segs)
526827bd09bSSatish Balay {
527a501084fSBarry Smith    int edge, type, dest, mask;
528a501084fSBarry Smith    int stage_n;
529827bd09bSSatish Balay   MPI_Status  status;
530827bd09bSSatish Balay 
531827bd09bSSatish Balay #ifdef SAFE
532827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
533827bd09bSSatish Balay   if (!p_init)
534827bd09bSSatish Balay     {comm_init();}
535827bd09bSSatish Balay #endif
536827bd09bSSatish Balay 
537827bd09bSSatish Balay 
538827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
539827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
540827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
541827bd09bSSatish Balay     {
542827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
543827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
544827bd09bSSatish Balay 	{
545827bd09bSSatish Balay 	  dest = edge_node[edge];
546827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
547827bd09bSSatish Balay 	  if (my_id>dest)
548a501084fSBarry Smith           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type,
549827bd09bSSatish Balay                       MPI_COMM_WORLD);}
550827bd09bSSatish Balay 	  else
551827bd09bSSatish Balay 	    {
552827bd09bSSatish Balay 	      type =  type - my_id + dest;
553a501084fSBarry Smith               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,
554827bd09bSSatish Balay                        MPI_COMM_WORLD,&status);
555827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
556827bd09bSSatish Balay 	    }
557827bd09bSSatish Balay 	}
558827bd09bSSatish Balay       mask <<= 1;
559827bd09bSSatish Balay     }
560827bd09bSSatish Balay   mask>>=1;
561827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
562827bd09bSSatish Balay     {
563827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
564827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
565827bd09bSSatish Balay 	{
566827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
567827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
568827bd09bSSatish Balay 	  if (my_id<dest)
569a501084fSBarry Smith             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,
570827bd09bSSatish Balay                       MPI_COMM_WORLD);}
571827bd09bSSatish Balay 	  else
572827bd09bSSatish Balay 	    {
573827bd09bSSatish Balay 	      type =  type - my_id + dest;
574a501084fSBarry Smith               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,
575827bd09bSSatish Balay                        MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
576827bd09bSSatish Balay 	    }
577827bd09bSSatish Balay 	}
578827bd09bSSatish Balay       mask >>= 1;
579827bd09bSSatish Balay     }
580827bd09bSSatish Balay }
581827bd09bSSatish Balay 
582827bd09bSSatish Balay 
583827bd09bSSatish Balay 
584827bd09bSSatish Balay /***********************************comm.c*************************************
585827bd09bSSatish Balay Function: grop_hc_vvl()
586827bd09bSSatish Balay 
587827bd09bSSatish Balay Input :
588827bd09bSSatish Balay Output:
589827bd09bSSatish Balay Return:
590827bd09bSSatish Balay Description: fan-in/out version
591827bd09bSSatish Balay 
592827bd09bSSatish Balay note good only for num_nodes=2^k!!!
593827bd09bSSatish Balay 
594827bd09bSSatish Balay ***********************************comm.c*************************************/
595827bd09bSSatish Balay void
596a501084fSBarry Smith grop_hc_vvl(PetscScalar *vals, PetscScalar *work, int *segs, int *oprs, int dim)
597827bd09bSSatish Balay {
598a501084fSBarry Smith    int mask, edge, n;
599827bd09bSSatish Balay   int type, dest;
600827bd09bSSatish Balay   vfp fp;
601827bd09bSSatish Balay   MPI_Status  status;
602827bd09bSSatish Balay 
603827bd09bSSatish Balay   error_msg_fatal("grop_hc_vvl() :: is not working!\n");
604827bd09bSSatish Balay 
605827bd09bSSatish Balay #ifdef SAFE
606827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
607827bd09bSSatish Balay   if (!vals||!work||!oprs||!segs)
60877431f27SBarry Smith     {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
609827bd09bSSatish Balay 
610827bd09bSSatish Balay   /* non-uniform should have at least two entries */
611827bd09bSSatish Balay 
612827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
613827bd09bSSatish Balay   if (!p_init)
614827bd09bSSatish Balay     {comm_init();}
615827bd09bSSatish Balay #endif
616827bd09bSSatish Balay 
617827bd09bSSatish Balay   /* if there's nothing to do return */
618827bd09bSSatish Balay   if ((num_nodes<2)||(dim<=0))
619827bd09bSSatish Balay     {return;}
620827bd09bSSatish Balay 
621827bd09bSSatish Balay   /* the error msg says it all!!! */
622827bd09bSSatish Balay   if (modfl_num_nodes)
623827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");}
624827bd09bSSatish Balay 
625827bd09bSSatish Balay   /* can't do more dimensions then exist */
62639945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
627827bd09bSSatish Balay 
628827bd09bSSatish Balay   /* advance to list of n operations for custom */
629827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
630827bd09bSSatish Balay     {oprs++;}
631827bd09bSSatish Balay 
632d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))){
633827bd09bSSatish Balay     error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n");
634827bd09bSSatish Balay     fp = (vfp) oprs;
635827bd09bSSatish Balay   }
636827bd09bSSatish Balay 
637827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
638827bd09bSSatish Balay     {
639827bd09bSSatish Balay       n = segs[dim]-segs[edge];
640827bd09bSSatish Balay       dest = my_id^mask;
641827bd09bSSatish Balay       if (my_id > dest)
642a501084fSBarry Smith 	{MPI_Send(vals+segs[edge],n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
643827bd09bSSatish Balay       else
644827bd09bSSatish Balay 	{
645a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,
646827bd09bSSatish Balay 		   MPI_COMM_WORLD, &status);
647827bd09bSSatish Balay 	  (*fp)(vals+segs[edge], work, n, oprs);
648827bd09bSSatish Balay 	}
649827bd09bSSatish Balay     }
650827bd09bSSatish Balay 
651827bd09bSSatish Balay   if (edge==dim)
652827bd09bSSatish Balay     {mask>>=1;}
653827bd09bSSatish Balay   else
654827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
655827bd09bSSatish Balay 
656827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
657827bd09bSSatish Balay     {
658827bd09bSSatish Balay       if (my_id%mask)
659827bd09bSSatish Balay 	{continue;}
660827bd09bSSatish Balay 
661827bd09bSSatish Balay       n = (segs[dim]-segs[dim-1-edge]);
662827bd09bSSatish Balay 
663827bd09bSSatish Balay       dest = my_id^mask;
664827bd09bSSatish Balay       if (my_id < dest)
665a501084fSBarry Smith 	{MPI_Send(vals+segs[dim-1-edge],n,MPIU_SCALAR,dest,MSGTAG4+my_id,
666827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
667827bd09bSSatish Balay       else
668827bd09bSSatish Balay 	{
669a501084fSBarry Smith 	  MPI_Recv(vals+segs[dim-1-edge],n,MPIU_SCALAR,MPI_ANY_SOURCE,
670827bd09bSSatish Balay 		   MSGTAG4+dest,MPI_COMM_WORLD, &status);
671827bd09bSSatish Balay 	}
672827bd09bSSatish Balay     }
673827bd09bSSatish Balay }
674827bd09bSSatish Balay 
675827bd09bSSatish Balay /******************************************************************************
676827bd09bSSatish Balay Function: giop()
677827bd09bSSatish Balay 
678827bd09bSSatish Balay Input :
679827bd09bSSatish Balay Output:
680827bd09bSSatish Balay Return:
681827bd09bSSatish Balay Description:
682827bd09bSSatish Balay 
683827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii
684827bd09bSSatish Balay 
685827bd09bSSatish Balay ******************************************************************************/
686a501084fSBarry Smith void new_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  int level, int *segs)
687827bd09bSSatish Balay {
688a501084fSBarry Smith    int edge, type, dest, mask;
689a501084fSBarry Smith    int stage_n;
690827bd09bSSatish Balay   MPI_Status  status;
691827bd09bSSatish Balay 
692827bd09bSSatish Balay #ifdef SAFE
693827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
694827bd09bSSatish Balay   if (!p_init)
695827bd09bSSatish Balay     {comm_init();}
696827bd09bSSatish Balay #endif
697827bd09bSSatish Balay 
698827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
699827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
700827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
701827bd09bSSatish Balay     {
702827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
703827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
704827bd09bSSatish Balay 	{
705827bd09bSSatish Balay 	  dest = edge_node[edge];
706827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
707827bd09bSSatish Balay 	  if (my_id>dest)
708a501084fSBarry Smith           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type,
709c700b191SBarry Smith                       MPI_COMM_WORLD);}
710827bd09bSSatish Balay 	  else
711827bd09bSSatish Balay 	    {
712827bd09bSSatish Balay 	      type =  type - my_id + dest;
713a501084fSBarry Smith               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,
714c700b191SBarry Smith                        MPI_COMM_WORLD,&status);
715827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
716827bd09bSSatish Balay 	    }
717827bd09bSSatish Balay 	}
718827bd09bSSatish Balay       mask <<= 1;
719827bd09bSSatish Balay     }
720827bd09bSSatish Balay   mask>>=1;
721827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
722827bd09bSSatish Balay     {
723827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
724827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
725827bd09bSSatish Balay 	{
726827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
727827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
728827bd09bSSatish Balay 	  if (my_id<dest)
729a501084fSBarry Smith             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,
730c700b191SBarry Smith                       MPI_COMM_WORLD);}
731827bd09bSSatish Balay 	  else
732827bd09bSSatish Balay 	    {
733827bd09bSSatish Balay 	      type =  type - my_id + dest;
734a501084fSBarry Smith               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,
735c700b191SBarry Smith                        MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
736827bd09bSSatish Balay 	    }
737827bd09bSSatish Balay 	}
738827bd09bSSatish Balay       mask >>= 1;
739827bd09bSSatish Balay     }
740827bd09bSSatish Balay }
741827bd09bSSatish Balay 
742827bd09bSSatish Balay 
743827bd09bSSatish Balay 
744827bd09bSSatish Balay /***********************************comm.c*************************************
745827bd09bSSatish Balay Function: giop()
746827bd09bSSatish Balay 
747827bd09bSSatish Balay Input :
748827bd09bSSatish Balay Output:
749827bd09bSSatish Balay Return:
750827bd09bSSatish Balay Description: fan-in/out version
751827bd09bSSatish Balay 
752827bd09bSSatish Balay note good only for num_nodes=2^k!!!
753827bd09bSSatish Balay 
754827bd09bSSatish Balay ***********************************comm.c*************************************/
755c700b191SBarry Smith void giop_hc(int *vals, int *work, int n, int *oprs, int dim)
756827bd09bSSatish Balay {
757a501084fSBarry Smith    int mask, edge;
758827bd09bSSatish Balay   int type, dest;
759827bd09bSSatish Balay   vfp fp;
760827bd09bSSatish Balay   MPI_Status  status;
761827bd09bSSatish Balay 
762827bd09bSSatish Balay #ifdef SAFE
763827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
764827bd09bSSatish Balay   if (!vals||!work||!oprs)
76577431f27SBarry Smith     {error_msg_fatal("giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
766827bd09bSSatish Balay 
767827bd09bSSatish Balay   /* non-uniform should have at least two entries */
768827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
769827bd09bSSatish Balay     {error_msg_fatal("giop_hc() :: non_uniform and n=0,1?");}
770827bd09bSSatish Balay 
771827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
772827bd09bSSatish Balay   if (!p_init)
773827bd09bSSatish Balay     {comm_init();}
774827bd09bSSatish Balay #endif
775827bd09bSSatish Balay 
776827bd09bSSatish Balay   /* if there's nothing to do return */
777827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
778827bd09bSSatish Balay     {return;}
779827bd09bSSatish Balay 
780827bd09bSSatish Balay   /* the error msg says it all!!! */
781827bd09bSSatish Balay   if (modfl_num_nodes)
782827bd09bSSatish Balay     {error_msg_fatal("giop_hc() :: num_nodes not a power of 2!?!");}
783827bd09bSSatish Balay 
784827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
785827bd09bSSatish Balay   if (n<0)
78677431f27SBarry Smith     {error_msg_fatal("giop_hc() :: n=%D<0?",n);}
787827bd09bSSatish Balay 
788827bd09bSSatish Balay   /* can't do more dimensions then exist */
78939945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
790827bd09bSSatish Balay 
791827bd09bSSatish Balay   /* advance to list of n operations for custom */
792827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
793827bd09bSSatish Balay     {oprs++;}
794827bd09bSSatish Balay 
795d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))){
796827bd09bSSatish Balay     error_msg_warning("giop_hc() :: hope you passed in a rbfp!\n");
797827bd09bSSatish Balay     fp = (vfp) oprs;
798827bd09bSSatish Balay   }
799827bd09bSSatish Balay 
800827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
801827bd09bSSatish Balay     {
802827bd09bSSatish Balay       dest = my_id^mask;
803827bd09bSSatish Balay       if (my_id > dest)
804a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
805827bd09bSSatish Balay       else
806827bd09bSSatish Balay 	{
807a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,
808827bd09bSSatish Balay 		   &status);
809827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
810827bd09bSSatish Balay 	}
811827bd09bSSatish Balay     }
812827bd09bSSatish Balay 
813827bd09bSSatish Balay   if (edge==dim)
814827bd09bSSatish Balay     {mask>>=1;}
815827bd09bSSatish Balay   else
816827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
817827bd09bSSatish Balay 
818827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
819827bd09bSSatish Balay     {
820827bd09bSSatish Balay       if (my_id%mask)
821827bd09bSSatish Balay 	{continue;}
822827bd09bSSatish Balay 
823827bd09bSSatish Balay       dest = my_id^mask;
824827bd09bSSatish Balay       if (my_id < dest)
825a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
826827bd09bSSatish Balay       else
827827bd09bSSatish Balay 	{
828a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,
829827bd09bSSatish Balay 		   &status);
830827bd09bSSatish Balay 	}
831827bd09bSSatish Balay     }
832827bd09bSSatish Balay }
833