xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision 7758a8cdabe51840ad55576361d16be55c4795b5)
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*************************************/
22*7758a8cdSBarry Smith #include "src/ksp/pc/impls/tfs/tfs.h"
23827bd09bSSatish Balay 
24827bd09bSSatish Balay 
25827bd09bSSatish Balay /* global program control variables - explicitly exported */
26827bd09bSSatish Balay int my_id            = 0;
27827bd09bSSatish Balay int num_nodes        = 1;
28827bd09bSSatish Balay int floor_num_nodes  = 0;
29827bd09bSSatish Balay int i_log2_num_nodes = 0;
30827bd09bSSatish Balay 
31827bd09bSSatish Balay /* global program control variables */
32827bd09bSSatish Balay static int p_init = 0;
33827bd09bSSatish Balay static int modfl_num_nodes;
34827bd09bSSatish Balay static int edge_not_pow_2;
35827bd09bSSatish Balay 
36a501084fSBarry Smith static unsigned int edge_node[sizeof(PetscInt)*32];
37827bd09bSSatish Balay 
38827bd09bSSatish Balay /***********************************comm.c*************************************
39827bd09bSSatish Balay Function: giop()
40827bd09bSSatish Balay 
41827bd09bSSatish Balay Input :
42827bd09bSSatish Balay Output:
43827bd09bSSatish Balay Return:
44827bd09bSSatish Balay Description:
45827bd09bSSatish Balay ***********************************comm.c*************************************/
46827bd09bSSatish Balay void
47827bd09bSSatish Balay comm_init (void)
48827bd09bSSatish Balay {
49827bd09bSSatish Balay 
50827bd09bSSatish Balay   if (p_init++) return;
51827bd09bSSatish Balay 
52827bd09bSSatish Balay #if PETSC_SIZEOF_INT != 4
53827bd09bSSatish Balay   error_msg_warning("Int != 4 Bytes!");
54827bd09bSSatish Balay #endif
55827bd09bSSatish Balay 
56827bd09bSSatish Balay 
57827bd09bSSatish Balay   MPI_Comm_size(MPI_COMM_WORLD,&num_nodes);
58827bd09bSSatish Balay   MPI_Comm_rank(MPI_COMM_WORLD,&my_id);
59827bd09bSSatish Balay 
60827bd09bSSatish Balay   if (num_nodes> (INT_MAX >> 1))
61827bd09bSSatish Balay   {error_msg_fatal("Can't have more then MAX_INT/2 nodes!!!");}
62827bd09bSSatish Balay 
63a501084fSBarry Smith   ivec_zero((int*)edge_node,sizeof(PetscInt)*32);
64827bd09bSSatish Balay 
65827bd09bSSatish Balay   floor_num_nodes = 1;
66827bd09bSSatish Balay   i_log2_num_nodes = modfl_num_nodes = 0;
67827bd09bSSatish Balay   while (floor_num_nodes <= num_nodes)
68827bd09bSSatish Balay     {
69827bd09bSSatish Balay       edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes;
70827bd09bSSatish Balay       floor_num_nodes <<= 1;
71827bd09bSSatish Balay       i_log2_num_nodes++;
72827bd09bSSatish Balay     }
73827bd09bSSatish Balay 
74827bd09bSSatish Balay   i_log2_num_nodes--;
75827bd09bSSatish Balay   floor_num_nodes >>= 1;
76827bd09bSSatish Balay   modfl_num_nodes = (num_nodes - floor_num_nodes);
77827bd09bSSatish Balay 
78827bd09bSSatish Balay   if ((my_id > 0) && (my_id <= modfl_num_nodes))
79827bd09bSSatish Balay     {edge_not_pow_2=((my_id|floor_num_nodes)-1);}
80827bd09bSSatish Balay   else if (my_id >= floor_num_nodes)
81827bd09bSSatish Balay     {edge_not_pow_2=((my_id^floor_num_nodes)+1);
82827bd09bSSatish Balay     }
83827bd09bSSatish Balay   else
84827bd09bSSatish Balay     {edge_not_pow_2 = 0;}
85827bd09bSSatish Balay 
86827bd09bSSatish Balay }
87827bd09bSSatish Balay 
88827bd09bSSatish Balay 
89827bd09bSSatish Balay 
90827bd09bSSatish Balay /***********************************comm.c*************************************
91827bd09bSSatish Balay Function: giop()
92827bd09bSSatish Balay 
93827bd09bSSatish Balay Input :
94827bd09bSSatish Balay Output:
95827bd09bSSatish Balay Return:
96827bd09bSSatish Balay Description: fan-in/out version
97827bd09bSSatish Balay ***********************************comm.c*************************************/
98827bd09bSSatish Balay void
99827bd09bSSatish Balay giop(int *vals, int *work, int n, int *oprs)
100827bd09bSSatish Balay {
101a501084fSBarry Smith    int mask, edge;
102827bd09bSSatish Balay   int type, dest;
103827bd09bSSatish Balay   vfp fp;
104827bd09bSSatish Balay   MPI_Status  status;
105827bd09bSSatish Balay 
106827bd09bSSatish Balay 
107827bd09bSSatish Balay #ifdef SAFE
108827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
109827bd09bSSatish Balay   if (!vals||!work||!oprs)
11077431f27SBarry Smith     {error_msg_fatal("giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
111827bd09bSSatish Balay 
112827bd09bSSatish Balay   /* non-uniform should have at least two entries */
113827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
114827bd09bSSatish Balay     {error_msg_fatal("giop() :: non_uniform and n=0,1?");}
115827bd09bSSatish Balay 
116827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
117827bd09bSSatish Balay   if (!p_init)
118827bd09bSSatish Balay     {comm_init();}
119827bd09bSSatish Balay #endif
120827bd09bSSatish Balay 
121827bd09bSSatish Balay   /* if there's nothing to do return */
122827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
123827bd09bSSatish Balay     {
124827bd09bSSatish Balay       return;
125827bd09bSSatish Balay     }
126827bd09bSSatish Balay 
127827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
128827bd09bSSatish Balay   if (n<0)
12977431f27SBarry Smith     {error_msg_fatal("giop() :: n=%D<0?",n);}
130827bd09bSSatish Balay 
131827bd09bSSatish Balay   /* advance to list of n operations for custom */
132827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
133827bd09bSSatish Balay     {oprs++;}
134827bd09bSSatish Balay 
135827bd09bSSatish Balay   /* major league hack */
136d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))) {
137827bd09bSSatish Balay     error_msg_warning("giop() :: hope you passed in a rbfp!\n");
138827bd09bSSatish Balay     fp = (vfp) oprs;
139827bd09bSSatish Balay   }
140827bd09bSSatish Balay 
141827bd09bSSatish Balay   /* all msgs will be of the same length */
142827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
143827bd09bSSatish Balay   if (edge_not_pow_2)
144827bd09bSSatish Balay     {
145827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
146827bd09bSSatish Balay 	{MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);}
147827bd09bSSatish Balay       else
148827bd09bSSatish Balay 	{
149827bd09bSSatish Balay 	  MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
150827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
151827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
152827bd09bSSatish Balay 	}
153827bd09bSSatish Balay     }
154827bd09bSSatish Balay 
155827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
156827bd09bSSatish Balay   if (my_id<floor_num_nodes)
157827bd09bSSatish Balay     {
158827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
159827bd09bSSatish Balay 	{
160827bd09bSSatish Balay 	  dest = my_id^mask;
161827bd09bSSatish Balay 	  if (my_id > dest)
162827bd09bSSatish Balay 	    {MPI_Send(vals,n,MPI_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
163827bd09bSSatish Balay 	  else
164827bd09bSSatish Balay 	    {
165827bd09bSSatish Balay 	      MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG2+dest,
166827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
167827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
168827bd09bSSatish Balay 	    }
169827bd09bSSatish Balay 	}
170827bd09bSSatish Balay 
171827bd09bSSatish Balay       mask=floor_num_nodes>>1;
172827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
173827bd09bSSatish Balay 	{
174827bd09bSSatish Balay 	  if (my_id%mask)
175827bd09bSSatish Balay 	    {continue;}
176827bd09bSSatish Balay 
177827bd09bSSatish Balay 	  dest = my_id^mask;
178827bd09bSSatish Balay 	  if (my_id < dest)
179827bd09bSSatish Balay 	    {MPI_Send(vals,n,MPI_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
180827bd09bSSatish Balay 	  else
181827bd09bSSatish Balay 	    {
182827bd09bSSatish Balay 	      MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG4+dest,
183827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
184827bd09bSSatish Balay 	    }
185827bd09bSSatish Balay 	}
186827bd09bSSatish Balay     }
187827bd09bSSatish Balay 
188827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
189827bd09bSSatish Balay   if (edge_not_pow_2)
190827bd09bSSatish Balay     {
191827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
192827bd09bSSatish Balay 	{
193827bd09bSSatish Balay 	  MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
194827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
195827bd09bSSatish Balay 	}
196827bd09bSSatish Balay       else
197827bd09bSSatish Balay 	{MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);}
198827bd09bSSatish Balay     }
199827bd09bSSatish Balay }
200827bd09bSSatish Balay 
201827bd09bSSatish Balay /***********************************comm.c*************************************
202827bd09bSSatish Balay Function: grop()
203827bd09bSSatish Balay 
204827bd09bSSatish Balay Input :
205827bd09bSSatish Balay Output:
206827bd09bSSatish Balay Return:
207827bd09bSSatish Balay Description: fan-in/out version
208827bd09bSSatish Balay ***********************************comm.c*************************************/
209827bd09bSSatish Balay void
210a501084fSBarry Smith grop(PetscScalar *vals, PetscScalar *work, int n, int *oprs)
211827bd09bSSatish Balay {
212a501084fSBarry Smith    int mask, edge;
213827bd09bSSatish Balay   int type, dest;
214827bd09bSSatish Balay   vfp fp;
215827bd09bSSatish Balay   MPI_Status  status;
216827bd09bSSatish Balay 
217827bd09bSSatish Balay 
218827bd09bSSatish Balay #ifdef SAFE
219827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
220827bd09bSSatish Balay   if (!vals||!work||!oprs)
22177431f27SBarry Smith     {error_msg_fatal("grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
222827bd09bSSatish Balay 
223827bd09bSSatish Balay   /* non-uniform should have at least two entries */
224827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
225827bd09bSSatish Balay     {error_msg_fatal("grop() :: non_uniform and n=0,1?");}
226827bd09bSSatish Balay 
227827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
228827bd09bSSatish Balay   if (!p_init)
229827bd09bSSatish Balay     {comm_init();}
230827bd09bSSatish Balay #endif
231827bd09bSSatish Balay 
232827bd09bSSatish Balay   /* if there's nothing to do return */
233827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
234827bd09bSSatish Balay     {return;}
235827bd09bSSatish Balay 
236827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
237827bd09bSSatish Balay   if (n<0)
23877431f27SBarry Smith     {error_msg_fatal("gdop() :: n=%D<0?",n);}
239827bd09bSSatish Balay 
240827bd09bSSatish Balay   /* advance to list of n operations for custom */
241827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
242827bd09bSSatish Balay     {oprs++;}
243827bd09bSSatish Balay 
244d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
245827bd09bSSatish Balay     error_msg_warning("grop() :: hope you passed in a rbfp!\n");
246827bd09bSSatish Balay     fp = (vfp) oprs;
247827bd09bSSatish Balay   }
248827bd09bSSatish Balay 
249827bd09bSSatish Balay   /* all msgs will be of the same length */
250827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
251827bd09bSSatish Balay   if (edge_not_pow_2)
252827bd09bSSatish Balay     {
253827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
254a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+my_id,
255827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
256827bd09bSSatish Balay       else
257827bd09bSSatish Balay 	{
258a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
259827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
260827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
261827bd09bSSatish Balay 	}
262827bd09bSSatish Balay     }
263827bd09bSSatish Balay 
264827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
265827bd09bSSatish Balay   if (my_id<floor_num_nodes)
266827bd09bSSatish Balay     {
267827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
268827bd09bSSatish Balay 	{
269827bd09bSSatish Balay 	  dest = my_id^mask;
270827bd09bSSatish Balay 	  if (my_id > dest)
271a501084fSBarry Smith 	    {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
272827bd09bSSatish Balay 	  else
273827bd09bSSatish Balay 	    {
274a501084fSBarry Smith 	      MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,
275827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
276827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
277827bd09bSSatish Balay 	    }
278827bd09bSSatish Balay 	}
279827bd09bSSatish Balay 
280827bd09bSSatish Balay       mask=floor_num_nodes>>1;
281827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
282827bd09bSSatish Balay 	{
283827bd09bSSatish Balay 	  if (my_id%mask)
284827bd09bSSatish Balay 	    {continue;}
285827bd09bSSatish Balay 
286827bd09bSSatish Balay 	  dest = my_id^mask;
287827bd09bSSatish Balay 	  if (my_id < dest)
288a501084fSBarry Smith 	    {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
289827bd09bSSatish Balay 	  else
290827bd09bSSatish Balay 	    {
291a501084fSBarry Smith 	      MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,
292827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
293827bd09bSSatish Balay 	    }
294827bd09bSSatish Balay 	}
295827bd09bSSatish Balay     }
296827bd09bSSatish Balay 
297827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
298827bd09bSSatish Balay   if (edge_not_pow_2)
299827bd09bSSatish Balay     {
300827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
301827bd09bSSatish Balay 	{
302a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
303827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
304827bd09bSSatish Balay 	}
305827bd09bSSatish Balay       else
306a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+my_id,
307827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
308827bd09bSSatish Balay     }
309827bd09bSSatish Balay }
310827bd09bSSatish Balay 
311827bd09bSSatish Balay 
312827bd09bSSatish Balay /***********************************comm.c*************************************
313827bd09bSSatish Balay Function: grop()
314827bd09bSSatish Balay 
315827bd09bSSatish Balay Input :
316827bd09bSSatish Balay Output:
317827bd09bSSatish Balay Return:
318827bd09bSSatish Balay Description: fan-in/out version
319827bd09bSSatish Balay 
320827bd09bSSatish Balay note good only for num_nodes=2^k!!!
321827bd09bSSatish Balay 
322827bd09bSSatish Balay ***********************************comm.c*************************************/
323827bd09bSSatish Balay void
324a501084fSBarry Smith grop_hc(PetscScalar *vals, PetscScalar *work, int n, int *oprs, int dim)
325827bd09bSSatish Balay {
326a501084fSBarry Smith    int mask, edge;
327827bd09bSSatish Balay   int type, dest;
328827bd09bSSatish Balay   vfp fp;
329827bd09bSSatish Balay   MPI_Status  status;
330827bd09bSSatish Balay 
331827bd09bSSatish Balay #ifdef SAFE
332827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
333827bd09bSSatish Balay   if (!vals||!work||!oprs)
33477431f27SBarry Smith     {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
335827bd09bSSatish Balay 
336827bd09bSSatish Balay   /* non-uniform should have at least two entries */
337827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
338827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: non_uniform and n=0,1?");}
339827bd09bSSatish Balay 
340827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
341827bd09bSSatish Balay   if (!p_init)
342827bd09bSSatish Balay     {comm_init();}
343827bd09bSSatish Balay #endif
344827bd09bSSatish Balay 
345827bd09bSSatish Balay   /* if there's nothing to do return */
346827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
347827bd09bSSatish Balay     {return;}
348827bd09bSSatish Balay 
349827bd09bSSatish Balay   /* the error msg says it all!!! */
350827bd09bSSatish Balay   if (modfl_num_nodes)
351827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");}
352827bd09bSSatish Balay 
353827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
354827bd09bSSatish Balay   if (n<0)
35577431f27SBarry Smith     {error_msg_fatal("grop_hc() :: n=%D<0?",n);}
356827bd09bSSatish Balay 
357827bd09bSSatish Balay   /* can't do more dimensions then exist */
35839945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
359827bd09bSSatish Balay 
360827bd09bSSatish Balay   /* advance to list of n operations for custom */
361827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
362827bd09bSSatish Balay     {oprs++;}
363827bd09bSSatish Balay 
364d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
365827bd09bSSatish Balay     error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n");
366827bd09bSSatish Balay     fp = (vfp) oprs;
367827bd09bSSatish Balay   }
368827bd09bSSatish Balay 
369827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
370827bd09bSSatish Balay     {
371827bd09bSSatish Balay       dest = my_id^mask;
372827bd09bSSatish Balay       if (my_id > dest)
373a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
374827bd09bSSatish Balay       else
375827bd09bSSatish Balay 	{
376a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,
377827bd09bSSatish Balay 		   &status);
378827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
379827bd09bSSatish Balay 	}
380827bd09bSSatish Balay     }
381827bd09bSSatish Balay 
382827bd09bSSatish Balay   if (edge==dim)
383827bd09bSSatish Balay     {mask>>=1;}
384827bd09bSSatish Balay   else
385827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
386827bd09bSSatish Balay 
387827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
388827bd09bSSatish Balay     {
389827bd09bSSatish Balay       if (my_id%mask)
390827bd09bSSatish Balay 	{continue;}
391827bd09bSSatish Balay 
392827bd09bSSatish Balay       dest = my_id^mask;
393827bd09bSSatish Balay       if (my_id < dest)
394a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
395827bd09bSSatish Balay       else
396827bd09bSSatish Balay 	{
397a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,
398827bd09bSSatish Balay 		   &status);
399827bd09bSSatish Balay 	}
400827bd09bSSatish Balay     }
401827bd09bSSatish Balay }
402827bd09bSSatish Balay 
403827bd09bSSatish Balay 
404827bd09bSSatish Balay /***********************************comm.c*************************************
405827bd09bSSatish Balay Function: gop()
406827bd09bSSatish Balay 
407827bd09bSSatish Balay Input :
408827bd09bSSatish Balay Output:
409827bd09bSSatish Balay Return:
410827bd09bSSatish Balay Description: fan-in/out version
411827bd09bSSatish Balay ***********************************comm.c*************************************/
412a501084fSBarry Smith void gfop(void *vals, void *work, int n, vbfp fp, MPI_Datatype dt, int comm_type)
413827bd09bSSatish Balay {
414a501084fSBarry Smith    int mask, edge;
415827bd09bSSatish Balay   int dest;
416827bd09bSSatish Balay   MPI_Status  status;
417827bd09bSSatish Balay   MPI_Op op;
418827bd09bSSatish Balay 
419827bd09bSSatish Balay #ifdef SAFE
420827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
421827bd09bSSatish Balay   if (!p_init)
422827bd09bSSatish Balay     {comm_init();}
423827bd09bSSatish Balay 
424827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
425827bd09bSSatish Balay   if (!vals||!work||!fp)
42677431f27SBarry Smith     {error_msg_fatal("gop() :: v=%D, w=%D, f=%D",vals,work,fp);}
427827bd09bSSatish Balay #endif
428827bd09bSSatish Balay 
429827bd09bSSatish Balay   /* if there's nothing to do return */
430827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
431827bd09bSSatish Balay     {return;}
432827bd09bSSatish Balay 
433827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
434827bd09bSSatish Balay   if (n<0)
43577431f27SBarry Smith     {error_msg_fatal("gop() :: n=%D<0?",n);}
436827bd09bSSatish Balay 
437827bd09bSSatish Balay   if (comm_type==MPI)
438827bd09bSSatish Balay     {
439827bd09bSSatish Balay       MPI_Op_create(fp,TRUE,&op);
440827bd09bSSatish Balay       MPI_Allreduce (vals, work, n, dt, op, MPI_COMM_WORLD);
441827bd09bSSatish Balay       MPI_Op_free(&op);
442827bd09bSSatish Balay       return;
443827bd09bSSatish Balay     }
444827bd09bSSatish Balay 
445827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
446827bd09bSSatish Balay   if (edge_not_pow_2)
447827bd09bSSatish Balay     {
448827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
449827bd09bSSatish Balay 	{MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG0+my_id,
450827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
451827bd09bSSatish Balay       else
452827bd09bSSatish Balay 	{
453827bd09bSSatish Balay 	  MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,
454827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
455827bd09bSSatish Balay 	  (*fp)(vals,work,&n,&dt);
456827bd09bSSatish Balay 	}
457827bd09bSSatish Balay     }
458827bd09bSSatish Balay 
459827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
460827bd09bSSatish Balay   if (my_id<floor_num_nodes)
461827bd09bSSatish Balay     {
462827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
463827bd09bSSatish Balay 	{
464827bd09bSSatish Balay 	  dest = my_id^mask;
465827bd09bSSatish Balay 	  if (my_id > dest)
466827bd09bSSatish Balay 	    {MPI_Send(vals,n,dt,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
467827bd09bSSatish Balay 	  else
468827bd09bSSatish Balay 	    {
469827bd09bSSatish Balay 	      MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG2+dest,
470827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
471827bd09bSSatish Balay 	      (*fp)(vals, work, &n, &dt);
472827bd09bSSatish Balay 	    }
473827bd09bSSatish Balay 	}
474827bd09bSSatish Balay 
475827bd09bSSatish Balay       mask=floor_num_nodes>>1;
476827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
477827bd09bSSatish Balay 	{
478827bd09bSSatish Balay 	  if (my_id%mask)
479827bd09bSSatish Balay 	    {continue;}
480827bd09bSSatish Balay 
481827bd09bSSatish Balay 	  dest = my_id^mask;
482827bd09bSSatish Balay 	  if (my_id < dest)
483827bd09bSSatish Balay 	    {MPI_Send(vals,n,dt,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
484827bd09bSSatish Balay 	  else
485827bd09bSSatish Balay 	    {
486827bd09bSSatish Balay 	      MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG4+dest,
487827bd09bSSatish Balay 		       MPI_COMM_WORLD, &status);
488827bd09bSSatish Balay 	    }
489827bd09bSSatish Balay 	}
490827bd09bSSatish Balay     }
491827bd09bSSatish Balay 
492827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
493827bd09bSSatish Balay   if (edge_not_pow_2)
494827bd09bSSatish Balay     {
495827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
496827bd09bSSatish Balay 	{
497827bd09bSSatish Balay 	  MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,
498827bd09bSSatish Balay 		   MPI_COMM_WORLD,&status);
499827bd09bSSatish Balay 	}
500827bd09bSSatish Balay       else
501827bd09bSSatish Balay 	{MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG5+my_id,
502827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
503827bd09bSSatish Balay     }
504827bd09bSSatish Balay }
505827bd09bSSatish Balay 
506827bd09bSSatish Balay 
507827bd09bSSatish Balay 
508827bd09bSSatish Balay 
509827bd09bSSatish Balay 
510827bd09bSSatish Balay 
511827bd09bSSatish Balay /******************************************************************************
512827bd09bSSatish Balay Function: giop()
513827bd09bSSatish Balay 
514827bd09bSSatish Balay Input :
515827bd09bSSatish Balay Output:
516827bd09bSSatish Balay Return:
517827bd09bSSatish Balay Description:
518827bd09bSSatish Balay 
519827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii
520827bd09bSSatish Balay 
521827bd09bSSatish Balay ******************************************************************************/
522827bd09bSSatish Balay void
523a501084fSBarry Smith ssgl_radd( PetscScalar *vals,  PetscScalar *work,  int level,
524a501084fSBarry Smith 	   int *segs)
525827bd09bSSatish Balay {
526a501084fSBarry Smith    int edge, type, dest, mask;
527a501084fSBarry Smith    int stage_n;
528827bd09bSSatish Balay   MPI_Status  status;
529827bd09bSSatish Balay 
530827bd09bSSatish Balay #ifdef SAFE
531827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
532827bd09bSSatish Balay   if (!p_init)
533827bd09bSSatish Balay     {comm_init();}
534827bd09bSSatish Balay #endif
535827bd09bSSatish Balay 
536827bd09bSSatish Balay 
537827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
538827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
539827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
540827bd09bSSatish Balay     {
541827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
542827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
543827bd09bSSatish Balay 	{
544827bd09bSSatish Balay 	  dest = edge_node[edge];
545827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
546827bd09bSSatish Balay 	  if (my_id>dest)
547a501084fSBarry Smith           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type,
548827bd09bSSatish Balay                       MPI_COMM_WORLD);}
549827bd09bSSatish Balay 	  else
550827bd09bSSatish Balay 	    {
551827bd09bSSatish Balay 	      type =  type - my_id + dest;
552a501084fSBarry Smith               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,
553827bd09bSSatish Balay                        MPI_COMM_WORLD,&status);
554827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
555827bd09bSSatish Balay 	    }
556827bd09bSSatish Balay 	}
557827bd09bSSatish Balay       mask <<= 1;
558827bd09bSSatish Balay     }
559827bd09bSSatish Balay   mask>>=1;
560827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
561827bd09bSSatish Balay     {
562827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
563827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
564827bd09bSSatish Balay 	{
565827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
566827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
567827bd09bSSatish Balay 	  if (my_id<dest)
568a501084fSBarry Smith             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,
569827bd09bSSatish Balay                       MPI_COMM_WORLD);}
570827bd09bSSatish Balay 	  else
571827bd09bSSatish Balay 	    {
572827bd09bSSatish Balay 	      type =  type - my_id + dest;
573a501084fSBarry Smith               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,
574827bd09bSSatish Balay                        MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
575827bd09bSSatish Balay 	    }
576827bd09bSSatish Balay 	}
577827bd09bSSatish Balay       mask >>= 1;
578827bd09bSSatish Balay     }
579827bd09bSSatish Balay }
580827bd09bSSatish Balay 
581827bd09bSSatish Balay 
582827bd09bSSatish Balay 
583827bd09bSSatish Balay /***********************************comm.c*************************************
584827bd09bSSatish Balay Function: grop_hc_vvl()
585827bd09bSSatish Balay 
586827bd09bSSatish Balay Input :
587827bd09bSSatish Balay Output:
588827bd09bSSatish Balay Return:
589827bd09bSSatish Balay Description: fan-in/out version
590827bd09bSSatish Balay 
591827bd09bSSatish Balay note good only for num_nodes=2^k!!!
592827bd09bSSatish Balay 
593827bd09bSSatish Balay ***********************************comm.c*************************************/
594827bd09bSSatish Balay void
595a501084fSBarry Smith grop_hc_vvl(PetscScalar *vals, PetscScalar *work, int *segs, int *oprs, int dim)
596827bd09bSSatish Balay {
597a501084fSBarry Smith    int mask, edge, n;
598827bd09bSSatish Balay   int type, dest;
599827bd09bSSatish Balay   vfp fp;
600827bd09bSSatish Balay   MPI_Status  status;
601827bd09bSSatish Balay 
602827bd09bSSatish Balay   error_msg_fatal("grop_hc_vvl() :: is not working!\n");
603827bd09bSSatish Balay 
604827bd09bSSatish Balay #ifdef SAFE
605827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
606827bd09bSSatish Balay   if (!vals||!work||!oprs||!segs)
60777431f27SBarry Smith     {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
608827bd09bSSatish Balay 
609827bd09bSSatish Balay   /* non-uniform should have at least two entries */
610827bd09bSSatish Balay 
611827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
612827bd09bSSatish Balay   if (!p_init)
613827bd09bSSatish Balay     {comm_init();}
614827bd09bSSatish Balay #endif
615827bd09bSSatish Balay 
616827bd09bSSatish Balay   /* if there's nothing to do return */
617827bd09bSSatish Balay   if ((num_nodes<2)||(dim<=0))
618827bd09bSSatish Balay     {return;}
619827bd09bSSatish Balay 
620827bd09bSSatish Balay   /* the error msg says it all!!! */
621827bd09bSSatish Balay   if (modfl_num_nodes)
622827bd09bSSatish Balay     {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");}
623827bd09bSSatish Balay 
624827bd09bSSatish Balay   /* can't do more dimensions then exist */
62539945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
626827bd09bSSatish Balay 
627827bd09bSSatish Balay   /* advance to list of n operations for custom */
628827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
629827bd09bSSatish Balay     {oprs++;}
630827bd09bSSatish Balay 
631d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))){
632827bd09bSSatish Balay     error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n");
633827bd09bSSatish Balay     fp = (vfp) oprs;
634827bd09bSSatish Balay   }
635827bd09bSSatish Balay 
636827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
637827bd09bSSatish Balay     {
638827bd09bSSatish Balay       n = segs[dim]-segs[edge];
639827bd09bSSatish Balay       dest = my_id^mask;
640827bd09bSSatish Balay       if (my_id > dest)
641a501084fSBarry Smith 	{MPI_Send(vals+segs[edge],n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
642827bd09bSSatish Balay       else
643827bd09bSSatish Balay 	{
644a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,
645827bd09bSSatish Balay 		   MPI_COMM_WORLD, &status);
646827bd09bSSatish Balay 	  (*fp)(vals+segs[edge], work, n, oprs);
647827bd09bSSatish Balay 	}
648827bd09bSSatish Balay     }
649827bd09bSSatish Balay 
650827bd09bSSatish Balay   if (edge==dim)
651827bd09bSSatish Balay     {mask>>=1;}
652827bd09bSSatish Balay   else
653827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
654827bd09bSSatish Balay 
655827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
656827bd09bSSatish Balay     {
657827bd09bSSatish Balay       if (my_id%mask)
658827bd09bSSatish Balay 	{continue;}
659827bd09bSSatish Balay 
660827bd09bSSatish Balay       n = (segs[dim]-segs[dim-1-edge]);
661827bd09bSSatish Balay 
662827bd09bSSatish Balay       dest = my_id^mask;
663827bd09bSSatish Balay       if (my_id < dest)
664a501084fSBarry Smith 	{MPI_Send(vals+segs[dim-1-edge],n,MPIU_SCALAR,dest,MSGTAG4+my_id,
665827bd09bSSatish Balay 		  MPI_COMM_WORLD);}
666827bd09bSSatish Balay       else
667827bd09bSSatish Balay 	{
668a501084fSBarry Smith 	  MPI_Recv(vals+segs[dim-1-edge],n,MPIU_SCALAR,MPI_ANY_SOURCE,
669827bd09bSSatish Balay 		   MSGTAG4+dest,MPI_COMM_WORLD, &status);
670827bd09bSSatish Balay 	}
671827bd09bSSatish Balay     }
672827bd09bSSatish Balay }
673827bd09bSSatish Balay 
674827bd09bSSatish Balay /******************************************************************************
675827bd09bSSatish Balay Function: giop()
676827bd09bSSatish Balay 
677827bd09bSSatish Balay Input :
678827bd09bSSatish Balay Output:
679827bd09bSSatish Balay Return:
680827bd09bSSatish Balay Description:
681827bd09bSSatish Balay 
682827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii
683827bd09bSSatish Balay 
684827bd09bSSatish Balay ******************************************************************************/
685a501084fSBarry Smith void new_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  int level, int *segs)
686827bd09bSSatish Balay {
687a501084fSBarry Smith    int edge, type, dest, mask;
688a501084fSBarry Smith    int stage_n;
689827bd09bSSatish Balay   MPI_Status  status;
690827bd09bSSatish Balay 
691827bd09bSSatish Balay #ifdef SAFE
692827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
693827bd09bSSatish Balay   if (!p_init)
694827bd09bSSatish Balay     {comm_init();}
695827bd09bSSatish Balay #endif
696827bd09bSSatish Balay 
697827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
698827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
699827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
700827bd09bSSatish Balay     {
701827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
702827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
703827bd09bSSatish Balay 	{
704827bd09bSSatish Balay 	  dest = edge_node[edge];
705827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
706827bd09bSSatish Balay 	  if (my_id>dest)
707a501084fSBarry Smith           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type,
708c700b191SBarry Smith                       MPI_COMM_WORLD);}
709827bd09bSSatish Balay 	  else
710827bd09bSSatish Balay 	    {
711827bd09bSSatish Balay 	      type =  type - my_id + dest;
712a501084fSBarry Smith               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,
713c700b191SBarry Smith                        MPI_COMM_WORLD,&status);
714827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
715827bd09bSSatish Balay 	    }
716827bd09bSSatish Balay 	}
717827bd09bSSatish Balay       mask <<= 1;
718827bd09bSSatish Balay     }
719827bd09bSSatish Balay   mask>>=1;
720827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
721827bd09bSSatish Balay     {
722827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
723827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
724827bd09bSSatish Balay 	{
725827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
726827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
727827bd09bSSatish Balay 	  if (my_id<dest)
728a501084fSBarry Smith             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,
729c700b191SBarry Smith                       MPI_COMM_WORLD);}
730827bd09bSSatish Balay 	  else
731827bd09bSSatish Balay 	    {
732827bd09bSSatish Balay 	      type =  type - my_id + dest;
733a501084fSBarry Smith               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,
734c700b191SBarry Smith                        MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
735827bd09bSSatish Balay 	    }
736827bd09bSSatish Balay 	}
737827bd09bSSatish Balay       mask >>= 1;
738827bd09bSSatish Balay     }
739827bd09bSSatish Balay }
740827bd09bSSatish Balay 
741827bd09bSSatish Balay 
742827bd09bSSatish Balay 
743827bd09bSSatish Balay /***********************************comm.c*************************************
744827bd09bSSatish Balay Function: giop()
745827bd09bSSatish Balay 
746827bd09bSSatish Balay Input :
747827bd09bSSatish Balay Output:
748827bd09bSSatish Balay Return:
749827bd09bSSatish Balay Description: fan-in/out version
750827bd09bSSatish Balay 
751827bd09bSSatish Balay note good only for num_nodes=2^k!!!
752827bd09bSSatish Balay 
753827bd09bSSatish Balay ***********************************comm.c*************************************/
754c700b191SBarry Smith void giop_hc(int *vals, int *work, int n, int *oprs, int dim)
755827bd09bSSatish Balay {
756a501084fSBarry Smith    int mask, edge;
757827bd09bSSatish Balay   int type, dest;
758827bd09bSSatish Balay   vfp fp;
759827bd09bSSatish Balay   MPI_Status  status;
760827bd09bSSatish Balay 
761827bd09bSSatish Balay #ifdef SAFE
762827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
763827bd09bSSatish Balay   if (!vals||!work||!oprs)
76477431f27SBarry Smith     {error_msg_fatal("giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
765827bd09bSSatish Balay 
766827bd09bSSatish Balay   /* non-uniform should have at least two entries */
767827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
768827bd09bSSatish Balay     {error_msg_fatal("giop_hc() :: non_uniform and n=0,1?");}
769827bd09bSSatish Balay 
770827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
771827bd09bSSatish Balay   if (!p_init)
772827bd09bSSatish Balay     {comm_init();}
773827bd09bSSatish Balay #endif
774827bd09bSSatish Balay 
775827bd09bSSatish Balay   /* if there's nothing to do return */
776827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
777827bd09bSSatish Balay     {return;}
778827bd09bSSatish Balay 
779827bd09bSSatish Balay   /* the error msg says it all!!! */
780827bd09bSSatish Balay   if (modfl_num_nodes)
781827bd09bSSatish Balay     {error_msg_fatal("giop_hc() :: num_nodes not a power of 2!?!");}
782827bd09bSSatish Balay 
783827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
784827bd09bSSatish Balay   if (n<0)
78577431f27SBarry Smith     {error_msg_fatal("giop_hc() :: n=%D<0?",n);}
786827bd09bSSatish Balay 
787827bd09bSSatish Balay   /* can't do more dimensions then exist */
78839945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
789827bd09bSSatish Balay 
790827bd09bSSatish Balay   /* advance to list of n operations for custom */
791827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
792827bd09bSSatish Balay     {oprs++;}
793827bd09bSSatish Balay 
794d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))){
795827bd09bSSatish Balay     error_msg_warning("giop_hc() :: hope you passed in a rbfp!\n");
796827bd09bSSatish Balay     fp = (vfp) oprs;
797827bd09bSSatish Balay   }
798827bd09bSSatish Balay 
799827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
800827bd09bSSatish Balay     {
801827bd09bSSatish Balay       dest = my_id^mask;
802827bd09bSSatish Balay       if (my_id > dest)
803a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);}
804827bd09bSSatish Balay       else
805827bd09bSSatish Balay 	{
806a501084fSBarry Smith 	  MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,
807827bd09bSSatish Balay 		   &status);
808827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
809827bd09bSSatish Balay 	}
810827bd09bSSatish Balay     }
811827bd09bSSatish Balay 
812827bd09bSSatish Balay   if (edge==dim)
813827bd09bSSatish Balay     {mask>>=1;}
814827bd09bSSatish Balay   else
815827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
816827bd09bSSatish Balay 
817827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
818827bd09bSSatish Balay     {
819827bd09bSSatish Balay       if (my_id%mask)
820827bd09bSSatish Balay 	{continue;}
821827bd09bSSatish Balay 
822827bd09bSSatish Balay       dest = my_id^mask;
823827bd09bSSatish Balay       if (my_id < dest)
824a501084fSBarry Smith 	{MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);}
825827bd09bSSatish Balay       else
826827bd09bSSatish Balay 	{
827a501084fSBarry Smith 	  MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,
828827bd09bSSatish Balay 		   &status);
829827bd09bSSatish Balay 	}
830827bd09bSSatish Balay     }
831827bd09bSSatish Balay }
832