xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
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*************************************/
16c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
17827bd09bSSatish Balay 
18827bd09bSSatish Balay 
19827bd09bSSatish Balay /* global program control variables - explicitly exported */
20b1c944f5SJed Brown PetscMPIInt PCTFS_my_id            = 0;
21b1c944f5SJed Brown PetscMPIInt PCTFS_num_nodes        = 1;
22b1c944f5SJed Brown PetscMPIInt PCTFS_floor_num_nodes  = 0;
23b1c944f5SJed Brown PetscMPIInt PCTFS_i_log2_num_nodes = 0;
24827bd09bSSatish Balay 
25827bd09bSSatish Balay /* global program control variables */
2652f87cdaSBarry Smith static PetscInt p_init = 0;
2752f87cdaSBarry Smith static PetscInt modfl_num_nodes;
2852f87cdaSBarry Smith static PetscInt edge_not_pow_2;
29827bd09bSSatish Balay 
3052f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt)*32];
31827bd09bSSatish Balay 
327b1ae94cSBarry Smith /***********************************comm.c*************************************/
33b1c944f5SJed Brown PetscErrorCode PCTFS_comm_init(void)
34827bd09bSSatish Balay {
35827bd09bSSatish Balay 
363fdc5746SBarry Smith   if (p_init++) PetscFunctionReturn(0);
37827bd09bSSatish Balay 
38b1c944f5SJed Brown   MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes);
39b1c944f5SJed Brown   MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id);
40827bd09bSSatish Balay 
41b1c944f5SJed Brown   if (PCTFS_num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!");
42827bd09bSSatish Balay 
43ca8e9878SJed Brown   PCTFS_ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32);
44827bd09bSSatish Balay 
45b1c944f5SJed Brown   PCTFS_floor_num_nodes  = 1;
46b1c944f5SJed Brown   PCTFS_i_log2_num_nodes = modfl_num_nodes = 0;
47db4deed7SKarl Rupp   while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) {
48b1c944f5SJed Brown     edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes;
49b1c944f5SJed Brown     PCTFS_floor_num_nodes           <<= 1;
50b1c944f5SJed Brown     PCTFS_i_log2_num_nodes++;
51827bd09bSSatish Balay   }
52827bd09bSSatish Balay 
53b1c944f5SJed Brown   PCTFS_i_log2_num_nodes--;
54b1c944f5SJed Brown   PCTFS_floor_num_nodes >>= 1;
55b1c944f5SJed Brown   modfl_num_nodes         = (PCTFS_num_nodes - PCTFS_floor_num_nodes);
56827bd09bSSatish Balay 
57*2fa5cd67SKarl Rupp   if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1);
58*2fa5cd67SKarl Rupp   else if (PCTFS_my_id >= PCTFS_floor_num_nodes) edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1);
59*2fa5cd67SKarl Rupp   else edge_not_pow_2 = 0;
603fdc5746SBarry Smith   PetscFunctionReturn(0);
61827bd09bSSatish Balay }
62827bd09bSSatish Balay 
637b1ae94cSBarry Smith /***********************************comm.c*************************************/
64b1c944f5SJed Brown PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
65827bd09bSSatish Balay {
663fdc5746SBarry Smith   PetscInt   mask, edge;
673fdc5746SBarry Smith   PetscInt   type, dest;
68827bd09bSSatish Balay   vfp        fp;
69827bd09bSSatish Balay   MPI_Status status;
703fdc5746SBarry Smith   PetscInt   ierr;
71827bd09bSSatish Balay 
723fdc5746SBarry Smith   PetscFunctionBegin;
73827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
74b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
75827bd09bSSatish Balay 
76827bd09bSSatish Balay   /* non-uniform should have at least two entries */
77b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?");
78827bd09bSSatish Balay 
79827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
80*2fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
81827bd09bSSatish Balay 
82827bd09bSSatish Balay   /* if there's nothing to do return */
83*2fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0);
8471a0148aSBarry Smith 
85827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
86b1c944f5SJed Brown   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n);
87827bd09bSSatish Balay 
88827bd09bSSatish Balay   /* advance to list of n operations for custom */
89*2fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
90827bd09bSSatish Balay 
91827bd09bSSatish Balay   /* major league hack */
92*2fa5cd67SKarl Rupp   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: Could not retrieve function pointer!\n");
93827bd09bSSatish Balay 
94827bd09bSSatish Balay   /* all msgs will be of the same length */
95827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
96db4deed7SKarl Rupp   if (edge_not_pow_2) {
97*2fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
98*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
99*2fa5cd67SKarl Rupp     } else {
1003fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
101827bd09bSSatish Balay       (*fp)(vals,work,n,oprs);
102827bd09bSSatish Balay     }
103827bd09bSSatish Balay   }
104827bd09bSSatish Balay 
105827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
106db4deed7SKarl Rupp   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
107db4deed7SKarl Rupp     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
108b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
109*2fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
110*2fa5cd67SKarl Rupp         ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
111*2fa5cd67SKarl Rupp       } else {
1123fdc5746SBarry Smith         ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
113827bd09bSSatish Balay         (*fp)(vals, work, n, oprs);
114827bd09bSSatish Balay       }
115827bd09bSSatish Balay     }
116827bd09bSSatish Balay 
117b1c944f5SJed Brown     mask=PCTFS_floor_num_nodes>>1;
118db4deed7SKarl Rupp     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
119*2fa5cd67SKarl Rupp       if (PCTFS_my_id%mask) continue;
120827bd09bSSatish Balay 
121b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
122*2fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
123*2fa5cd67SKarl Rupp         ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
124*2fa5cd67SKarl Rupp       } else {
1253fdc5746SBarry Smith         ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
126827bd09bSSatish Balay       }
127827bd09bSSatish Balay     }
128827bd09bSSatish Balay   }
129827bd09bSSatish Balay 
130827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
131db4deed7SKarl Rupp   if (edge_not_pow_2) {
132*2fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
1333fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
134*2fa5cd67SKarl Rupp     } else {
135*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
136827bd09bSSatish Balay     }
137827bd09bSSatish Balay   }
1383fdc5746SBarry Smith   PetscFunctionReturn(0);
139827bd09bSSatish Balay }
140827bd09bSSatish Balay 
1417b1ae94cSBarry Smith /***********************************comm.c*************************************/
142b1c944f5SJed Brown PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
143827bd09bSSatish Balay {
1443fdc5746SBarry Smith   PetscInt       mask, edge;
1453fdc5746SBarry Smith   PetscInt       type, dest;
146827bd09bSSatish Balay   vfp            fp;
147827bd09bSSatish Balay   MPI_Status     status;
1483fdc5746SBarry Smith   PetscErrorCode ierr;
149827bd09bSSatish Balay 
1503fdc5746SBarry Smith   PetscFunctionBegin;
151827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
152b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
153827bd09bSSatish Balay 
154827bd09bSSatish Balay   /* non-uniform should have at least two entries */
155b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?");
156827bd09bSSatish Balay 
157827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
158*2fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
159827bd09bSSatish Balay 
160827bd09bSSatish Balay   /* if there's nothing to do return */
161*2fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0);
162827bd09bSSatish Balay 
163827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
164c1235816SBarry Smith   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);
165827bd09bSSatish Balay 
166827bd09bSSatish Balay   /* advance to list of n operations for custom */
167*2fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
168827bd09bSSatish Balay 
169*2fa5cd67SKarl Rupp   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: Could not retrieve function pointer!\n");
170827bd09bSSatish Balay 
171827bd09bSSatish Balay   /* all msgs will be of the same length */
172827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
173*2fa5cd67SKarl Rupp   if (edge_not_pow_2) {
174*2fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
175*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
176*2fa5cd67SKarl Rupp     } else {
1773fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
178827bd09bSSatish Balay       (*fp)(vals,work,n,oprs);
179827bd09bSSatish Balay     }
180827bd09bSSatish Balay   }
181827bd09bSSatish Balay 
182827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
183db4deed7SKarl Rupp   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
184db4deed7SKarl Rupp     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
185b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
186*2fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
187*2fa5cd67SKarl Rupp         ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
188*2fa5cd67SKarl Rupp       } else {
1893fdc5746SBarry Smith         ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
190827bd09bSSatish Balay         (*fp)(vals, work, n, oprs);
191827bd09bSSatish Balay       }
192827bd09bSSatish Balay     }
193827bd09bSSatish Balay 
194b1c944f5SJed Brown     mask=PCTFS_floor_num_nodes>>1;
195db4deed7SKarl Rupp     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
196*2fa5cd67SKarl Rupp       if (PCTFS_my_id%mask) continue;
197827bd09bSSatish Balay 
198b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
199*2fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
200*2fa5cd67SKarl Rupp         ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
201*2fa5cd67SKarl Rupp       } else {
2023fdc5746SBarry Smith         ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
203827bd09bSSatish Balay       }
204827bd09bSSatish Balay     }
205827bd09bSSatish Balay   }
206827bd09bSSatish Balay 
207827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
208db4deed7SKarl Rupp   if (edge_not_pow_2) {
209db4deed7SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
2103fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
211db4deed7SKarl Rupp     } else {
212db4deed7SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
213827bd09bSSatish Balay     }
214827bd09bSSatish Balay   }
2153fdc5746SBarry Smith   PetscFunctionReturn(0);
216827bd09bSSatish Balay }
217827bd09bSSatish Balay 
2187b1ae94cSBarry Smith /***********************************comm.c*************************************/
219b1c944f5SJed Brown PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
220827bd09bSSatish Balay {
2213fdc5746SBarry Smith   PetscInt       mask, edge;
2223fdc5746SBarry Smith   PetscInt       type, dest;
223827bd09bSSatish Balay   vfp            fp;
224827bd09bSSatish Balay   MPI_Status     status;
2253fdc5746SBarry Smith   PetscErrorCode ierr;
226827bd09bSSatish Balay 
2273fdc5746SBarry Smith   PetscFunctionBegin;
228827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
229b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
230827bd09bSSatish Balay 
231827bd09bSSatish Balay   /* non-uniform should have at least two entries */
232b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?");
233827bd09bSSatish Balay 
234827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
235*2fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
236827bd09bSSatish Balay 
237827bd09bSSatish Balay   /* if there's nothing to do return */
238*2fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0);
239827bd09bSSatish Balay 
240827bd09bSSatish Balay   /* the error msg says it all!!! */
241b1c944f5SJed Brown   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!");
242827bd09bSSatish Balay 
243827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
244b1c944f5SJed Brown   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n);
245827bd09bSSatish Balay 
246827bd09bSSatish Balay   /* can't do more dimensions then exist */
247b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
248827bd09bSSatish Balay 
249827bd09bSSatish Balay   /* advance to list of n operations for custom */
250*2fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
251827bd09bSSatish Balay 
252*2fa5cd67SKarl Rupp   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: Could not retrieve function pointer!\n");
253827bd09bSSatish Balay 
254db4deed7SKarl Rupp   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
255b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
256*2fa5cd67SKarl Rupp     if (PCTFS_my_id > dest) {
257*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
258*2fa5cd67SKarl Rupp     } else {
2593fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
260827bd09bSSatish Balay       (*fp)(vals, work, n, oprs);
261827bd09bSSatish Balay     }
262827bd09bSSatish Balay   }
263827bd09bSSatish Balay 
264*2fa5cd67SKarl Rupp   if (edge==dim) mask>>=1;
265db4deed7SKarl Rupp   else {
266*2fa5cd67SKarl Rupp     while (++edge<dim) mask<<=1;
267db4deed7SKarl Rupp   }
268827bd09bSSatish Balay 
269db4deed7SKarl Rupp   for (edge=0; edge<dim; edge++,mask>>=1) {
270*2fa5cd67SKarl Rupp     if (PCTFS_my_id%mask) continue;
271827bd09bSSatish Balay 
272b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
273*2fa5cd67SKarl Rupp     if (PCTFS_my_id < dest) {
274*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
275*2fa5cd67SKarl Rupp     } else {
2763fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
277827bd09bSSatish Balay     }
278827bd09bSSatish Balay   }
2793fdc5746SBarry Smith   PetscFunctionReturn(0);
280827bd09bSSatish Balay }
281827bd09bSSatish Balay 
2827b1ae94cSBarry Smith /******************************************************************************/
283b1c944f5SJed Brown PetscErrorCode PCTFS_ssgl_radd(PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
284827bd09bSSatish Balay {
2853fdc5746SBarry Smith   PetscInt       edge, type, dest, mask;
2863fdc5746SBarry Smith   PetscInt       stage_n;
287827bd09bSSatish Balay   MPI_Status     status;
2883fdc5746SBarry Smith   PetscErrorCode ierr;
289827bd09bSSatish Balay 
2903fdc5746SBarry Smith   PetscFunctionBegin;
291827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
292*2fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
293827bd09bSSatish Balay 
294827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
295827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
296db4deed7SKarl Rupp   for (mask=0, edge=0; edge<level; edge++, mask++) {
297827bd09bSSatish Balay     stage_n = (segs[level] - segs[edge]);
298db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
299827bd09bSSatish Balay       dest = edge_node[edge];
300b1c944f5SJed Brown       type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge);
301*2fa5cd67SKarl Rupp       if (PCTFS_my_id>dest) {
302*2fa5cd67SKarl Rupp         ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);
303*2fa5cd67SKarl Rupp       } else {
304b1c944f5SJed Brown         type =  type - PCTFS_my_id + dest;
3053fdc5746SBarry Smith         ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
306ca8e9878SJed Brown         PCTFS_rvec_add(vals+segs[edge], work, stage_n);
307827bd09bSSatish Balay       }
308827bd09bSSatish Balay     }
309827bd09bSSatish Balay     mask <<= 1;
310827bd09bSSatish Balay   }
311827bd09bSSatish Balay   mask>>=1;
312db4deed7SKarl Rupp   for (edge=0; edge<level; edge++) {
313827bd09bSSatish Balay     stage_n = (segs[level] - segs[level-1-edge]);
314db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
315827bd09bSSatish Balay       dest = edge_node[level-edge-1];
316b1c944f5SJed Brown       type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge);
317*2fa5cd67SKarl Rupp       if (PCTFS_my_id<dest) {
318*2fa5cd67SKarl Rupp         ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);
319*2fa5cd67SKarl Rupp       } else {
320b1c944f5SJed Brown         type =  type - PCTFS_my_id + dest;
3213fdc5746SBarry Smith         ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
322827bd09bSSatish Balay       }
323827bd09bSSatish Balay     }
324827bd09bSSatish Balay     mask >>= 1;
325827bd09bSSatish Balay   }
3263fdc5746SBarry Smith   PetscFunctionReturn(0);
327827bd09bSSatish Balay }
328827bd09bSSatish Balay 
3297b1ae94cSBarry Smith /***********************************comm.c*************************************/
330b1c944f5SJed Brown PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
331827bd09bSSatish Balay {
33252f87cdaSBarry Smith   PetscInt       mask, edge;
33352f87cdaSBarry Smith   PetscInt       type, dest;
334827bd09bSSatish Balay   vfp            fp;
335827bd09bSSatish Balay   MPI_Status     status;
3363fdc5746SBarry Smith   PetscErrorCode ierr;
337827bd09bSSatish Balay 
3383fdc5746SBarry Smith   PetscFunctionBegin;
339827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
340b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
341827bd09bSSatish Balay 
342827bd09bSSatish Balay   /* non-uniform should have at least two entries */
343b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?");
344827bd09bSSatish Balay 
345827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
346*2fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
347827bd09bSSatish Balay 
348827bd09bSSatish Balay   /* if there's nothing to do return */
349*2fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0);
350827bd09bSSatish Balay 
351827bd09bSSatish Balay   /* the error msg says it all!!! */
352b1c944f5SJed Brown   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!");
353827bd09bSSatish Balay 
354827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
355b1c944f5SJed Brown   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n);
356827bd09bSSatish Balay 
357827bd09bSSatish Balay   /* can't do more dimensions then exist */
358b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
359827bd09bSSatish Balay 
360827bd09bSSatish Balay   /* advance to list of n operations for custom */
361*2fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
362827bd09bSSatish Balay 
363*2fa5cd67SKarl Rupp   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: Could not retrieve function pointer!\n");
364827bd09bSSatish Balay 
365db4deed7SKarl Rupp   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
366b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
367*2fa5cd67SKarl Rupp     if (PCTFS_my_id > dest) {
368*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
369*2fa5cd67SKarl Rupp     } else {
3703fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
371827bd09bSSatish Balay       (*fp)(vals, work, n, oprs);
372827bd09bSSatish Balay     }
373827bd09bSSatish Balay   }
374827bd09bSSatish Balay 
375*2fa5cd67SKarl Rupp   if (edge==dim) mask>>=1;
376*2fa5cd67SKarl Rupp   else {
377*2fa5cd67SKarl Rupp     while (++edge<dim) mask<<=1;
378*2fa5cd67SKarl Rupp   }
379827bd09bSSatish Balay 
380db4deed7SKarl Rupp   for (edge=0; edge<dim; edge++,mask>>=1) {
381*2fa5cd67SKarl Rupp     if (PCTFS_my_id%mask) continue;
382827bd09bSSatish Balay 
383b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
384*2fa5cd67SKarl Rupp     if (PCTFS_my_id < dest) {
385*2fa5cd67SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
386*2fa5cd67SKarl Rupp     } else {
3873fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
388827bd09bSSatish Balay     }
389827bd09bSSatish Balay   }
3903fdc5746SBarry Smith   PetscFunctionReturn(0);
391827bd09bSSatish Balay }
392