xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision e7e72b3d0edcd0d15e7f68c03be08666507fc872)
1dba47a55SKris 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*************************************/
177c4f633dSBarry Smith #include "../src/ksp/pc/impls/tfs/tfs.h"
18827bd09bSSatish Balay 
19827bd09bSSatish Balay 
20827bd09bSSatish Balay /* global program control variables - explicitly exported */
2152f87cdaSBarry Smith PetscMPIInt my_id            = 0;
2252f87cdaSBarry Smith PetscMPIInt num_nodes        = 1;
236e4f4d19SBarry Smith PetscMPIInt floor_num_nodes  = 0;
246e4f4d19SBarry Smith PetscMPIInt i_log2_num_nodes = 0;
25827bd09bSSatish Balay 
26827bd09bSSatish Balay /* global program control variables */
2752f87cdaSBarry Smith static PetscInt p_init = 0;
2852f87cdaSBarry Smith static PetscInt modfl_num_nodes;
2952f87cdaSBarry Smith static PetscInt edge_not_pow_2;
30827bd09bSSatish Balay 
3152f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt)*32];
32827bd09bSSatish Balay 
337b1ae94cSBarry Smith /***********************************comm.c*************************************/
340924e98cSBarry Smith PetscErrorCode comm_init (void)
35827bd09bSSatish Balay {
36827bd09bSSatish Balay 
373fdc5746SBarry Smith   if (p_init++)   PetscFunctionReturn(0);
38827bd09bSSatish Balay 
39827bd09bSSatish Balay   MPI_Comm_size(MPI_COMM_WORLD,&num_nodes);
40827bd09bSSatish Balay   MPI_Comm_rank(MPI_COMM_WORLD,&my_id);
41827bd09bSSatish Balay 
42*e7e72b3dSBarry Smith   if (num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!");
43827bd09bSSatish Balay 
443fdc5746SBarry Smith   ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32);
45827bd09bSSatish Balay 
46827bd09bSSatish Balay   floor_num_nodes = 1;
47827bd09bSSatish Balay   i_log2_num_nodes = modfl_num_nodes = 0;
48827bd09bSSatish Balay   while (floor_num_nodes <= num_nodes)
49827bd09bSSatish Balay     {
50827bd09bSSatish Balay       edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes;
51827bd09bSSatish Balay       floor_num_nodes <<= 1;
52827bd09bSSatish Balay       i_log2_num_nodes++;
53827bd09bSSatish Balay     }
54827bd09bSSatish Balay 
55827bd09bSSatish Balay   i_log2_num_nodes--;
56827bd09bSSatish Balay   floor_num_nodes >>= 1;
57827bd09bSSatish Balay   modfl_num_nodes = (num_nodes - floor_num_nodes);
58827bd09bSSatish Balay 
59827bd09bSSatish Balay   if ((my_id > 0) && (my_id <= modfl_num_nodes))
60827bd09bSSatish Balay     {edge_not_pow_2=((my_id|floor_num_nodes)-1);}
61827bd09bSSatish Balay   else if (my_id >= floor_num_nodes)
62827bd09bSSatish Balay     {edge_not_pow_2=((my_id^floor_num_nodes)+1);
63827bd09bSSatish Balay     }
64827bd09bSSatish Balay   else
65827bd09bSSatish Balay     {edge_not_pow_2 = 0;}
663fdc5746SBarry Smith   PetscFunctionReturn(0);
67827bd09bSSatish Balay }
68827bd09bSSatish Balay 
697b1ae94cSBarry Smith /***********************************comm.c*************************************/
700924e98cSBarry Smith PetscErrorCode giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
71827bd09bSSatish Balay {
723fdc5746SBarry Smith   PetscInt   mask, edge;
733fdc5746SBarry Smith   PetscInt    type, dest;
74827bd09bSSatish Balay   vfp         fp;
75827bd09bSSatish Balay   MPI_Status  status;
763fdc5746SBarry Smith   PetscInt    ierr;
77827bd09bSSatish Balay 
783fdc5746SBarry Smith    PetscFunctionBegin;
79827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
80827bd09bSSatish Balay   if (!vals||!work||!oprs)
81e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
82827bd09bSSatish Balay 
83827bd09bSSatish Balay   /* non-uniform should have at least two entries */
84*e7e72b3dSBarry Smith   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: non_uniform and n=0,1?");
85827bd09bSSatish Balay 
86827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
87827bd09bSSatish Balay   if (!p_init)
88827bd09bSSatish Balay     {comm_init();}
89827bd09bSSatish Balay 
90827bd09bSSatish Balay   /* if there's nothing to do return */
91827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
92827bd09bSSatish Balay     {
933fdc5746SBarry Smith         PetscFunctionReturn(0);
94827bd09bSSatish Balay     }
95827bd09bSSatish Balay 
9671a0148aSBarry Smith 
97827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
98827bd09bSSatish Balay   if (n<0)
99e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: n=%D<0?",n);}
100827bd09bSSatish Balay 
101827bd09bSSatish Balay   /* advance to list of n operations for custom */
102827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
103827bd09bSSatish Balay     {oprs++;}
104827bd09bSSatish Balay 
105827bd09bSSatish Balay   /* major league hack */
106d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))) {
107f1ed62a8SBarry Smith     ierr = PetscInfo(0,"giop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
108827bd09bSSatish Balay     fp = (vfp) oprs;
109827bd09bSSatish Balay   }
110827bd09bSSatish Balay 
111827bd09bSSatish Balay   /* all msgs will be of the same length */
112827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
113827bd09bSSatish Balay   if (edge_not_pow_2)
114827bd09bSSatish Balay     {
115827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
1163fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
117827bd09bSSatish Balay       else
118827bd09bSSatish Balay 	{
1193fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
120827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
121827bd09bSSatish Balay 	}
122827bd09bSSatish Balay     }
123827bd09bSSatish Balay 
124827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
125827bd09bSSatish Balay   if (my_id<floor_num_nodes)
126827bd09bSSatish Balay     {
127827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
128827bd09bSSatish Balay 	{
129827bd09bSSatish Balay 	  dest = my_id^mask;
130827bd09bSSatish Balay 	  if (my_id > dest)
1313fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
132827bd09bSSatish Balay 	  else
133827bd09bSSatish Balay 	    {
1343fdc5746SBarry Smith 	      ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
135827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
136827bd09bSSatish Balay 	    }
137827bd09bSSatish Balay 	}
138827bd09bSSatish Balay 
139827bd09bSSatish Balay       mask=floor_num_nodes>>1;
140827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
141827bd09bSSatish Balay 	{
142827bd09bSSatish Balay 	  if (my_id%mask)
143827bd09bSSatish Balay 	    {continue;}
144827bd09bSSatish Balay 
145827bd09bSSatish Balay 	  dest = my_id^mask;
146827bd09bSSatish Balay 	  if (my_id < dest)
1473fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
148827bd09bSSatish Balay 	  else
149827bd09bSSatish Balay 	    {
1503fdc5746SBarry Smith 	      ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
151827bd09bSSatish Balay 	    }
152827bd09bSSatish Balay 	}
153827bd09bSSatish Balay     }
154827bd09bSSatish Balay 
155827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
156827bd09bSSatish Balay   if (edge_not_pow_2)
157827bd09bSSatish Balay     {
158827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
159827bd09bSSatish Balay 	{
1603fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
161827bd09bSSatish Balay 	}
162827bd09bSSatish Balay       else
1633fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
164827bd09bSSatish Balay     }
1653fdc5746SBarry Smith         PetscFunctionReturn(0);
166827bd09bSSatish Balay }
167827bd09bSSatish Balay 
1687b1ae94cSBarry Smith /***********************************comm.c*************************************/
1696e4f4d19SBarry Smith PetscErrorCode grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
170827bd09bSSatish Balay {
1713fdc5746SBarry Smith   PetscInt       mask, edge;
1723fdc5746SBarry Smith   PetscInt       type, dest;
173827bd09bSSatish Balay   vfp            fp;
174827bd09bSSatish Balay   MPI_Status     status;
1753fdc5746SBarry Smith   PetscErrorCode ierr;
176827bd09bSSatish Balay 
1773fdc5746SBarry Smith    PetscFunctionBegin;
178827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
179827bd09bSSatish Balay   if (!vals||!work||!oprs)
180e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
181827bd09bSSatish Balay 
182827bd09bSSatish Balay   /* non-uniform should have at least two entries */
183*e7e72b3dSBarry Smith   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop() :: non_uniform and n=0,1?");
184827bd09bSSatish Balay 
185827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
186827bd09bSSatish Balay   if (!p_init)
187827bd09bSSatish Balay     {comm_init();}
188827bd09bSSatish Balay 
189827bd09bSSatish Balay   /* if there's nothing to do return */
190827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
1913fdc5746SBarry Smith     {        PetscFunctionReturn(0);}
192827bd09bSSatish Balay 
193827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
194827bd09bSSatish Balay   if (n<0)
195e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);}
196827bd09bSSatish Balay 
197827bd09bSSatish Balay   /* advance to list of n operations for custom */
198827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
199827bd09bSSatish Balay     {oprs++;}
200827bd09bSSatish Balay 
201d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
202f1ed62a8SBarry Smith     ierr = PetscInfo(0,"grop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
203827bd09bSSatish Balay     fp = (vfp) oprs;
204827bd09bSSatish Balay   }
205827bd09bSSatish Balay 
206827bd09bSSatish Balay   /* all msgs will be of the same length */
207827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
208827bd09bSSatish Balay   if (edge_not_pow_2)
209827bd09bSSatish Balay     {
210827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
2113fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
212827bd09bSSatish Balay       else
213827bd09bSSatish Balay 	{
2143fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
215827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
216827bd09bSSatish Balay 	}
217827bd09bSSatish Balay     }
218827bd09bSSatish Balay 
219827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
220827bd09bSSatish Balay   if (my_id<floor_num_nodes)
221827bd09bSSatish Balay     {
222827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
223827bd09bSSatish Balay 	{
224827bd09bSSatish Balay 	  dest = my_id^mask;
225827bd09bSSatish Balay 	  if (my_id > dest)
2263fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
227827bd09bSSatish Balay 	  else
228827bd09bSSatish Balay 	    {
2293fdc5746SBarry Smith 	      ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
230827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
231827bd09bSSatish Balay 	    }
232827bd09bSSatish Balay 	}
233827bd09bSSatish Balay 
234827bd09bSSatish Balay       mask=floor_num_nodes>>1;
235827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
236827bd09bSSatish Balay 	{
237827bd09bSSatish Balay 	  if (my_id%mask)
238827bd09bSSatish Balay 	    {continue;}
239827bd09bSSatish Balay 
240827bd09bSSatish Balay 	  dest = my_id^mask;
241827bd09bSSatish Balay 	  if (my_id < dest)
2423fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
243827bd09bSSatish Balay 	  else
244827bd09bSSatish Balay 	    {
2453fdc5746SBarry Smith 	      ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
246827bd09bSSatish Balay 	    }
247827bd09bSSatish Balay 	}
248827bd09bSSatish Balay     }
249827bd09bSSatish Balay 
250827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
251827bd09bSSatish Balay   if (edge_not_pow_2)
252827bd09bSSatish Balay     {
253827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
254827bd09bSSatish Balay 	{
2553fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
256827bd09bSSatish Balay 	}
257827bd09bSSatish Balay       else
2583fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
259827bd09bSSatish Balay     }
2603fdc5746SBarry Smith         PetscFunctionReturn(0);
261827bd09bSSatish Balay }
262827bd09bSSatish Balay 
2637b1ae94cSBarry Smith /***********************************comm.c*************************************/
26452f87cdaSBarry Smith PetscErrorCode grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
265827bd09bSSatish Balay {
2663fdc5746SBarry Smith   PetscInt       mask, edge;
2673fdc5746SBarry Smith   PetscInt       type, dest;
268827bd09bSSatish Balay   vfp            fp;
269827bd09bSSatish Balay   MPI_Status     status;
2703fdc5746SBarry Smith   PetscErrorCode ierr;
271827bd09bSSatish Balay 
2723fdc5746SBarry Smith    PetscFunctionBegin;
273827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
274827bd09bSSatish Balay   if (!vals||!work||!oprs)
275e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
276827bd09bSSatish Balay 
277827bd09bSSatish Balay   /* non-uniform should have at least two entries */
278*e7e72b3dSBarry Smith   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: non_uniform and n=0,1?");
279827bd09bSSatish Balay 
280827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
281827bd09bSSatish Balay   if (!p_init)
282827bd09bSSatish Balay     {comm_init();}
283827bd09bSSatish Balay 
284827bd09bSSatish Balay   /* if there's nothing to do return */
285827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
2860924e98cSBarry Smith     {PetscFunctionReturn(0);}
287827bd09bSSatish Balay 
288827bd09bSSatish Balay   /* the error msg says it all!!! */
289*e7e72b3dSBarry Smith   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: num_nodes not a power of 2!?!");
290827bd09bSSatish Balay 
291827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
292827bd09bSSatish Balay   if (n<0)
293e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: n=%D<0?",n);}
294827bd09bSSatish Balay 
295827bd09bSSatish Balay   /* can't do more dimensions then exist */
29639945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
297827bd09bSSatish Balay 
298827bd09bSSatish Balay   /* advance to list of n operations for custom */
299827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
300827bd09bSSatish Balay     {oprs++;}
301827bd09bSSatish Balay 
302d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
303f1ed62a8SBarry Smith     ierr = PetscInfo(0,"grop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
304827bd09bSSatish Balay     fp = (vfp) oprs;
305827bd09bSSatish Balay   }
306827bd09bSSatish Balay 
307827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
308827bd09bSSatish Balay     {
309827bd09bSSatish Balay       dest = my_id^mask;
310827bd09bSSatish Balay       if (my_id > dest)
3113fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
312827bd09bSSatish Balay       else
313827bd09bSSatish Balay 	{
3143fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
315827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
316827bd09bSSatish Balay 	}
317827bd09bSSatish Balay     }
318827bd09bSSatish Balay 
319827bd09bSSatish Balay   if (edge==dim)
320827bd09bSSatish Balay     {mask>>=1;}
321827bd09bSSatish Balay   else
322827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
323827bd09bSSatish Balay 
324827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
325827bd09bSSatish Balay     {
326827bd09bSSatish Balay       if (my_id%mask)
327827bd09bSSatish Balay 	{continue;}
328827bd09bSSatish Balay 
329827bd09bSSatish Balay       dest = my_id^mask;
330827bd09bSSatish Balay       if (my_id < dest)
3313fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
332827bd09bSSatish Balay       else
333827bd09bSSatish Balay 	{
3343fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
335827bd09bSSatish Balay 	}
336827bd09bSSatish Balay     }
3373fdc5746SBarry Smith         PetscFunctionReturn(0);
338827bd09bSSatish Balay }
339827bd09bSSatish Balay 
3407b1ae94cSBarry Smith /******************************************************************************/
3410924e98cSBarry Smith PetscErrorCode ssgl_radd( PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
342827bd09bSSatish Balay {
3433fdc5746SBarry Smith   PetscInt       edge, type, dest, mask;
3443fdc5746SBarry Smith   PetscInt       stage_n;
345827bd09bSSatish Balay   MPI_Status     status;
3463fdc5746SBarry Smith   PetscErrorCode ierr;
347827bd09bSSatish Balay 
3483fdc5746SBarry Smith    PetscFunctionBegin;
349827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
350827bd09bSSatish Balay   if (!p_init)
351827bd09bSSatish Balay     {comm_init();}
352827bd09bSSatish Balay 
353827bd09bSSatish Balay 
354827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
355827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
356827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
357827bd09bSSatish Balay     {
358827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
359827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
360827bd09bSSatish Balay 	{
361827bd09bSSatish Balay 	  dest = edge_node[edge];
362827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
363827bd09bSSatish Balay 	  if (my_id>dest)
3643fdc5746SBarry Smith           {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);}
365827bd09bSSatish Balay 	  else
366827bd09bSSatish Balay 	    {
367827bd09bSSatish Balay 	      type =  type - my_id + dest;
3683fdc5746SBarry Smith               ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
369827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
370827bd09bSSatish Balay 	    }
371827bd09bSSatish Balay 	}
372827bd09bSSatish Balay       mask <<= 1;
373827bd09bSSatish Balay     }
374827bd09bSSatish Balay   mask>>=1;
375827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
376827bd09bSSatish Balay     {
377827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
378827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
379827bd09bSSatish Balay 	{
380827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
381827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
382827bd09bSSatish Balay 	  if (my_id<dest)
3833fdc5746SBarry Smith             {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);}
384827bd09bSSatish Balay 	  else
385827bd09bSSatish Balay 	    {
386827bd09bSSatish Balay 	      type =  type - my_id + dest;
3873fdc5746SBarry Smith               ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
388827bd09bSSatish Balay 	    }
389827bd09bSSatish Balay 	}
390827bd09bSSatish Balay       mask >>= 1;
391827bd09bSSatish Balay     }
3923fdc5746SBarry Smith   PetscFunctionReturn(0);
393827bd09bSSatish Balay }
394827bd09bSSatish Balay 
3957b1ae94cSBarry Smith /******************************************************************************/
39652f87cdaSBarry Smith PetscErrorCode new_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
397827bd09bSSatish Balay {
39852f87cdaSBarry Smith   PetscInt            edge, type, dest, mask;
39952f87cdaSBarry Smith   PetscInt            stage_n;
400827bd09bSSatish Balay   MPI_Status     status;
4013fdc5746SBarry Smith   PetscErrorCode ierr;
402827bd09bSSatish Balay 
4033fdc5746SBarry Smith    PetscFunctionBegin;
404827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
405827bd09bSSatish Balay   if (!p_init)
406827bd09bSSatish Balay     {comm_init();}
407827bd09bSSatish Balay 
408827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
409827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
410827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
411827bd09bSSatish Balay     {
412827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
413827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
414827bd09bSSatish Balay 	{
415827bd09bSSatish Balay 	  dest = edge_node[edge];
416827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
417827bd09bSSatish Balay 	  if (my_id>dest)
4183fdc5746SBarry Smith           {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);}
419827bd09bSSatish Balay 	  else
420827bd09bSSatish Balay 	    {
421827bd09bSSatish Balay 	      type =  type - my_id + dest;
4223fdc5746SBarry Smith               ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
423827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
424827bd09bSSatish Balay 	    }
425827bd09bSSatish Balay 	}
426827bd09bSSatish Balay       mask <<= 1;
427827bd09bSSatish Balay     }
428827bd09bSSatish Balay   mask>>=1;
429827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
430827bd09bSSatish Balay     {
431827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
432827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
433827bd09bSSatish Balay 	{
434827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
435827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
436827bd09bSSatish Balay 	  if (my_id<dest)
4373fdc5746SBarry Smith             {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);}
438827bd09bSSatish Balay 	  else
439827bd09bSSatish Balay 	    {
440827bd09bSSatish Balay 	      type =  type - my_id + dest;
4413fdc5746SBarry Smith               ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
442827bd09bSSatish Balay 	    }
443827bd09bSSatish Balay 	}
444827bd09bSSatish Balay       mask >>= 1;
445827bd09bSSatish Balay     }
4463fdc5746SBarry Smith   PetscFunctionReturn(0);
447827bd09bSSatish Balay }
448827bd09bSSatish Balay 
4497b1ae94cSBarry Smith /***********************************comm.c*************************************/
45052f87cdaSBarry Smith PetscErrorCode giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
451827bd09bSSatish Balay {
45252f87cdaSBarry Smith   PetscInt            mask, edge;
45352f87cdaSBarry Smith   PetscInt            type, dest;
454827bd09bSSatish Balay   vfp            fp;
455827bd09bSSatish Balay   MPI_Status     status;
4563fdc5746SBarry Smith   PetscErrorCode ierr;
457827bd09bSSatish Balay 
4583fdc5746SBarry Smith    PetscFunctionBegin;
459827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
460827bd09bSSatish Balay   if (!vals||!work||!oprs)
461e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
462827bd09bSSatish Balay 
463827bd09bSSatish Balay   /* non-uniform should have at least two entries */
464*e7e72b3dSBarry Smith   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: non_uniform and n=0,1?");
465827bd09bSSatish Balay 
466827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
467827bd09bSSatish Balay   if (!p_init)
468827bd09bSSatish Balay     {comm_init();}
469827bd09bSSatish Balay 
470827bd09bSSatish Balay   /* if there's nothing to do return */
471827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
4723fdc5746SBarry Smith     {  PetscFunctionReturn(0);}
473827bd09bSSatish Balay 
474827bd09bSSatish Balay   /* the error msg says it all!!! */
475*e7e72b3dSBarry Smith   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: num_nodes not a power of 2!?!");
476827bd09bSSatish Balay 
477827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
478827bd09bSSatish Balay   if (n<0)
479e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: n=%D<0?",n);}
480827bd09bSSatish Balay 
481827bd09bSSatish Balay   /* can't do more dimensions then exist */
48239945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
483827bd09bSSatish Balay 
484827bd09bSSatish Balay   /* advance to list of n operations for custom */
485827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
486827bd09bSSatish Balay     {oprs++;}
487827bd09bSSatish Balay 
488d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))){
489f1ed62a8SBarry Smith     ierr = PetscInfo(0,"giop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
490827bd09bSSatish Balay     fp = (vfp) oprs;
491827bd09bSSatish Balay   }
492827bd09bSSatish Balay 
493827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
494827bd09bSSatish Balay     {
495827bd09bSSatish Balay       dest = my_id^mask;
496827bd09bSSatish Balay       if (my_id > dest)
4973fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
498827bd09bSSatish Balay       else
499827bd09bSSatish Balay 	{
5003fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
501827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
502827bd09bSSatish Balay 	}
503827bd09bSSatish Balay     }
504827bd09bSSatish Balay 
505827bd09bSSatish Balay   if (edge==dim)
506827bd09bSSatish Balay     {mask>>=1;}
507827bd09bSSatish Balay   else
508827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
509827bd09bSSatish Balay 
510827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
511827bd09bSSatish Balay     {
512827bd09bSSatish Balay       if (my_id%mask)
513827bd09bSSatish Balay 	{continue;}
514827bd09bSSatish Balay 
515827bd09bSSatish Balay       dest = my_id^mask;
516827bd09bSSatish Balay       if (my_id < dest)
5173fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
518827bd09bSSatish Balay       else
519827bd09bSSatish Balay 	{
5203fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
521827bd09bSSatish Balay 	}
522827bd09bSSatish Balay     }
5233fdc5746SBarry Smith   PetscFunctionReturn(0);
524827bd09bSSatish Balay }
525