xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision 7827d75ba736e00c19d105c058b1c2ddcca945c7)
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 /* global program control variables - explicitly exported */
19b1c944f5SJed Brown PetscMPIInt PCTFS_my_id            = 0;
20b1c944f5SJed Brown PetscMPIInt PCTFS_num_nodes        = 1;
21b1c944f5SJed Brown PetscMPIInt PCTFS_floor_num_nodes  = 0;
22b1c944f5SJed Brown PetscMPIInt PCTFS_i_log2_num_nodes = 0;
23827bd09bSSatish Balay 
24827bd09bSSatish Balay /* global program control variables */
2552f87cdaSBarry Smith static PetscInt p_init = 0;
2652f87cdaSBarry Smith static PetscInt modfl_num_nodes;
2752f87cdaSBarry Smith static PetscInt edge_not_pow_2;
28827bd09bSSatish Balay 
2952f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt)*32];
30827bd09bSSatish Balay 
317b1ae94cSBarry Smith /***********************************comm.c*************************************/
32b1c944f5SJed Brown PetscErrorCode PCTFS_comm_init(void)
33827bd09bSSatish Balay {
34362febeeSStefano Zampini   PetscFunctionBegin;
353fdc5746SBarry Smith   if (p_init++) PetscFunctionReturn(0);
36827bd09bSSatish Balay 
37b1c944f5SJed Brown   MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes);
38b1c944f5SJed Brown   MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id);
39827bd09bSSatish Balay 
402c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PCTFS_num_nodes> (INT_MAX >> 1),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!");
41827bd09bSSatish Balay 
42ca8e9878SJed Brown   PCTFS_ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32);
43827bd09bSSatish Balay 
44b1c944f5SJed Brown   PCTFS_floor_num_nodes  = 1;
45b1c944f5SJed Brown   PCTFS_i_log2_num_nodes = modfl_num_nodes = 0;
46db4deed7SKarl Rupp   while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) {
47b1c944f5SJed Brown     edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes;
48b1c944f5SJed Brown     PCTFS_floor_num_nodes           <<= 1;
49b1c944f5SJed Brown     PCTFS_i_log2_num_nodes++;
50827bd09bSSatish Balay   }
51827bd09bSSatish Balay 
52b1c944f5SJed Brown   PCTFS_i_log2_num_nodes--;
53b1c944f5SJed Brown   PCTFS_floor_num_nodes >>= 1;
54b1c944f5SJed Brown   modfl_num_nodes         = (PCTFS_num_nodes - PCTFS_floor_num_nodes);
55827bd09bSSatish Balay 
562fa5cd67SKarl Rupp   if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1);
572fa5cd67SKarl Rupp   else if (PCTFS_my_id >= PCTFS_floor_num_nodes) edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1);
582fa5cd67SKarl Rupp   else edge_not_pow_2 = 0;
593fdc5746SBarry Smith   PetscFunctionReturn(0);
60827bd09bSSatish Balay }
61827bd09bSSatish Balay 
627b1ae94cSBarry Smith /***********************************comm.c*************************************/
63b1c944f5SJed Brown PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
64827bd09bSSatish Balay {
653fdc5746SBarry Smith   PetscInt   mask, edge;
663fdc5746SBarry Smith   PetscInt   type, dest;
67827bd09bSSatish Balay   vfp        fp;
68827bd09bSSatish Balay   MPI_Status status;
69827bd09bSSatish Balay 
703fdc5746SBarry Smith   PetscFunctionBegin;
71827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
72*7827d75bSBarry Smith   PetscCheck(vals && work && oprs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
73827bd09bSSatish Balay 
74827bd09bSSatish Balay   /* non-uniform should have at least two entries */
752c71b3e2SJacob Faibussowitsch   PetscCheckFalse((oprs[0] == NON_UNIFORM)&&(n<2),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?");
76827bd09bSSatish Balay 
77827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
782fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
79827bd09bSSatish Balay 
80827bd09bSSatish Balay   /* if there's nothing to do return */
812fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0);
8271a0148aSBarry Smith 
83827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
842c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n<0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n);
85827bd09bSSatish Balay 
86827bd09bSSatish Balay   /* advance to list of n operations for custom */
872fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
88827bd09bSSatish Balay 
89827bd09bSSatish Balay   /* major league hack */
90*7827d75bSBarry Smith   PetscCheck(fp = (vfp) PCTFS_ivec_fct_addr(type),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: Could not retrieve function pointer!");
91827bd09bSSatish Balay 
92827bd09bSSatish Balay   /* all msgs will be of the same length */
93827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
94db4deed7SKarl Rupp   if (edge_not_pow_2) {
952fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
969566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD));
972fa5cd67SKarl Rupp     } else {
989566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status));
99827bd09bSSatish Balay       (*fp)(vals,work,n,oprs);
100827bd09bSSatish Balay     }
101827bd09bSSatish Balay   }
102827bd09bSSatish Balay 
103827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
104db4deed7SKarl Rupp   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
105db4deed7SKarl Rupp     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
106b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
1072fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
1089566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD));
1092fa5cd67SKarl Rupp       } else {
1109566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status));
111827bd09bSSatish Balay         (*fp)(vals, work, n, oprs);
112827bd09bSSatish Balay       }
113827bd09bSSatish Balay     }
114827bd09bSSatish Balay 
115b1c944f5SJed Brown     mask=PCTFS_floor_num_nodes>>1;
116db4deed7SKarl Rupp     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
1172fa5cd67SKarl Rupp       if (PCTFS_my_id%mask) continue;
118827bd09bSSatish Balay 
119b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
1202fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
1219566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD));
1222fa5cd67SKarl Rupp       } else {
1239566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status));
124827bd09bSSatish Balay       }
125827bd09bSSatish Balay     }
126827bd09bSSatish Balay   }
127827bd09bSSatish Balay 
128827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
129db4deed7SKarl Rupp   if (edge_not_pow_2) {
1302fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
1319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status));
1322fa5cd67SKarl Rupp     } else {
1339566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD));
134827bd09bSSatish Balay     }
135827bd09bSSatish Balay   }
1363fdc5746SBarry Smith   PetscFunctionReturn(0);
137827bd09bSSatish Balay }
138827bd09bSSatish Balay 
1397b1ae94cSBarry Smith /***********************************comm.c*************************************/
140b1c944f5SJed Brown PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
141827bd09bSSatish Balay {
1423fdc5746SBarry Smith   PetscInt       mask, edge;
1433fdc5746SBarry Smith   PetscInt       type, dest;
144827bd09bSSatish Balay   vfp            fp;
145827bd09bSSatish Balay   MPI_Status     status;
146827bd09bSSatish Balay 
1473fdc5746SBarry Smith   PetscFunctionBegin;
148827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
149*7827d75bSBarry Smith   PetscCheck(vals && work && oprs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
150827bd09bSSatish Balay 
151827bd09bSSatish Balay   /* non-uniform should have at least two entries */
1522c71b3e2SJacob Faibussowitsch   PetscCheckFalse((oprs[0] == NON_UNIFORM)&&(n<2),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?");
153827bd09bSSatish Balay 
154827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
1552fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
156827bd09bSSatish Balay 
157827bd09bSSatish Balay   /* if there's nothing to do return */
1582fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0);
159827bd09bSSatish Balay 
160827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
1612c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n<0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);
162827bd09bSSatish Balay 
163827bd09bSSatish Balay   /* advance to list of n operations for custom */
1642fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
165827bd09bSSatish Balay 
166*7827d75bSBarry Smith   PetscCheck(fp = (vfp) PCTFS_rvec_fct_addr(type),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: Could not retrieve function pointer!");
167827bd09bSSatish Balay 
168827bd09bSSatish Balay   /* all msgs will be of the same length */
169827bd09bSSatish Balay   /* if not a hypercube must colapse partial dim */
1702fa5cd67SKarl Rupp   if (edge_not_pow_2) {
1712fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
1729566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD));
1732fa5cd67SKarl Rupp     } else {
1749566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status));
175827bd09bSSatish Balay       (*fp)(vals,work,n,oprs);
176827bd09bSSatish Balay     }
177827bd09bSSatish Balay   }
178827bd09bSSatish Balay 
179827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
180db4deed7SKarl Rupp   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
181db4deed7SKarl Rupp     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
182b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
1832fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
1849566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD));
1852fa5cd67SKarl Rupp       } else {
1869566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status));
187827bd09bSSatish Balay         (*fp)(vals, work, n, oprs);
188827bd09bSSatish Balay       }
189827bd09bSSatish Balay     }
190827bd09bSSatish Balay 
191b1c944f5SJed Brown     mask=PCTFS_floor_num_nodes>>1;
192db4deed7SKarl Rupp     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
1932fa5cd67SKarl Rupp       if (PCTFS_my_id%mask) continue;
194827bd09bSSatish Balay 
195b1c944f5SJed Brown       dest = PCTFS_my_id^mask;
1962fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
1979566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD));
1982fa5cd67SKarl Rupp       } else {
1999566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status));
200827bd09bSSatish Balay       }
201827bd09bSSatish Balay     }
202827bd09bSSatish Balay   }
203827bd09bSSatish Balay 
204827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
205db4deed7SKarl Rupp   if (edge_not_pow_2) {
206db4deed7SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
2079566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status));
208db4deed7SKarl Rupp     } else {
2099566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD));
210827bd09bSSatish Balay     }
211827bd09bSSatish Balay   }
2123fdc5746SBarry Smith   PetscFunctionReturn(0);
213827bd09bSSatish Balay }
214827bd09bSSatish Balay 
2157b1ae94cSBarry Smith /***********************************comm.c*************************************/
216b1c944f5SJed Brown PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
217827bd09bSSatish Balay {
2183fdc5746SBarry Smith   PetscInt       mask, edge;
2193fdc5746SBarry Smith   PetscInt       type, dest;
220827bd09bSSatish Balay   vfp            fp;
221827bd09bSSatish Balay   MPI_Status     status;
222827bd09bSSatish Balay 
2233fdc5746SBarry Smith   PetscFunctionBegin;
224827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
225*7827d75bSBarry Smith   PetscCheck(vals && work && oprs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
226827bd09bSSatish Balay 
227827bd09bSSatish Balay   /* non-uniform should have at least two entries */
2282c71b3e2SJacob Faibussowitsch   PetscCheckFalse((oprs[0] == NON_UNIFORM)&&(n<2),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?");
229827bd09bSSatish Balay 
230827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
2312fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
232827bd09bSSatish Balay 
233827bd09bSSatish Balay   /* if there's nothing to do return */
2342fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0);
235827bd09bSSatish Balay 
236827bd09bSSatish Balay   /* the error msg says it all!!! */
23728b400f6SJacob Faibussowitsch   PetscCheck(!modfl_num_nodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!");
238827bd09bSSatish Balay 
239827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
2402c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n<0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n);
241827bd09bSSatish Balay 
242827bd09bSSatish Balay   /* can't do more dimensions then exist */
243b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
244827bd09bSSatish Balay 
245827bd09bSSatish Balay   /* advance to list of n operations for custom */
2462fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
247827bd09bSSatish Balay 
248*7827d75bSBarry Smith   PetscCheck(fp = (vfp) PCTFS_rvec_fct_addr(type),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: Could not retrieve function pointer!");
249827bd09bSSatish Balay 
250db4deed7SKarl Rupp   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
251b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
2522fa5cd67SKarl Rupp     if (PCTFS_my_id > dest) {
2539566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD));
2542fa5cd67SKarl Rupp     } else {
2559566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status));
256827bd09bSSatish Balay       (*fp)(vals, work, n, oprs);
257827bd09bSSatish Balay     }
258827bd09bSSatish Balay   }
259827bd09bSSatish Balay 
2602fa5cd67SKarl Rupp   if (edge==dim) mask>>=1;
261db4deed7SKarl Rupp   else {
2622fa5cd67SKarl Rupp     while (++edge<dim) mask<<=1;
263db4deed7SKarl Rupp   }
264827bd09bSSatish Balay 
265db4deed7SKarl Rupp   for (edge=0; edge<dim; edge++,mask>>=1) {
2662fa5cd67SKarl Rupp     if (PCTFS_my_id%mask) continue;
267827bd09bSSatish Balay 
268b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
2692fa5cd67SKarl Rupp     if (PCTFS_my_id < dest) {
2709566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD));
2712fa5cd67SKarl Rupp     } else {
2729566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status));
273827bd09bSSatish Balay     }
274827bd09bSSatish Balay   }
2753fdc5746SBarry Smith   PetscFunctionReturn(0);
276827bd09bSSatish Balay }
277827bd09bSSatish Balay 
2787b1ae94cSBarry Smith /******************************************************************************/
279b1c944f5SJed Brown PetscErrorCode PCTFS_ssgl_radd(PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
280827bd09bSSatish Balay {
2813fdc5746SBarry Smith   PetscInt       edge, type, dest, mask;
2823fdc5746SBarry Smith   PetscInt       stage_n;
283827bd09bSSatish Balay   MPI_Status     status;
2840912c85aSBarry Smith   PetscMPIInt    *maxval,flg;
285827bd09bSSatish Balay 
2863fdc5746SBarry Smith   PetscFunctionBegin;
287827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
2882fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
289827bd09bSSatish Balay 
290827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
291827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
292db4deed7SKarl Rupp   for (mask=0, edge=0; edge<level; edge++, mask++) {
293827bd09bSSatish Balay     stage_n = (segs[level] - segs[edge]);
294db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
295827bd09bSSatish Balay       dest = edge_node[edge];
296b1c944f5SJed Brown       type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge);
2972fa5cd67SKarl Rupp       if (PCTFS_my_id>dest) {
2989566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD));
2992fa5cd67SKarl Rupp       } else {
300b1c944f5SJed Brown         type =  type - PCTFS_my_id + dest;
3019566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status));
302ca8e9878SJed Brown         PCTFS_rvec_add(vals+segs[edge], work, stage_n);
303827bd09bSSatish Balay       }
304827bd09bSSatish Balay     }
305827bd09bSSatish Balay     mask <<= 1;
306827bd09bSSatish Balay   }
307827bd09bSSatish Balay   mask>>=1;
308db4deed7SKarl Rupp   for (edge=0; edge<level; edge++) {
309827bd09bSSatish Balay     stage_n = (segs[level] - segs[level-1-edge]);
310db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
311827bd09bSSatish Balay       dest = edge_node[level-edge-1];
312b1c944f5SJed Brown       type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge);
3139566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_get_attr(MPI_COMM_WORLD,MPI_TAG_UB,&maxval,&flg));
31428b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_LIB,"MPI error: MPI_Comm_get_attr() is not returning a MPI_TAG_UB");
3152c71b3e2SJacob Faibussowitsch       PetscCheckFalse(*maxval <= type,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MPI_TAG_UB for your current MPI implementation is not large enough to use PCTFS");
3162fa5cd67SKarl Rupp       if (PCTFS_my_id<dest) {
3179566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD));
3182fa5cd67SKarl Rupp       } else {
319b1c944f5SJed Brown         type =  type - PCTFS_my_id + dest;
3209566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status));
321827bd09bSSatish Balay       }
322827bd09bSSatish Balay     }
323827bd09bSSatish Balay     mask >>= 1;
324827bd09bSSatish Balay   }
3253fdc5746SBarry Smith   PetscFunctionReturn(0);
326827bd09bSSatish Balay }
327827bd09bSSatish Balay 
3287b1ae94cSBarry Smith /***********************************comm.c*************************************/
329b1c944f5SJed Brown PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
330827bd09bSSatish Balay {
33152f87cdaSBarry Smith   PetscInt       mask, edge;
33252f87cdaSBarry Smith   PetscInt       type, dest;
333827bd09bSSatish Balay   vfp            fp;
334827bd09bSSatish Balay   MPI_Status     status;
335827bd09bSSatish Balay 
3363fdc5746SBarry Smith   PetscFunctionBegin;
337827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
338*7827d75bSBarry Smith   PetscCheck(vals && work && oprs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
339827bd09bSSatish Balay 
340827bd09bSSatish Balay   /* non-uniform should have at least two entries */
3412c71b3e2SJacob Faibussowitsch   PetscCheckFalse((oprs[0] == NON_UNIFORM)&&(n<2),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?");
342827bd09bSSatish Balay 
343827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
3442fa5cd67SKarl Rupp   if (!p_init) PCTFS_comm_init();
345827bd09bSSatish Balay 
346827bd09bSSatish Balay   /* if there's nothing to do return */
3472fa5cd67SKarl Rupp   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0);
348827bd09bSSatish Balay 
349827bd09bSSatish Balay   /* the error msg says it all!!! */
35028b400f6SJacob Faibussowitsch   PetscCheck(!modfl_num_nodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!");
351827bd09bSSatish Balay 
352827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
3532c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n<0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n);
354827bd09bSSatish Balay 
355827bd09bSSatish Balay   /* can't do more dimensions then exist */
356b1c944f5SJed Brown   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
357827bd09bSSatish Balay 
358827bd09bSSatish Balay   /* advance to list of n operations for custom */
3592fa5cd67SKarl Rupp   if ((type=oprs[0])==NON_UNIFORM) oprs++;
360827bd09bSSatish Balay 
361*7827d75bSBarry Smith   PetscCheck(fp = (vfp) PCTFS_ivec_fct_addr(type),PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: Could not retrieve function pointer!");
362827bd09bSSatish Balay 
363db4deed7SKarl Rupp   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
364b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
3652fa5cd67SKarl Rupp     if (PCTFS_my_id > dest) {
3669566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD));
3672fa5cd67SKarl Rupp     } else {
3689566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status));
369827bd09bSSatish Balay       (*fp)(vals, work, n, oprs);
370827bd09bSSatish Balay     }
371827bd09bSSatish Balay   }
372827bd09bSSatish Balay 
3732fa5cd67SKarl Rupp   if (edge==dim) mask>>=1;
3742fa5cd67SKarl Rupp   else {
3752fa5cd67SKarl Rupp     while (++edge<dim) mask<<=1;
3762fa5cd67SKarl Rupp   }
377827bd09bSSatish Balay 
378db4deed7SKarl Rupp   for (edge=0; edge<dim; edge++,mask>>=1) {
3792fa5cd67SKarl Rupp     if (PCTFS_my_id%mask) continue;
380827bd09bSSatish Balay 
381b1c944f5SJed Brown     dest = PCTFS_my_id^mask;
3822fa5cd67SKarl Rupp     if (PCTFS_my_id < dest) {
3839566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD));
3842fa5cd67SKarl Rupp     } else {
3859566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status));
386827bd09bSSatish Balay     }
387827bd09bSSatish Balay   }
3883fdc5746SBarry Smith   PetscFunctionReturn(0);
389827bd09bSSatish Balay }
390