xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision db4deed73f06ae556a0f2e02a0fd6e016e156849)
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;
47*db4deed7SKarl 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*db4deed7SKarl 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*db4deed7SKarl Rupp   else if (PCTFS_my_id >= PCTFS_floor_num_nodes) { edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1); }
59*db4deed7SKarl 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*db4deed7SKarl Rupp   if (!p_init) {PCTFS_comm_init();}
81827bd09bSSatish Balay 
82827bd09bSSatish Balay   /* if there's nothing to do return */
83*db4deed7SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)) {
843fdc5746SBarry Smith     PetscFunctionReturn(0);
85827bd09bSSatish Balay   }
86827bd09bSSatish Balay 
8771a0148aSBarry Smith 
88827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
89b1c944f5SJed Brown   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n);
90827bd09bSSatish Balay 
91827bd09bSSatish Balay   /* advance to list of n operations for custom */
92*db4deed7SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) {oprs++;}
93827bd09bSSatish Balay 
94827bd09bSSatish Balay   /* major league hack */
95ca8e9878SJed Brown   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) {
96b1c944f5SJed Brown     ierr = PetscInfo(0,"PCTFS_giop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
97827bd09bSSatish Balay     fp = (vfp) oprs;
98827bd09bSSatish Balay   }
99827bd09bSSatish Balay 
100827bd09bSSatish Balay   /* all msgs will be of the same length */
101827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
102*db4deed7SKarl Rupp   if (edge_not_pow_2) {
103*db4deed7SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
104*db4deed7SKarl Rupp     else {
1053fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
106827bd09bSSatish Balay       (*fp)(vals,work,n,oprs);
107827bd09bSSatish Balay     }
108827bd09bSSatish Balay   }
109827bd09bSSatish Balay 
110827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
111*db4deed7SKarl Rupp   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
112*db4deed7SKarl Rupp     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
113b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
114*db4deed7SKarl Rupp       if (PCTFS_my_id > dest) {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
115*db4deed7SKarl Rupp       else {
1163fdc5746SBarry Smith         ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
117827bd09bSSatish Balay         (*fp)(vals, work, n, oprs);
118827bd09bSSatish Balay       }
119827bd09bSSatish Balay     }
120827bd09bSSatish Balay 
121b1c944f5SJed Brown     mask=PCTFS_floor_num_nodes>>1;
122*db4deed7SKarl Rupp     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
123*db4deed7SKarl Rupp       if (PCTFS_my_id%mask) {continue;}
124827bd09bSSatish Balay 
125b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
126*db4deed7SKarl Rupp       if (PCTFS_my_id < dest) {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
127*db4deed7SKarl Rupp       else {
1283fdc5746SBarry Smith         ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
129827bd09bSSatish Balay       }
130827bd09bSSatish Balay     }
131827bd09bSSatish Balay   }
132827bd09bSSatish Balay 
133827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
134*db4deed7SKarl Rupp   if (edge_not_pow_2) {
135b1c944f5SJed Brown     if (PCTFS_my_id >= PCTFS_floor_num_nodes)
136827bd09bSSatish Balay     {
1373fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
138827bd09bSSatish Balay     }
139*db4deed7SKarl Rupp     else {ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
140827bd09bSSatish Balay   }
1413fdc5746SBarry Smith   PetscFunctionReturn(0);
142827bd09bSSatish Balay }
143827bd09bSSatish Balay 
1447b1ae94cSBarry Smith /***********************************comm.c*************************************/
145b1c944f5SJed Brown PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
146827bd09bSSatish Balay {
1473fdc5746SBarry Smith   PetscInt       mask, edge;
1483fdc5746SBarry Smith   PetscInt       type, dest;
149827bd09bSSatish Balay   vfp            fp;
150827bd09bSSatish Balay   MPI_Status     status;
1513fdc5746SBarry Smith   PetscErrorCode ierr;
152827bd09bSSatish Balay 
1533fdc5746SBarry Smith    PetscFunctionBegin;
154827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
155b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
156827bd09bSSatish Balay 
157827bd09bSSatish Balay   /* non-uniform should have at least two entries */
158b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?");
159827bd09bSSatish Balay 
160827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
161*db4deed7SKarl Rupp   if (!p_init) {PCTFS_comm_init();}
162827bd09bSSatish Balay 
163827bd09bSSatish Balay   /* if there's nothing to do return */
164*db4deed7SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)) { PetscFunctionReturn(0); }
165827bd09bSSatish Balay 
166827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
167c1235816SBarry Smith   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);
168827bd09bSSatish Balay 
169827bd09bSSatish Balay   /* advance to list of n operations for custom */
170*db4deed7SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) {oprs++;}
171827bd09bSSatish Balay 
172ca8e9878SJed Brown   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) {
173b1c944f5SJed Brown     ierr = PetscInfo(0,"PCTFS_grop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
174827bd09bSSatish Balay     fp = (vfp) oprs;
175827bd09bSSatish Balay   }
176827bd09bSSatish Balay 
177827bd09bSSatish Balay   /* all msgs will be of the same length */
178827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
179827bd09bSSatish Balay   if (edge_not_pow_2)
180827bd09bSSatish Balay   {
181*db4deed7SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
182*db4deed7SKarl Rupp     else {
1833fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
184827bd09bSSatish Balay       (*fp)(vals,work,n,oprs);
185827bd09bSSatish Balay     }
186827bd09bSSatish Balay   }
187827bd09bSSatish Balay 
188827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
189*db4deed7SKarl Rupp   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
190*db4deed7SKarl Rupp     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
191b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
192*db4deed7SKarl Rupp       if (PCTFS_my_id > dest) {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
193*db4deed7SKarl Rupp       else {
1943fdc5746SBarry Smith         ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
195827bd09bSSatish Balay         (*fp)(vals, work, n, oprs);
196827bd09bSSatish Balay       }
197827bd09bSSatish Balay     }
198827bd09bSSatish Balay 
199b1c944f5SJed Brown     mask=PCTFS_floor_num_nodes>>1;
200*db4deed7SKarl Rupp     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
201*db4deed7SKarl Rupp       if (PCTFS_my_id%mask) { continue; }
202827bd09bSSatish Balay 
203b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
204*db4deed7SKarl Rupp       if (PCTFS_my_id < dest) { ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr); }
205*db4deed7SKarl Rupp       else {
2063fdc5746SBarry Smith         ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
207827bd09bSSatish Balay       }
208827bd09bSSatish Balay     }
209827bd09bSSatish Balay   }
210827bd09bSSatish Balay 
211827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
212*db4deed7SKarl Rupp   if (edge_not_pow_2) {
213*db4deed7SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
2143fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
215*db4deed7SKarl Rupp     } else {
216*db4deed7SKarl Rupp       ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
217827bd09bSSatish Balay     }
218827bd09bSSatish Balay   }
2193fdc5746SBarry Smith   PetscFunctionReturn(0);
220827bd09bSSatish Balay }
221827bd09bSSatish Balay 
2227b1ae94cSBarry Smith /***********************************comm.c*************************************/
223b1c944f5SJed Brown PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
224827bd09bSSatish Balay {
2253fdc5746SBarry Smith   PetscInt       mask, edge;
2263fdc5746SBarry Smith   PetscInt       type, dest;
227827bd09bSSatish Balay   vfp            fp;
228827bd09bSSatish Balay   MPI_Status     status;
2293fdc5746SBarry Smith   PetscErrorCode ierr;
230827bd09bSSatish Balay 
2313fdc5746SBarry Smith    PetscFunctionBegin;
232827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
233b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
234827bd09bSSatish Balay 
235827bd09bSSatish Balay   /* non-uniform should have at least two entries */
236b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?");
237827bd09bSSatish Balay 
238827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
239*db4deed7SKarl Rupp   if (!p_init) {PCTFS_comm_init();}
240827bd09bSSatish Balay 
241827bd09bSSatish Balay   /* if there's nothing to do return */
242*db4deed7SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) {PetscFunctionReturn(0);}
243827bd09bSSatish Balay 
244827bd09bSSatish Balay   /* the error msg says it all!!! */
245b1c944f5SJed Brown   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!");
246827bd09bSSatish Balay 
247827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
248b1c944f5SJed Brown   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n);
249827bd09bSSatish Balay 
250827bd09bSSatish Balay   /* can't do more dimensions then exist */
251b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
252827bd09bSSatish Balay 
253827bd09bSSatish Balay   /* advance to list of n operations for custom */
254*db4deed7SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) {oprs++;}
255827bd09bSSatish Balay 
256ca8e9878SJed Brown   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) {
257b1c944f5SJed Brown     ierr = PetscInfo(0,"PCTFS_grop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
258827bd09bSSatish Balay     fp = (vfp) oprs;
259827bd09bSSatish Balay   }
260827bd09bSSatish Balay 
261*db4deed7SKarl Rupp   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
262b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
263*db4deed7SKarl Rupp     if (PCTFS_my_id > dest) { ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr); }
264*db4deed7SKarl Rupp     else {
2653fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
266827bd09bSSatish Balay       (*fp)(vals, work, n, oprs);
267827bd09bSSatish Balay     }
268827bd09bSSatish Balay   }
269827bd09bSSatish Balay 
270*db4deed7SKarl Rupp   if (edge==dim) { mask>>=1; }
271*db4deed7SKarl Rupp   else {
272*db4deed7SKarl Rupp     while (++edge<dim) {mask<<=1;}
273*db4deed7SKarl Rupp   }
274827bd09bSSatish Balay 
275*db4deed7SKarl Rupp   for (edge=0; edge<dim; edge++,mask>>=1) {
276*db4deed7SKarl Rupp     if (PCTFS_my_id%mask) { continue; }
277827bd09bSSatish Balay 
278b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
279*db4deed7SKarl Rupp     if (PCTFS_my_id < dest) { ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr); }
280*db4deed7SKarl Rupp     else {
2813fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
282827bd09bSSatish Balay     }
283827bd09bSSatish Balay   }
2843fdc5746SBarry Smith   PetscFunctionReturn(0);
285827bd09bSSatish Balay }
286827bd09bSSatish Balay 
2877b1ae94cSBarry Smith /******************************************************************************/
288b1c944f5SJed Brown PetscErrorCode PCTFS_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
289827bd09bSSatish Balay {
2903fdc5746SBarry Smith   PetscInt       edge, type, dest, mask;
2913fdc5746SBarry Smith   PetscInt       stage_n;
292827bd09bSSatish Balay   MPI_Status     status;
2933fdc5746SBarry Smith   PetscErrorCode ierr;
294827bd09bSSatish Balay 
2953fdc5746SBarry Smith    PetscFunctionBegin;
296827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
297*db4deed7SKarl Rupp   if (!p_init) { PCTFS_comm_init(); }
298827bd09bSSatish Balay 
299827bd09bSSatish Balay 
300827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
301827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
302*db4deed7SKarl Rupp   for (mask=0, edge=0; edge<level; edge++, mask++) {
303827bd09bSSatish Balay     stage_n = (segs[level] - segs[edge]);
304*db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
305827bd09bSSatish Balay       dest = edge_node[edge];
306b1c944f5SJed Brown       type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge);
307*db4deed7SKarl Rupp       if (PCTFS_my_id>dest) { ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr); }
308*db4deed7SKarl Rupp       else {
309b1c944f5SJed Brown         type =  type - PCTFS_my_id + dest;
3103fdc5746SBarry Smith         ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
311ca8e9878SJed Brown         PCTFS_rvec_add(vals+segs[edge], work, stage_n);
312827bd09bSSatish Balay       }
313827bd09bSSatish Balay     }
314827bd09bSSatish Balay     mask <<= 1;
315827bd09bSSatish Balay   }
316827bd09bSSatish Balay   mask>>=1;
317*db4deed7SKarl Rupp   for (edge=0; edge<level; edge++) {
318827bd09bSSatish Balay     stage_n = (segs[level] - segs[level-1-edge]);
319*db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
320827bd09bSSatish Balay       dest = edge_node[level-edge-1];
321b1c944f5SJed Brown       type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge);
322*db4deed7SKarl Rupp       if (PCTFS_my_id<dest) { ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr); }
323*db4deed7SKarl Rupp       else {
324b1c944f5SJed Brown         type =  type - PCTFS_my_id + dest;
3253fdc5746SBarry Smith         ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
326827bd09bSSatish Balay       }
327827bd09bSSatish Balay     }
328827bd09bSSatish Balay     mask >>= 1;
329827bd09bSSatish Balay   }
3303fdc5746SBarry Smith   PetscFunctionReturn(0);
331827bd09bSSatish Balay }
332827bd09bSSatish Balay 
3337b1ae94cSBarry Smith /***********************************comm.c*************************************/
334b1c944f5SJed Brown PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
335827bd09bSSatish Balay {
33652f87cdaSBarry Smith   PetscInt            mask, edge;
33752f87cdaSBarry Smith   PetscInt            type, dest;
338827bd09bSSatish Balay   vfp            fp;
339827bd09bSSatish Balay   MPI_Status     status;
3403fdc5746SBarry Smith   PetscErrorCode ierr;
341827bd09bSSatish Balay 
3423fdc5746SBarry Smith   PetscFunctionBegin;
343827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
344b1c944f5SJed Brown   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
345827bd09bSSatish Balay 
346827bd09bSSatish Balay   /* non-uniform should have at least two entries */
347b1c944f5SJed Brown   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?");
348827bd09bSSatish Balay 
349827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
350*db4deed7SKarl Rupp   if (!p_init) { PCTFS_comm_init(); }
351827bd09bSSatish Balay 
352827bd09bSSatish Balay   /* if there's nothing to do return */
353*db4deed7SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) { PetscFunctionReturn(0); }
354827bd09bSSatish Balay 
355827bd09bSSatish Balay   /* the error msg says it all!!! */
356b1c944f5SJed Brown   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!");
357827bd09bSSatish Balay 
358827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
359b1c944f5SJed Brown   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n);
360827bd09bSSatish Balay 
361827bd09bSSatish Balay   /* can't do more dimensions then exist */
362b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
363827bd09bSSatish Balay 
364827bd09bSSatish Balay   /* advance to list of n operations for custom */
365*db4deed7SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) { oprs++; }
366827bd09bSSatish Balay 
367ca8e9878SJed Brown   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) {
368b1c944f5SJed Brown     ierr = PetscInfo(0,"PCTFS_giop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr);
369827bd09bSSatish Balay     fp = (vfp) oprs;
370827bd09bSSatish Balay   }
371827bd09bSSatish Balay 
372*db4deed7SKarl Rupp   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
373b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
374*db4deed7SKarl Rupp     if (PCTFS_my_id > dest) {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
375*db4deed7SKarl Rupp     else {
3763fdc5746SBarry Smith       ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
377827bd09bSSatish Balay       (*fp)(vals, work, n, oprs);
378827bd09bSSatish Balay     }
379827bd09bSSatish Balay   }
380827bd09bSSatish Balay 
381*db4deed7SKarl Rupp   if (edge==dim) {mask>>=1;}
382*db4deed7SKarl Rupp   else {while (++edge<dim) {mask<<=1;}}
383827bd09bSSatish Balay 
384*db4deed7SKarl Rupp   for (edge=0; edge<dim; edge++,mask>>=1) {
385*db4deed7SKarl Rupp     if (PCTFS_my_id%mask) {continue;}
386827bd09bSSatish Balay 
387b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
388*db4deed7SKarl Rupp     if (PCTFS_my_id < dest) {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);}
389*db4deed7SKarl Rupp     else {
3903fdc5746SBarry Smith       ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
391827bd09bSSatish Balay     }
392827bd09bSSatish Balay   }
3933fdc5746SBarry Smith   PetscFunctionReturn(0);
394827bd09bSSatish Balay }
395