xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
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 
42827bd09bSSatish Balay   if (num_nodes> (INT_MAX >> 1))
43*e32f2f54SBarry Smith   {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!");}
44827bd09bSSatish Balay 
453fdc5746SBarry Smith   ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32);
46827bd09bSSatish Balay 
47827bd09bSSatish Balay   floor_num_nodes = 1;
48827bd09bSSatish Balay   i_log2_num_nodes = modfl_num_nodes = 0;
49827bd09bSSatish Balay   while (floor_num_nodes <= num_nodes)
50827bd09bSSatish Balay     {
51827bd09bSSatish Balay       edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes;
52827bd09bSSatish Balay       floor_num_nodes <<= 1;
53827bd09bSSatish Balay       i_log2_num_nodes++;
54827bd09bSSatish Balay     }
55827bd09bSSatish Balay 
56827bd09bSSatish Balay   i_log2_num_nodes--;
57827bd09bSSatish Balay   floor_num_nodes >>= 1;
58827bd09bSSatish Balay   modfl_num_nodes = (num_nodes - floor_num_nodes);
59827bd09bSSatish Balay 
60827bd09bSSatish Balay   if ((my_id > 0) && (my_id <= modfl_num_nodes))
61827bd09bSSatish Balay     {edge_not_pow_2=((my_id|floor_num_nodes)-1);}
62827bd09bSSatish Balay   else if (my_id >= floor_num_nodes)
63827bd09bSSatish Balay     {edge_not_pow_2=((my_id^floor_num_nodes)+1);
64827bd09bSSatish Balay     }
65827bd09bSSatish Balay   else
66827bd09bSSatish Balay     {edge_not_pow_2 = 0;}
673fdc5746SBarry Smith   PetscFunctionReturn(0);
68827bd09bSSatish Balay }
69827bd09bSSatish Balay 
707b1ae94cSBarry Smith /***********************************comm.c*************************************/
710924e98cSBarry Smith PetscErrorCode giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
72827bd09bSSatish Balay {
733fdc5746SBarry Smith   PetscInt   mask, edge;
743fdc5746SBarry Smith   PetscInt    type, dest;
75827bd09bSSatish Balay   vfp         fp;
76827bd09bSSatish Balay   MPI_Status  status;
773fdc5746SBarry Smith   PetscInt    ierr;
78827bd09bSSatish Balay 
793fdc5746SBarry Smith    PetscFunctionBegin;
80827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
81827bd09bSSatish Balay   if (!vals||!work||!oprs)
82*e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
83827bd09bSSatish Balay 
84827bd09bSSatish Balay   /* non-uniform should have at least two entries */
85827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
86*e32f2f54SBarry Smith     {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: non_uniform and n=0,1?");}
87827bd09bSSatish Balay 
88827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
89827bd09bSSatish Balay   if (!p_init)
90827bd09bSSatish Balay     {comm_init();}
91827bd09bSSatish Balay 
92827bd09bSSatish Balay   /* if there's nothing to do return */
93827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
94827bd09bSSatish Balay     {
953fdc5746SBarry Smith         PetscFunctionReturn(0);
96827bd09bSSatish Balay     }
97827bd09bSSatish Balay 
9871a0148aSBarry Smith 
99827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
100827bd09bSSatish Balay   if (n<0)
101*e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: n=%D<0?",n);}
102827bd09bSSatish Balay 
103827bd09bSSatish Balay   /* advance to list of n operations for custom */
104827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
105827bd09bSSatish Balay     {oprs++;}
106827bd09bSSatish Balay 
107827bd09bSSatish Balay   /* major league hack */
108d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))) {
109f1ed62a8SBarry Smith     ierr = PetscInfo(0,"giop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
110827bd09bSSatish Balay     fp = (vfp) oprs;
111827bd09bSSatish Balay   }
112827bd09bSSatish Balay 
113827bd09bSSatish Balay   /* all msgs will be of the same length */
114827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
115827bd09bSSatish Balay   if (edge_not_pow_2)
116827bd09bSSatish Balay     {
117827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
1183fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
119827bd09bSSatish Balay       else
120827bd09bSSatish Balay 	{
1213fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
122827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
123827bd09bSSatish Balay 	}
124827bd09bSSatish Balay     }
125827bd09bSSatish Balay 
126827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
127827bd09bSSatish Balay   if (my_id<floor_num_nodes)
128827bd09bSSatish Balay     {
129827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
130827bd09bSSatish Balay 	{
131827bd09bSSatish Balay 	  dest = my_id^mask;
132827bd09bSSatish Balay 	  if (my_id > dest)
1333fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
134827bd09bSSatish Balay 	  else
135827bd09bSSatish Balay 	    {
1363fdc5746SBarry Smith 	      ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
137827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
138827bd09bSSatish Balay 	    }
139827bd09bSSatish Balay 	}
140827bd09bSSatish Balay 
141827bd09bSSatish Balay       mask=floor_num_nodes>>1;
142827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
143827bd09bSSatish Balay 	{
144827bd09bSSatish Balay 	  if (my_id%mask)
145827bd09bSSatish Balay 	    {continue;}
146827bd09bSSatish Balay 
147827bd09bSSatish Balay 	  dest = my_id^mask;
148827bd09bSSatish Balay 	  if (my_id < dest)
1493fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
150827bd09bSSatish Balay 	  else
151827bd09bSSatish Balay 	    {
1523fdc5746SBarry Smith 	      ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
153827bd09bSSatish Balay 	    }
154827bd09bSSatish Balay 	}
155827bd09bSSatish Balay     }
156827bd09bSSatish Balay 
157827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
158827bd09bSSatish Balay   if (edge_not_pow_2)
159827bd09bSSatish Balay     {
160827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
161827bd09bSSatish Balay 	{
1623fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
163827bd09bSSatish Balay 	}
164827bd09bSSatish Balay       else
1653fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
166827bd09bSSatish Balay     }
1673fdc5746SBarry Smith         PetscFunctionReturn(0);
168827bd09bSSatish Balay }
169827bd09bSSatish Balay 
1707b1ae94cSBarry Smith /***********************************comm.c*************************************/
1716e4f4d19SBarry Smith PetscErrorCode grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
172827bd09bSSatish Balay {
1733fdc5746SBarry Smith   PetscInt       mask, edge;
1743fdc5746SBarry Smith   PetscInt       type, dest;
175827bd09bSSatish Balay   vfp            fp;
176827bd09bSSatish Balay   MPI_Status     status;
1773fdc5746SBarry Smith   PetscErrorCode ierr;
178827bd09bSSatish Balay 
1793fdc5746SBarry Smith    PetscFunctionBegin;
180827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
181827bd09bSSatish Balay   if (!vals||!work||!oprs)
182*e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
183827bd09bSSatish Balay 
184827bd09bSSatish Balay   /* non-uniform should have at least two entries */
185827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
186*e32f2f54SBarry Smith     {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop() :: non_uniform and n=0,1?");}
187827bd09bSSatish Balay 
188827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
189827bd09bSSatish Balay   if (!p_init)
190827bd09bSSatish Balay     {comm_init();}
191827bd09bSSatish Balay 
192827bd09bSSatish Balay   /* if there's nothing to do return */
193827bd09bSSatish Balay   if ((num_nodes<2)||(!n))
1943fdc5746SBarry Smith     {        PetscFunctionReturn(0);}
195827bd09bSSatish Balay 
196827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
197827bd09bSSatish Balay   if (n<0)
198*e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);}
199827bd09bSSatish Balay 
200827bd09bSSatish Balay   /* advance to list of n operations for custom */
201827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
202827bd09bSSatish Balay     {oprs++;}
203827bd09bSSatish Balay 
204d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
205f1ed62a8SBarry Smith     ierr = PetscInfo(0,"grop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
206827bd09bSSatish Balay     fp = (vfp) oprs;
207827bd09bSSatish Balay   }
208827bd09bSSatish Balay 
209827bd09bSSatish Balay   /* all msgs will be of the same length */
210827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
211827bd09bSSatish Balay   if (edge_not_pow_2)
212827bd09bSSatish Balay     {
213827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
2143fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
215827bd09bSSatish Balay       else
216827bd09bSSatish Balay 	{
2173fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
218827bd09bSSatish Balay 	  (*fp)(vals,work,n,oprs);
219827bd09bSSatish Balay 	}
220827bd09bSSatish Balay     }
221827bd09bSSatish Balay 
222827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
223827bd09bSSatish Balay   if (my_id<floor_num_nodes)
224827bd09bSSatish Balay     {
225827bd09bSSatish Balay       for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1)
226827bd09bSSatish Balay 	{
227827bd09bSSatish Balay 	  dest = my_id^mask;
228827bd09bSSatish Balay 	  if (my_id > dest)
2293fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
230827bd09bSSatish Balay 	  else
231827bd09bSSatish Balay 	    {
2323fdc5746SBarry Smith 	      ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
233827bd09bSSatish Balay 	      (*fp)(vals, work, n, oprs);
234827bd09bSSatish Balay 	    }
235827bd09bSSatish Balay 	}
236827bd09bSSatish Balay 
237827bd09bSSatish Balay       mask=floor_num_nodes>>1;
238827bd09bSSatish Balay       for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1)
239827bd09bSSatish Balay 	{
240827bd09bSSatish Balay 	  if (my_id%mask)
241827bd09bSSatish Balay 	    {continue;}
242827bd09bSSatish Balay 
243827bd09bSSatish Balay 	  dest = my_id^mask;
244827bd09bSSatish Balay 	  if (my_id < dest)
2453fdc5746SBarry Smith 	    {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
246827bd09bSSatish Balay 	  else
247827bd09bSSatish Balay 	    {
2483fdc5746SBarry Smith 	      ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
249827bd09bSSatish Balay 	    }
250827bd09bSSatish Balay 	}
251827bd09bSSatish Balay     }
252827bd09bSSatish Balay 
253827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
254827bd09bSSatish Balay   if (edge_not_pow_2)
255827bd09bSSatish Balay     {
256827bd09bSSatish Balay       if (my_id >= floor_num_nodes)
257827bd09bSSatish Balay 	{
2583fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
259827bd09bSSatish Balay 	}
260827bd09bSSatish Balay       else
2613fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
262827bd09bSSatish Balay     }
2633fdc5746SBarry Smith         PetscFunctionReturn(0);
264827bd09bSSatish Balay }
265827bd09bSSatish Balay 
2667b1ae94cSBarry Smith /***********************************comm.c*************************************/
26752f87cdaSBarry Smith PetscErrorCode grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
268827bd09bSSatish Balay {
2693fdc5746SBarry Smith   PetscInt       mask, edge;
2703fdc5746SBarry Smith   PetscInt       type, dest;
271827bd09bSSatish Balay   vfp            fp;
272827bd09bSSatish Balay   MPI_Status     status;
2733fdc5746SBarry Smith   PetscErrorCode ierr;
274827bd09bSSatish Balay 
2753fdc5746SBarry Smith    PetscFunctionBegin;
276827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
277827bd09bSSatish Balay   if (!vals||!work||!oprs)
278*e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
279827bd09bSSatish Balay 
280827bd09bSSatish Balay   /* non-uniform should have at least two entries */
281827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
282*e32f2f54SBarry Smith     {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: non_uniform and n=0,1?");}
283827bd09bSSatish Balay 
284827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
285827bd09bSSatish Balay   if (!p_init)
286827bd09bSSatish Balay     {comm_init();}
287827bd09bSSatish Balay 
288827bd09bSSatish Balay   /* if there's nothing to do return */
289827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
2900924e98cSBarry Smith     {PetscFunctionReturn(0);}
291827bd09bSSatish Balay 
292827bd09bSSatish Balay   /* the error msg says it all!!! */
293827bd09bSSatish Balay   if (modfl_num_nodes)
294*e32f2f54SBarry Smith     {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: num_nodes not a power of 2!?!");}
295827bd09bSSatish Balay 
296827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
297827bd09bSSatish Balay   if (n<0)
298*e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: n=%D<0?",n);}
299827bd09bSSatish Balay 
300827bd09bSSatish Balay   /* can't do more dimensions then exist */
30139945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
302827bd09bSSatish Balay 
303827bd09bSSatish Balay   /* advance to list of n operations for custom */
304827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
305827bd09bSSatish Balay     {oprs++;}
306827bd09bSSatish Balay 
307d890fc11SSatish Balay   if (!(fp = (vfp) rvec_fct_addr(type))) {
308f1ed62a8SBarry Smith     ierr = PetscInfo(0,"grop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
309827bd09bSSatish Balay     fp = (vfp) oprs;
310827bd09bSSatish Balay   }
311827bd09bSSatish Balay 
312827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
313827bd09bSSatish Balay     {
314827bd09bSSatish Balay       dest = my_id^mask;
315827bd09bSSatish Balay       if (my_id > dest)
3163fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
317827bd09bSSatish Balay       else
318827bd09bSSatish Balay 	{
3193fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
320827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
321827bd09bSSatish Balay 	}
322827bd09bSSatish Balay     }
323827bd09bSSatish Balay 
324827bd09bSSatish Balay   if (edge==dim)
325827bd09bSSatish Balay     {mask>>=1;}
326827bd09bSSatish Balay   else
327827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
328827bd09bSSatish Balay 
329827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
330827bd09bSSatish Balay     {
331827bd09bSSatish Balay       if (my_id%mask)
332827bd09bSSatish Balay 	{continue;}
333827bd09bSSatish Balay 
334827bd09bSSatish Balay       dest = my_id^mask;
335827bd09bSSatish Balay       if (my_id < dest)
3363fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
337827bd09bSSatish Balay       else
338827bd09bSSatish Balay 	{
3393fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
340827bd09bSSatish Balay 	}
341827bd09bSSatish Balay     }
3423fdc5746SBarry Smith         PetscFunctionReturn(0);
343827bd09bSSatish Balay }
344827bd09bSSatish Balay 
3457b1ae94cSBarry Smith /******************************************************************************/
3460924e98cSBarry Smith PetscErrorCode ssgl_radd( PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
347827bd09bSSatish Balay {
3483fdc5746SBarry Smith   PetscInt       edge, type, dest, mask;
3493fdc5746SBarry Smith   PetscInt       stage_n;
350827bd09bSSatish Balay   MPI_Status     status;
3513fdc5746SBarry Smith   PetscErrorCode ierr;
352827bd09bSSatish Balay 
3533fdc5746SBarry Smith    PetscFunctionBegin;
354827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
355827bd09bSSatish Balay   if (!p_init)
356827bd09bSSatish Balay     {comm_init();}
357827bd09bSSatish Balay 
358827bd09bSSatish Balay 
359827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
360827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
361827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
362827bd09bSSatish Balay     {
363827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
364827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
365827bd09bSSatish Balay 	{
366827bd09bSSatish Balay 	  dest = edge_node[edge];
367827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
368827bd09bSSatish Balay 	  if (my_id>dest)
3693fdc5746SBarry Smith           {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);}
370827bd09bSSatish Balay 	  else
371827bd09bSSatish Balay 	    {
372827bd09bSSatish Balay 	      type =  type - my_id + dest;
3733fdc5746SBarry Smith               ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
374827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
375827bd09bSSatish Balay 	    }
376827bd09bSSatish Balay 	}
377827bd09bSSatish Balay       mask <<= 1;
378827bd09bSSatish Balay     }
379827bd09bSSatish Balay   mask>>=1;
380827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
381827bd09bSSatish Balay     {
382827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
383827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
384827bd09bSSatish Balay 	{
385827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
386827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
387827bd09bSSatish Balay 	  if (my_id<dest)
3883fdc5746SBarry Smith             {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);}
389827bd09bSSatish Balay 	  else
390827bd09bSSatish Balay 	    {
391827bd09bSSatish Balay 	      type =  type - my_id + dest;
3923fdc5746SBarry Smith               ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
393827bd09bSSatish Balay 	    }
394827bd09bSSatish Balay 	}
395827bd09bSSatish Balay       mask >>= 1;
396827bd09bSSatish Balay     }
3973fdc5746SBarry Smith   PetscFunctionReturn(0);
398827bd09bSSatish Balay }
399827bd09bSSatish Balay 
4007b1ae94cSBarry Smith /******************************************************************************/
40152f87cdaSBarry Smith PetscErrorCode new_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
402827bd09bSSatish Balay {
40352f87cdaSBarry Smith   PetscInt            edge, type, dest, mask;
40452f87cdaSBarry Smith   PetscInt            stage_n;
405827bd09bSSatish Balay   MPI_Status     status;
4063fdc5746SBarry Smith   PetscErrorCode ierr;
407827bd09bSSatish Balay 
4083fdc5746SBarry Smith    PetscFunctionBegin;
409827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
410827bd09bSSatish Balay   if (!p_init)
411827bd09bSSatish Balay     {comm_init();}
412827bd09bSSatish Balay 
413827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
414827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
415827bd09bSSatish Balay   for (mask=0, edge=0; edge<level; edge++, mask++)
416827bd09bSSatish Balay     {
417827bd09bSSatish Balay       stage_n = (segs[level] - segs[edge]);
418827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
419827bd09bSSatish Balay 	{
420827bd09bSSatish Balay 	  dest = edge_node[edge];
421827bd09bSSatish Balay 	  type = MSGTAG3 + my_id + (num_nodes*edge);
422827bd09bSSatish Balay 	  if (my_id>dest)
4233fdc5746SBarry Smith           {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);}
424827bd09bSSatish Balay 	  else
425827bd09bSSatish Balay 	    {
426827bd09bSSatish Balay 	      type =  type - my_id + dest;
4273fdc5746SBarry Smith               ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
428827bd09bSSatish Balay 	      rvec_add(vals+segs[edge], work, stage_n);
429827bd09bSSatish Balay 	    }
430827bd09bSSatish Balay 	}
431827bd09bSSatish Balay       mask <<= 1;
432827bd09bSSatish Balay     }
433827bd09bSSatish Balay   mask>>=1;
434827bd09bSSatish Balay   for (edge=0; edge<level; edge++)
435827bd09bSSatish Balay     {
436827bd09bSSatish Balay       stage_n = (segs[level] - segs[level-1-edge]);
437827bd09bSSatish Balay       if (stage_n && !(my_id & mask))
438827bd09bSSatish Balay 	{
439827bd09bSSatish Balay 	  dest = edge_node[level-edge-1];
440827bd09bSSatish Balay 	  type = MSGTAG6 + my_id + (num_nodes*edge);
441827bd09bSSatish Balay 	  if (my_id<dest)
4423fdc5746SBarry Smith             {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);}
443827bd09bSSatish Balay 	  else
444827bd09bSSatish Balay 	    {
445827bd09bSSatish Balay 	      type =  type - my_id + dest;
4463fdc5746SBarry Smith               ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
447827bd09bSSatish Balay 	    }
448827bd09bSSatish Balay 	}
449827bd09bSSatish Balay       mask >>= 1;
450827bd09bSSatish Balay     }
4513fdc5746SBarry Smith   PetscFunctionReturn(0);
452827bd09bSSatish Balay }
453827bd09bSSatish Balay 
4547b1ae94cSBarry Smith /***********************************comm.c*************************************/
45552f87cdaSBarry Smith PetscErrorCode giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
456827bd09bSSatish Balay {
45752f87cdaSBarry Smith   PetscInt            mask, edge;
45852f87cdaSBarry Smith   PetscInt            type, dest;
459827bd09bSSatish Balay   vfp            fp;
460827bd09bSSatish Balay   MPI_Status     status;
4613fdc5746SBarry Smith   PetscErrorCode ierr;
462827bd09bSSatish Balay 
4633fdc5746SBarry Smith    PetscFunctionBegin;
464827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
465827bd09bSSatish Balay   if (!vals||!work||!oprs)
466*e32f2f54SBarry Smith     {SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);}
467827bd09bSSatish Balay 
468827bd09bSSatish Balay   /* non-uniform should have at least two entries */
469827bd09bSSatish Balay   if ((oprs[0] == NON_UNIFORM)&&(n<2))
470*e32f2f54SBarry Smith     {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: non_uniform and n=0,1?");}
471827bd09bSSatish Balay 
472827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
473827bd09bSSatish Balay   if (!p_init)
474827bd09bSSatish Balay     {comm_init();}
475827bd09bSSatish Balay 
476827bd09bSSatish Balay   /* if there's nothing to do return */
477827bd09bSSatish Balay   if ((num_nodes<2)||(!n)||(dim<=0))
4783fdc5746SBarry Smith     {  PetscFunctionReturn(0);}
479827bd09bSSatish Balay 
480827bd09bSSatish Balay   /* the error msg says it all!!! */
481827bd09bSSatish Balay   if (modfl_num_nodes)
482*e32f2f54SBarry Smith     {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: num_nodes not a power of 2!?!");}
483827bd09bSSatish Balay 
484827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
485827bd09bSSatish Balay   if (n<0)
486*e32f2f54SBarry Smith     {SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: n=%D<0?",n);}
487827bd09bSSatish Balay 
488827bd09bSSatish Balay   /* can't do more dimensions then exist */
48939945688SSatish Balay   dim = PetscMin(dim,i_log2_num_nodes);
490827bd09bSSatish Balay 
491827bd09bSSatish Balay   /* advance to list of n operations for custom */
492827bd09bSSatish Balay   if ((type=oprs[0])==NON_UNIFORM)
493827bd09bSSatish Balay     {oprs++;}
494827bd09bSSatish Balay 
495d890fc11SSatish Balay   if (!(fp = (vfp) ivec_fct_addr(type))){
496f1ed62a8SBarry Smith     ierr = PetscInfo(0,"giop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
497827bd09bSSatish Balay     fp = (vfp) oprs;
498827bd09bSSatish Balay   }
499827bd09bSSatish Balay 
500827bd09bSSatish Balay   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
501827bd09bSSatish Balay     {
502827bd09bSSatish Balay       dest = my_id^mask;
503827bd09bSSatish Balay       if (my_id > dest)
5043fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
505827bd09bSSatish Balay       else
506827bd09bSSatish Balay 	{
5073fdc5746SBarry Smith 	  ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
508827bd09bSSatish Balay 	  (*fp)(vals, work, n, oprs);
509827bd09bSSatish Balay 	}
510827bd09bSSatish Balay     }
511827bd09bSSatish Balay 
512827bd09bSSatish Balay   if (edge==dim)
513827bd09bSSatish Balay     {mask>>=1;}
514827bd09bSSatish Balay   else
515827bd09bSSatish Balay     {while (++edge<dim) {mask<<=1;}}
516827bd09bSSatish Balay 
517827bd09bSSatish Balay   for (edge=0; edge<dim; edge++,mask>>=1)
518827bd09bSSatish Balay     {
519827bd09bSSatish Balay       if (my_id%mask)
520827bd09bSSatish Balay 	{continue;}
521827bd09bSSatish Balay 
522827bd09bSSatish Balay       dest = my_id^mask;
523827bd09bSSatish Balay       if (my_id < dest)
5243fdc5746SBarry Smith 	{ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
525827bd09bSSatish Balay       else
526827bd09bSSatish Balay 	{
5273fdc5746SBarry Smith 	  ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
528827bd09bSSatish Balay 	}
529827bd09bSSatish Balay     }
5303fdc5746SBarry Smith   PetscFunctionReturn(0);
531827bd09bSSatish Balay }
532