xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision a501084f0ec6f311bc0b878b7b6f8a6508165e40)
1827bd09bSSatish Balay 
2827bd09bSSatish Balay /***********************************comm.c*************************************
3827bd09bSSatish Balay 
4827bd09bSSatish Balay Author: Henry M. Tufo III
5827bd09bSSatish Balay 
6827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
7827bd09bSSatish Balay 
8827bd09bSSatish Balay snail-mail:
9827bd09bSSatish Balay Division of Applied Mathematics
10827bd09bSSatish Balay Brown University
11827bd09bSSatish Balay Providence, RI 02912
12827bd09bSSatish Balay 
13827bd09bSSatish Balay Last Modification:
14827bd09bSSatish Balay 11.21.97
15827bd09bSSatish Balay ***********************************comm.c*************************************/
16827bd09bSSatish Balay 
17827bd09bSSatish Balay /***********************************comm.c*************************************
18827bd09bSSatish Balay File Description:
19827bd09bSSatish Balay -----------------
20827bd09bSSatish Balay 
21827bd09bSSatish Balay ***********************************comm.c*************************************/
22d890fc11SSatish Balay #include "petsc.h"
23827bd09bSSatish Balay 
24827bd09bSSatish Balay #include "const.h"
25827bd09bSSatish Balay #include "types.h"
26827bd09bSSatish Balay #include "error.h"
27827bd09bSSatish Balay #include "ivec.h"
28827bd09bSSatish Balay #include "comm.h"
29827bd09bSSatish Balay 
30827bd09bSSatish Balay 
31827bd09bSSatish Balay /* global program control variables - explicitly exported */
32827bd09bSSatish Balay int my_id            = 0;
33827bd09bSSatish Balay int num_nodes        = 1;
34827bd09bSSatish Balay int floor_num_nodes  = 0;
35827bd09bSSatish Balay int i_log2_num_nodes = 0;
36827bd09bSSatish Balay 
37827bd09bSSatish Balay /* global program control variables */
38827bd09bSSatish Balay static int p_init = 0;
39827bd09bSSatish Balay static int modfl_num_nodes;
40827bd09bSSatish Balay static int edge_not_pow_2;
41827bd09bSSatish Balay 
42*a501084fSBarry Smith static unsigned int edge_node[sizeof(PetscInt)*32];
43827bd09bSSatish Balay 
44827bd09bSSatish Balay /***********************************comm.c*************************************
45827bd09bSSatish Balay Function: giop()
46827bd09bSSatish Balay 
47827bd09bSSatish Balay Input :
48827bd09bSSatish Balay Output:
49827bd09bSSatish Balay Return:
50827bd09bSSatish Balay Description:
51827bd09bSSatish Balay ***********************************comm.c*************************************/
52827bd09bSSatish Balay void
53827bd09bSSatish Balay comm_init (void)
54827bd09bSSatish Balay {
55827bd09bSSatish Balay 
56827bd09bSSatish Balay   if (p_init++) return;
57827bd09bSSatish Balay 
58827bd09bSSatish Balay #if PETSC_SIZEOF_INT != 4
59827bd09bSSatish Balay   error_msg_warning("Int != 4 Bytes!");
60827bd09bSSatish Balay #endif
61827bd09bSSatish Balay 
62827bd09bSSatish Balay 
63827bd09bSSatish Balay   MPI_Comm_size(MPI_COMM_WORLD,&num_nodes);
64827bd09bSSatish Balay   MPI_Comm_rank(MPI_COMM_WORLD,&my_id);
65827bd09bSSatish Balay 
66827bd09bSSatish Balay   if (num_nodes> (INT_MAX >> 1))
67827bd09bSSatish Balay   {error_msg_fatal("Can't have more then MAX_INT/2 nodes!!!");}
68827bd09bSSatish Balay 
69*a501084fSBarry Smith   ivec_zero((int*)edge_node,sizeof(PetscInt)*32);
70827bd09bSSatish Balay 
71827bd09bSSatish Balay   floor_num_nodes = 1;
72827bd09bSSatish Balay   i_log2_num_nodes = modfl_num_nodes = 0;
73827bd09bSSatish Balay   while (floor_num_nodes <= num_nodes)
74827bd09bSSatish Balay     {
75827bd09bSSatish Balay       edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes;
76827bd09bSSatish Balay       floor_num_nodes <<= 1;
77827bd09bSSatish Balay       i_log2_num_nodes++;
78827bd09bSSatish Balay     }
79827bd09bSSatish Balay 
80827bd09bSSatish Balay   i_log2_num_nodes--;
81827bd09bSSatish Balay   floor_num_nodes >>= 1;
82827bd09bSSatish Balay   modfl_num_nodes = (num_nodes - floor_num_nodes);
83827bd09bSSatish Balay 
84827bd09bSSatish Balay   if ((my_id > 0) && (my_id <= modfl_num_nodes))
85827bd09bSSatish Balay     {edge_not_pow_2=((my_id|floor_num_nodes)-1);}
86827bd09bSSatish Balay   else if (my_id >= floor_num_nodes)
87827bd09bSSatish Balay     {edge_not_pow_2=((my_id^floor_num_nodes)+1);
88827bd09bSSatish Balay     }
89827bd09bSSatish Balay   else
90827bd09bSSatish Balay     {edge_not_pow_2 = 0;}
91827bd09bSSatish Balay 
92827bd09bSSatish Balay }
93827bd09bSSatish Balay 
94827bd09bSSatish Balay 
95827bd09bSSatish Balay 
96827bd09bSSatish Balay /***********************************comm.c*************************************
97827bd09bSSatish Balay Function: giop()
98827bd09bSSatish Balay 
99827bd09bSSatish Balay Input :
100827bd09bSSatish Balay Output:
101827bd09bSSatish Balay Return:
102827bd09bSSatish Balay Description: fan-in/out version
103827bd09bSSatish Balay ***********************************comm.c*************************************/
104827bd09bSSatish Balay void
105827bd09bSSatish Balay giop(int *vals, int *work, int n, int *oprs)
106827bd09bSSatish Balay {
107*a501084fSBarry Smith    int mask, edge;
108827bd09bSSatish Balay   int type, dest;
109827bd09bSSatish Balay   vfp fp;
110827bd09bSSatish Balay   MPI_Status  status;
111827bd09bSSatish Balay 
112827bd09bSSatish Balay 
113827bd09bSSatish Balay #ifdef SAFE
114827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
115827bd09bSSatish Balay   if (!vals||!work||!oprs)
11677431f27SBarry Smith     {error_msg_fatal("giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
117827bd09bSSatish Balay 
118827bd09bSSatish Balay   /* non-uniform should have at least two entries */
119827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
120827bd09bSSatish Balay     {error_msg_fatal("giop() :: non_uniform and n=0,1?");}
121827bd09bSSatish Balay 
122827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
123827bd09bSSatish Balay   if (!p_init)
124827bd09bSSatish Balay     {comm_init();}
125827bd09bSSatish Balay #endif
126827bd09bSSatish Balay 
127827bd09bSSatish Balay   /* if there's nothing to do return */
128827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
129827bd09bSSatish Balay     {
130827bd09bSSatish Balay       return;
131827bd09bSSatish Balay     }
132827bd09bSSatish Balay 
133827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
134827bd09bSSatish Balay   if (n<0)
13577431f27SBarry Smith     {error_msg_fatal("giop() :: n=%D<0?",n);}
136827bd09bSSatish Balay 
137827bd09bSSatish Balay   /* advance to list of n operations for custom */
138827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
139827bd09bSSatish Balay     {oprs++;}
140827bd09bSSatish Balay 
141827bd09bSSatish Balay   /* major league hack */
142d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))) {
143827bd09bSSatish Balay     error_msg_warning("giop() :: hope you passed in a rbfp!\n");
144827bd09bSSatish Balay     fp = (vfp) oprs;
145827bd09bSSatish Balay   }
146827bd09bSSatish Balay 
147827bd09bSSatish Balay   /* all msgs will be of the same length */
148827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
149827bd09bSSatish Balay   if (edge_not_pow_2)
150827bd09bSSatish Balay     {
151827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
152827bd09bSSatish Balay 	{MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);}
153827bd09bSSatish Balay       else
154827bd09bSSatish Balay 	{
155827bd09bSSatish Balay 	  MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
156827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
157827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
158827bd09bSSatish Balay 	}
159827bd09bSSatish Balay     }
160827bd09bSSatish Balay 
161827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
162827bd09bSSatish Balay   if (my_id<floor_num_nodes)
163827bd09bSSatish Balay     {
164827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
165827bd09bSSatish Balay 	{
166827bd09bSSatish Balay 	  dest = my_id^mask;
167827bd09bSSatish Balay 	  if (my_id > dest)
168827bd09bSSatish Balay 	    {MPI_Send(vals,n,MPI_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
169827bd09bSSatish Balay 	  else
170827bd09bSSatish Balay 	    {
171827bd09bSSatish Balay 	      MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG2+dest,
172827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
173827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
174827bd09bSSatish Balay 	    }
175827bd09bSSatish Balay 	}
176827bd09bSSatish Balay 
177827bd09bSSatish Balay       mask=floor_num_nodes>>1;
178827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
179827bd09bSSatish Balay 	{
180827bd09bSSatish Balay 	  if (my_id%mask)
181827bd09bSSatish Balay 	    {continue;}
182827bd09bSSatish Balay 
183827bd09bSSatish Balay 	  dest = my_id^mask;
184827bd09bSSatish Balay 	  if (my_id < dest)
185827bd09bSSatish Balay 	    {MPI_Send(vals,n,MPI_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
186827bd09bSSatish Balay 	  else
187827bd09bSSatish Balay 	    {
188827bd09bSSatish Balay 	      MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG4+dest,
189827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
190827bd09bSSatish Balay 	    }
191827bd09bSSatish Balay 	}
192827bd09bSSatish Balay     }
193827bd09bSSatish Balay 
194827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
195827bd09bSSatish Balay   if (edge_not_pow_2)
196827bd09bSSatish Balay     {
197827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
198827bd09bSSatish Balay 	{
199827bd09bSSatish Balay 	  MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
200827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
201827bd09bSSatish Balay 	}
202827bd09bSSatish Balay       else
203827bd09bSSatish Balay 	{MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);}
204827bd09bSSatish Balay     }
205827bd09bSSatish Balay }
206827bd09bSSatish Balay 
207827bd09bSSatish Balay /***********************************comm.c*************************************
208827bd09bSSatish Balay Function: grop()
209827bd09bSSatish Balay 
210827bd09bSSatish Balay Input :
211827bd09bSSatish Balay Output:
212827bd09bSSatish Balay Return:
213827bd09bSSatish Balay Description: fan-in/out version
214827bd09bSSatish Balay ***********************************comm.c*************************************/
215827bd09bSSatish Balay void
216*a501084fSBarry Smith grop(PetscScalar *vals, PetscScalar *work, int n, int *oprs)
217827bd09bSSatish Balay {
218*a501084fSBarry Smith    int mask, edge;
219827bd09bSSatish Balay   int type, dest;
220827bd09bSSatish Balay   vfp fp;
221827bd09bSSatish Balay   MPI_Status  status;
222827bd09bSSatish Balay 
223827bd09bSSatish Balay 
224827bd09bSSatish Balay #ifdef SAFE
225827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
226827bd09bSSatish Balay   if (!vals||!work||!oprs)
22777431f27SBarry Smith     {error_msg_fatal("grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
228827bd09bSSatish Balay 
229827bd09bSSatish Balay   /* non-uniform should have at least two entries */
230827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
231827bd09bSSatish Balay     {error_msg_fatal("grop() :: non_uniform and n=0,1?");}
232827bd09bSSatish Balay 
233827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
234827bd09bSSatish Balay   if (!p_init)
235827bd09bSSatish Balay     {comm_init();}
236827bd09bSSatish Balay #endif
237827bd09bSSatish Balay 
238827bd09bSSatish Balay   /* if there's nothing to do return */
239827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
240827bd09bSSatish Balay     {return;}
241827bd09bSSatish Balay 
242827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
243827bd09bSSatish Balay   if (n<0)
24477431f27SBarry Smith     {error_msg_fatal("gdop() :: n=%D<0?",n);}
245827bd09bSSatish Balay 
246827bd09bSSatish Balay   /* advance to list of n operations for custom */
247827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
248827bd09bSSatish Balay     {oprs++;}
249827bd09bSSatish Balay 
250d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
251827bd09bSSatish Balay     error_msg_warning("grop() :: hope you passed in a rbfp!\n");
252827bd09bSSatish Balay     fp = (vfp) oprs;
253827bd09bSSatish Balay   }
254827bd09bSSatish Balay 
255827bd09bSSatish Balay   /* all msgs will be of the same length */
256827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
257827bd09bSSatish Balay   if (edge_not_pow_2)
258827bd09bSSatish Balay     {
259827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
260*a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+my_id,
261827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
262827bd09bSSatish Balay       else
263827bd09bSSatish Balay 	{
264*a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
265827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
266827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
267827bd09bSSatish Balay 	}
268827bd09bSSatish Balay     }
269827bd09bSSatish Balay 
270827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
271827bd09bSSatish Balay   if (my_id<floor_num_nodes)
272827bd09bSSatish Balay     {
273827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
274827bd09bSSatish Balay 	{
275827bd09bSSatish Balay 	  dest = my_id^mask;
276827bd09bSSatish Balay 	  if (my_id > dest)
277*a501084fSBarry Smith 	    {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
278827bd09bSSatish Balay 	  else
279827bd09bSSatish Balay 	    {
280*a501084fSBarry Smith 	      MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,
281827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
282827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
283827bd09bSSatish Balay 	    }
284827bd09bSSatish Balay 	}
285827bd09bSSatish Balay 
286827bd09bSSatish Balay       mask=floor_num_nodes>>1;
287827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
288827bd09bSSatish Balay 	{
289827bd09bSSatish Balay 	  if (my_id%mask)
290827bd09bSSatish Balay 	    {continue;}
291827bd09bSSatish Balay 
292827bd09bSSatish Balay 	  dest = my_id^mask;
293827bd09bSSatish Balay 	  if (my_id < dest)
294*a501084fSBarry Smith 	    {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
295827bd09bSSatish Balay 	  else
296827bd09bSSatish Balay 	    {
297*a501084fSBarry Smith 	      MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,
298827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
299827bd09bSSatish Balay 	    }
300827bd09bSSatish Balay 	}
301827bd09bSSatish Balay     }
302827bd09bSSatish Balay 
303827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
304827bd09bSSatish Balay   if (edge_not_pow_2)
305827bd09bSSatish Balay     {
306827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
307827bd09bSSatish Balay 	{
308*a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
309827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
310827bd09bSSatish Balay 	}
311827bd09bSSatish Balay       else
312*a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+my_id,
313827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
314827bd09bSSatish Balay     }
315827bd09bSSatish Balay }
316827bd09bSSatish Balay 
317827bd09bSSatish Balay 
318827bd09bSSatish Balay /***********************************comm.c*************************************
319827bd09bSSatish Balay Function: grop()
320827bd09bSSatish Balay 
321827bd09bSSatish Balay Input :
322827bd09bSSatish Balay Output:
323827bd09bSSatish Balay Return:
324827bd09bSSatish Balay Description: fan-in/out version
325827bd09bSSatish Balay 
326827bd09bSSatish Balay note good only for num_nodes=2^k!!!
327827bd09bSSatish Balay 
328827bd09bSSatish Balay ***********************************comm.c*************************************/
329827bd09bSSatish Balay void
330*a501084fSBarry Smith grop_hc(PetscScalar *vals, PetscScalar *work, int n, int *oprs, int dim)
331827bd09bSSatish Balay {
332*a501084fSBarry Smith    int mask, edge;
333827bd09bSSatish Balay   int type, dest;
334827bd09bSSatish Balay   vfp fp;
335827bd09bSSatish Balay   MPI_Status  status;
336827bd09bSSatish Balay 
337827bd09bSSatish Balay #ifdef SAFE
338827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
339827bd09bSSatish Balay   if (!vals||!work||!oprs)
34077431f27SBarry Smith     {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
341827bd09bSSatish Balay 
342827bd09bSSatish Balay   /* non-uniform should have at least two entries */
343827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
344827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: non_uniform and n=0,1?");}
345827bd09bSSatish Balay 
346827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
347827bd09bSSatish Balay   if (!p_init)
348827bd09bSSatish Balay     {comm_init();}
349827bd09bSSatish Balay #endif
350827bd09bSSatish Balay 
351827bd09bSSatish Balay   /* if there's nothing to do return */
352827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
353827bd09bSSatish Balay     {return;}
354827bd09bSSatish Balay 
355827bd09bSSatish Balay   /* the error msg says it all!!! */
356827bd09bSSatish Balay   if (modfl_num_nodes)
357827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");}
358827bd09bSSatish Balay 
359827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
360827bd09bSSatish Balay   if (n<0)
36177431f27SBarry Smith     {error_msg_fatal("grop_hc() :: n=%D<0?",n);}
362827bd09bSSatish Balay 
363827bd09bSSatish Balay   /* can't do more dimensions then exist */
36439945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
365827bd09bSSatish Balay 
366827bd09bSSatish Balay   /* advance to list of n operations for custom */
367827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
368827bd09bSSatish Balay     {oprs++;}
369827bd09bSSatish Balay 
370d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
371827bd09bSSatish Balay     error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n");
372827bd09bSSatish Balay     fp = (vfp) oprs;
373827bd09bSSatish Balay   }
374827bd09bSSatish Balay 
375827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
376827bd09bSSatish Balay     {
377827bd09bSSatish Balay       dest = my_id^mask;
378827bd09bSSatish Balay       if (my_id > dest)
379*a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
380827bd09bSSatish Balay       else
381827bd09bSSatish Balay 	{
382*a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,
383827bd09bSSatish Balay 		   &status);
384827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
385827bd09bSSatish Balay 	}
386827bd09bSSatish Balay     }
387827bd09bSSatish Balay 
388827bd09bSSatish Balay   if (edge==dim)
389827bd09bSSatish Balay     {mask>>=1;}
390827bd09bSSatish Balay   else
391827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
392827bd09bSSatish Balay 
393827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
394827bd09bSSatish Balay     {
395827bd09bSSatish Balay       if (my_id%mask)
396827bd09bSSatish Balay 	{continue;}
397827bd09bSSatish Balay 
398827bd09bSSatish Balay       dest = my_id^mask;
399827bd09bSSatish Balay       if (my_id < dest)
400*a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
401827bd09bSSatish Balay       else
402827bd09bSSatish Balay 	{
403*a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,
404827bd09bSSatish Balay 		   &status);
405827bd09bSSatish Balay 	}
406827bd09bSSatish Balay     }
407827bd09bSSatish Balay }
408827bd09bSSatish Balay 
409827bd09bSSatish Balay 
410827bd09bSSatish Balay /***********************************comm.c*************************************
411827bd09bSSatish Balay Function: gop()
412827bd09bSSatish Balay 
413827bd09bSSatish Balay Input :
414827bd09bSSatish Balay Output:
415827bd09bSSatish Balay Return:
416827bd09bSSatish Balay Description: fan-in/out version
417827bd09bSSatish Balay ***********************************comm.c*************************************/
418*a501084fSBarry Smith void gfop(void *vals, void *work, int n, vbfp fp, MPI_Datatype dt, int comm_type)
419827bd09bSSatish Balay {
420*a501084fSBarry Smith    int mask, edge;
421827bd09bSSatish Balay   int dest;
422827bd09bSSatish Balay   MPI_Status  status;
423827bd09bSSatish Balay   MPI_Op op;
424827bd09bSSatish Balay 
425827bd09bSSatish Balay #ifdef SAFE
426827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
427827bd09bSSatish Balay   if (!p_init)
428827bd09bSSatish Balay     {comm_init();}
429827bd09bSSatish Balay 
430827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
431827bd09bSSatish Balay   if (!vals||!work||!fp)
43277431f27SBarry Smith     {error_msg_fatal("gop() :: v=%D, w=%D, f=%D",vals,work,fp);}
433827bd09bSSatish Balay #endif
434827bd09bSSatish Balay 
435827bd09bSSatish Balay   /* if there's nothing to do return */
436827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
437827bd09bSSatish Balay     {return;}
438827bd09bSSatish Balay 
439827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
440827bd09bSSatish Balay   if (n<0)
44177431f27SBarry Smith     {error_msg_fatal("gop() :: n=%D<0?",n);}
442827bd09bSSatish Balay 
443827bd09bSSatish Balay   if (comm_type==MPI)
444827bd09bSSatish Balay     {
445827bd09bSSatish Balay       MPI_Op_create(fp,TRUE,&op);
446827bd09bSSatish Balay       MPI_Allreduce (vals, work, n, dt, op, MPI_COMM_WORLD);
447827bd09bSSatish Balay       MPI_Op_free(&op);
448827bd09bSSatish Balay       return;
449827bd09bSSatish Balay     }
450827bd09bSSatish Balay 
451827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
452827bd09bSSatish Balay   if (edge_not_pow_2)
453827bd09bSSatish Balay     {
454827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
455827bd09bSSatish Balay 	{MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG0+my_id,
456827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
457827bd09bSSatish Balay       else
458827bd09bSSatish Balay 	{
459827bd09bSSatish Balay 	  MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
460827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
461827bd09bSSatish Balay 	  (*fp)(vals,work,&n,&dt);
462827bd09bSSatish Balay 	}
463827bd09bSSatish Balay     }
464827bd09bSSatish Balay 
465827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
466827bd09bSSatish Balay   if (my_id<floor_num_nodes)
467827bd09bSSatish Balay     {
468827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
469827bd09bSSatish Balay 	{
470827bd09bSSatish Balay 	  dest = my_id^mask;
471827bd09bSSatish Balay 	  if (my_id > dest)
472827bd09bSSatish Balay 	    {MPI_Send(vals,n,dt,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
473827bd09bSSatish Balay 	  else
474827bd09bSSatish Balay 	    {
475827bd09bSSatish Balay 	      MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG2+dest,
476827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
477827bd09bSSatish Balay 	      (*fp)(vals, work, &n, &dt);
478827bd09bSSatish Balay 	    }
479827bd09bSSatish Balay 	}
480827bd09bSSatish Balay 
481827bd09bSSatish Balay       mask=floor_num_nodes>>1;
482827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
483827bd09bSSatish Balay 	{
484827bd09bSSatish Balay 	  if (my_id%mask)
485827bd09bSSatish Balay 	    {continue;}
486827bd09bSSatish Balay 
487827bd09bSSatish Balay 	  dest = my_id^mask;
488827bd09bSSatish Balay 	  if (my_id < dest)
489827bd09bSSatish Balay 	    {MPI_Send(vals,n,dt,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
490827bd09bSSatish Balay 	  else
491827bd09bSSatish Balay 	    {
492827bd09bSSatish Balay 	      MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG4+dest,
493827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
494827bd09bSSatish Balay 	    }
495827bd09bSSatish Balay 	}
496827bd09bSSatish Balay     }
497827bd09bSSatish Balay 
498827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
499827bd09bSSatish Balay   if (edge_not_pow_2)
500827bd09bSSatish Balay     {
501827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
502827bd09bSSatish Balay 	{
503827bd09bSSatish Balay 	  MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
504827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
505827bd09bSSatish Balay 	}
506827bd09bSSatish Balay       else
507827bd09bSSatish Balay 	{MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG5+my_id,
508827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
509827bd09bSSatish Balay     }
510827bd09bSSatish Balay }
511827bd09bSSatish Balay 
512827bd09bSSatish Balay 
513827bd09bSSatish Balay 
514827bd09bSSatish Balay 
515827bd09bSSatish Balay 
516827bd09bSSatish Balay 
517827bd09bSSatish Balay /******************************************************************************
518827bd09bSSatish Balay Function: giop()
519827bd09bSSatish Balay 
520827bd09bSSatish Balay Input :
521827bd09bSSatish Balay Output:
522827bd09bSSatish Balay Return:
523827bd09bSSatish Balay Description:
524827bd09bSSatish Balay 
525827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii
526827bd09bSSatish Balay 
527827bd09bSSatish Balay ******************************************************************************/
528827bd09bSSatish Balay void
529*a501084fSBarry Smith ssgl_radd( PetscScalar *vals,  PetscScalar *work,  int level,
530*a501084fSBarry Smith 	   int *segs)
531827bd09bSSatish Balay {
532*a501084fSBarry Smith    int edge, type, dest, mask;
533*a501084fSBarry Smith    int stage_n;
534827bd09bSSatish Balay   MPI_Status  status;
535827bd09bSSatish Balay 
536827bd09bSSatish Balay #ifdef SAFE
537827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
538827bd09bSSatish Balay   if (!p_init)
539827bd09bSSatish Balay     {comm_init();}
540827bd09bSSatish Balay #endif
541827bd09bSSatish Balay 
542827bd09bSSatish Balay 
543827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
544827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
545827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
546827bd09bSSatish Balay     {
547827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
548827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
549827bd09bSSatish Balay 	{
550827bd09bSSatish Balay 	  dest = edge_node[edge];
551827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
552827bd09bSSatish Balay 	  if (my_id>dest)
553*a501084fSBarry Smith           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type,
554827bd09bSSatish Balay                       MPI_COMM_WORLD);}
555827bd09bSSatish Balay 	  else
556827bd09bSSatish Balay 	    {
557827bd09bSSatish Balay 	      type =  type - my_id + dest;
558*a501084fSBarry Smith               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,
559827bd09bSSatish Balay                        MPI_COMM_WORLD,&status);
560827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
561827bd09bSSatish Balay 	    }
562827bd09bSSatish Balay 	}
563827bd09bSSatish Balay       mask <<= 1;
564827bd09bSSatish Balay     }
565827bd09bSSatish Balay   mask>>=1;
566827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
567827bd09bSSatish Balay     {
568827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
569827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
570827bd09bSSatish Balay 	{
571827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
572827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
573827bd09bSSatish Balay 	  if (my_id<dest)
574*a501084fSBarry Smith             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,
575827bd09bSSatish Balay                       MPI_COMM_WORLD);}
576827bd09bSSatish Balay 	  else
577827bd09bSSatish Balay 	    {
578827bd09bSSatish Balay 	      type =  type - my_id + dest;
579*a501084fSBarry Smith               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,
580827bd09bSSatish Balay                        MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
581827bd09bSSatish Balay 	    }
582827bd09bSSatish Balay 	}
583827bd09bSSatish Balay       mask >>= 1;
584827bd09bSSatish Balay     }
585827bd09bSSatish Balay }
586827bd09bSSatish Balay 
587827bd09bSSatish Balay 
588827bd09bSSatish Balay 
589827bd09bSSatish Balay /***********************************comm.c*************************************
590827bd09bSSatish Balay Function: grop_hc_vvl()
591827bd09bSSatish Balay 
592827bd09bSSatish Balay Input :
593827bd09bSSatish Balay Output:
594827bd09bSSatish Balay Return:
595827bd09bSSatish Balay Description: fan-in/out version
596827bd09bSSatish Balay 
597827bd09bSSatish Balay note good only for num_nodes=2^k!!!
598827bd09bSSatish Balay 
599827bd09bSSatish Balay ***********************************comm.c*************************************/
600827bd09bSSatish Balay void
601*a501084fSBarry Smith grop_hc_vvl(PetscScalar *vals, PetscScalar *work, int *segs, int *oprs, int dim)
602827bd09bSSatish Balay {
603*a501084fSBarry Smith    int mask, edge, n;
604827bd09bSSatish Balay   int type, dest;
605827bd09bSSatish Balay   vfp fp;
606827bd09bSSatish Balay   MPI_Status  status;
607827bd09bSSatish Balay 
608827bd09bSSatish Balay   error_msg_fatal("grop_hc_vvl() :: is not working!\n");
609827bd09bSSatish Balay 
610827bd09bSSatish Balay #ifdef SAFE
611827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
612827bd09bSSatish Balay   if (!vals||!work||!oprs||!segs)
61377431f27SBarry Smith     {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
614827bd09bSSatish Balay 
615827bd09bSSatish Balay   /* non-uniform should have at least two entries */
616827bd09bSSatish Balay 
617827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
618827bd09bSSatish Balay   if (!p_init)
619827bd09bSSatish Balay     {comm_init();}
620827bd09bSSatish Balay #endif
621827bd09bSSatish Balay 
622827bd09bSSatish Balay   /* if there's nothing to do return */
623827bd09bSSatish Balay   if ((num_nodes<2)||(dim<=0))
624827bd09bSSatish Balay     {return;}
625827bd09bSSatish Balay 
626827bd09bSSatish Balay   /* the error msg says it all!!! */
627827bd09bSSatish Balay   if (modfl_num_nodes)
628827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");}
629827bd09bSSatish Balay 
630827bd09bSSatish Balay   /* can't do more dimensions then exist */
63139945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
632827bd09bSSatish Balay 
633827bd09bSSatish Balay   /* advance to list of n operations for custom */
634827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
635827bd09bSSatish Balay     {oprs++;}
636827bd09bSSatish Balay 
637d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))){
638827bd09bSSatish Balay     error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n");
639827bd09bSSatish Balay     fp = (vfp) oprs;
640827bd09bSSatish Balay   }
641827bd09bSSatish Balay 
642827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
643827bd09bSSatish Balay     {
644827bd09bSSatish Balay       n = segs[dim]-segs[edge];
645827bd09bSSatish Balay       dest = my_id^mask;
646827bd09bSSatish Balay       if (my_id > dest)
647*a501084fSBarry Smith 	{MPI_Send(vals+segs[edge],n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
648827bd09bSSatish Balay       else
649827bd09bSSatish Balay 	{
650*a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,
651827bd09bSSatish Balay 		   MPI_COMM_WORLD, &status);
652827bd09bSSatish Balay 	  (*fp)(vals+segs[edge], work, n, oprs);
653827bd09bSSatish Balay 	}
654827bd09bSSatish Balay     }
655827bd09bSSatish Balay 
656827bd09bSSatish Balay   if (edge==dim)
657827bd09bSSatish Balay     {mask>>=1;}
658827bd09bSSatish Balay   else
659827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
660827bd09bSSatish Balay 
661827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
662827bd09bSSatish Balay     {
663827bd09bSSatish Balay       if (my_id%mask)
664827bd09bSSatish Balay 	{continue;}
665827bd09bSSatish Balay 
666827bd09bSSatish Balay       n = (segs[dim]-segs[dim-1-edge]);
667827bd09bSSatish Balay 
668827bd09bSSatish Balay       dest = my_id^mask;
669827bd09bSSatish Balay       if (my_id < dest)
670*a501084fSBarry Smith 	{MPI_Send(vals+segs[dim-1-edge],n,MPIU_SCALAR,dest,MSGTAG4+my_id,
671827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
672827bd09bSSatish Balay       else
673827bd09bSSatish Balay 	{
674*a501084fSBarry Smith 	  MPI_Recv(vals+segs[dim-1-edge],n,MPIU_SCALAR,MPI_ANY_SOURCE,
675827bd09bSSatish Balay 		   MSGTAG4+dest,MPI_COMM_WORLD, &status);
676827bd09bSSatish Balay 	}
677827bd09bSSatish Balay     }
678827bd09bSSatish Balay }
679827bd09bSSatish Balay 
680827bd09bSSatish Balay /******************************************************************************
681827bd09bSSatish Balay Function: giop()
682827bd09bSSatish Balay 
683827bd09bSSatish Balay Input :
684827bd09bSSatish Balay Output:
685827bd09bSSatish Balay Return:
686827bd09bSSatish Balay Description:
687827bd09bSSatish Balay 
688827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii
689827bd09bSSatish Balay 
690827bd09bSSatish Balay ******************************************************************************/
691*a501084fSBarry Smith void new_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  int level, int *segs)
692827bd09bSSatish Balay {
693*a501084fSBarry Smith    int edge, type, dest, mask;
694*a501084fSBarry Smith    int stage_n;
695827bd09bSSatish Balay   MPI_Status  status;
696827bd09bSSatish Balay 
697827bd09bSSatish Balay #ifdef SAFE
698827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
699827bd09bSSatish Balay   if (!p_init)
700827bd09bSSatish Balay     {comm_init();}
701827bd09bSSatish Balay #endif
702827bd09bSSatish Balay 
703827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
704827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
705827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
706827bd09bSSatish Balay     {
707827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
708827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
709827bd09bSSatish Balay 	{
710827bd09bSSatish Balay 	  dest = edge_node[edge];
711827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
712827bd09bSSatish Balay 	  if (my_id>dest)
713*a501084fSBarry Smith           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type,
714c700b191SBarry Smith                       MPI_COMM_WORLD);}
715827bd09bSSatish Balay 	  else
716827bd09bSSatish Balay 	    {
717827bd09bSSatish Balay 	      type =  type - my_id + dest;
718*a501084fSBarry Smith               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,
719c700b191SBarry Smith                        MPI_COMM_WORLD,&status);
720827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
721827bd09bSSatish Balay 	    }
722827bd09bSSatish Balay 	}
723827bd09bSSatish Balay       mask <<= 1;
724827bd09bSSatish Balay     }
725827bd09bSSatish Balay   mask>>=1;
726827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
727827bd09bSSatish Balay     {
728827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
729827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
730827bd09bSSatish Balay 	{
731827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
732827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
733827bd09bSSatish Balay 	  if (my_id<dest)
734*a501084fSBarry Smith             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,
735c700b191SBarry Smith                       MPI_COMM_WORLD);}
736827bd09bSSatish Balay 	  else
737827bd09bSSatish Balay 	    {
738827bd09bSSatish Balay 	      type =  type - my_id + dest;
739*a501084fSBarry Smith               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,
740c700b191SBarry Smith                        MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
741827bd09bSSatish Balay 	    }
742827bd09bSSatish Balay 	}
743827bd09bSSatish Balay       mask >>= 1;
744827bd09bSSatish Balay     }
745827bd09bSSatish Balay }
746827bd09bSSatish Balay 
747827bd09bSSatish Balay 
748827bd09bSSatish Balay 
749827bd09bSSatish Balay /***********************************comm.c*************************************
750827bd09bSSatish Balay Function: giop()
751827bd09bSSatish Balay 
752827bd09bSSatish Balay Input :
753827bd09bSSatish Balay Output:
754827bd09bSSatish Balay Return:
755827bd09bSSatish Balay Description: fan-in/out version
756827bd09bSSatish Balay 
757827bd09bSSatish Balay note good only for num_nodes=2^k!!!
758827bd09bSSatish Balay 
759827bd09bSSatish Balay ***********************************comm.c*************************************/
760c700b191SBarry Smith void giop_hc(int *vals, int *work, int n, int *oprs, int dim)
761827bd09bSSatish Balay {
762*a501084fSBarry Smith    int mask, edge;
763827bd09bSSatish Balay   int type, dest;
764827bd09bSSatish Balay   vfp fp;
765827bd09bSSatish Balay   MPI_Status  status;
766827bd09bSSatish Balay 
767827bd09bSSatish Balay #ifdef SAFE
768827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
769827bd09bSSatish Balay   if (!vals||!work||!oprs)
77077431f27SBarry Smith     {error_msg_fatal("giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
771827bd09bSSatish Balay 
772827bd09bSSatish Balay   /* non-uniform should have at least two entries */
773827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
774827bd09bSSatish Balay     {error_msg_fatal("giop_hc() :: non_uniform and n=0,1?");}
775827bd09bSSatish Balay 
776827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
777827bd09bSSatish Balay   if (!p_init)
778827bd09bSSatish Balay     {comm_init();}
779827bd09bSSatish Balay #endif
780827bd09bSSatish Balay 
781827bd09bSSatish Balay   /* if there's nothing to do return */
782827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
783827bd09bSSatish Balay     {return;}
784827bd09bSSatish Balay 
785827bd09bSSatish Balay   /* the error msg says it all!!! */
786827bd09bSSatish Balay   if (modfl_num_nodes)
787827bd09bSSatish Balay     {error_msg_fatal("giop_hc() :: num_nodes not a power of 2!?!");}
788827bd09bSSatish Balay 
789827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
790827bd09bSSatish Balay   if (n<0)
79177431f27SBarry Smith     {error_msg_fatal("giop_hc() :: n=%D<0?",n);}
792827bd09bSSatish Balay 
793827bd09bSSatish Balay   /* can't do more dimensions then exist */
79439945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
795827bd09bSSatish Balay 
796827bd09bSSatish Balay   /* advance to list of n operations for custom */
797827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
798827bd09bSSatish Balay     {oprs++;}
799827bd09bSSatish Balay 
800d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))){
801827bd09bSSatish Balay     error_msg_warning("giop_hc() :: hope you passed in a rbfp!\n");
802827bd09bSSatish Balay     fp = (vfp) oprs;
803827bd09bSSatish Balay   }
804827bd09bSSatish Balay 
805827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
806827bd09bSSatish Balay     {
807827bd09bSSatish Balay       dest = my_id^mask;
808827bd09bSSatish Balay       if (my_id > dest)
809*a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
810827bd09bSSatish Balay       else
811827bd09bSSatish Balay 	{
812*a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,
813827bd09bSSatish Balay 		   &status);
814827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
815827bd09bSSatish Balay 	}
816827bd09bSSatish Balay     }
817827bd09bSSatish Balay 
818827bd09bSSatish Balay   if (edge==dim)
819827bd09bSSatish Balay     {mask>>=1;}
820827bd09bSSatish Balay   else
821827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
822827bd09bSSatish Balay 
823827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
824827bd09bSSatish Balay     {
825827bd09bSSatish Balay       if (my_id%mask)
826827bd09bSSatish Balay 	{continue;}
827827bd09bSSatish Balay 
828827bd09bSSatish Balay       dest = my_id^mask;
829827bd09bSSatish Balay       if (my_id < dest)
830*a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
831827bd09bSSatish Balay       else
832827bd09bSSatish Balay 	{
833*a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,
834827bd09bSSatish Balay 		   &status);
835827bd09bSSatish Balay 	}
836827bd09bSSatish Balay     }
837827bd09bSSatish Balay }
838