xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision 51b144c619aff302b570817d6f78637b8418d403)
1827bd09bSSatish Balay /***********************************comm.c*************************************
2827bd09bSSatish Balay 
3827bd09bSSatish Balay Author: Henry M. Tufo III
4827bd09bSSatish Balay 
5827bd09bSSatish Balay e-mail: hmt@cs.brown.edu
6827bd09bSSatish Balay 
7827bd09bSSatish Balay snail-mail:
8827bd09bSSatish Balay Division of Applied Mathematics
9827bd09bSSatish Balay Brown University
10827bd09bSSatish Balay Providence, RI 02912
11827bd09bSSatish Balay 
12827bd09bSSatish Balay Last Modification:
13827bd09bSSatish Balay 11.21.97
14827bd09bSSatish Balay ***********************************comm.c*************************************/
15c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h>
16827bd09bSSatish Balay 
17827bd09bSSatish Balay /* global program control variables - explicitly exported */
18b1c944f5SJed Brown PetscMPIInt PCTFS_my_id            = 0;
19b1c944f5SJed Brown PetscMPIInt PCTFS_num_nodes        = 1;
20b1c944f5SJed Brown PetscMPIInt PCTFS_floor_num_nodes  = 0;
21b1c944f5SJed Brown PetscMPIInt PCTFS_i_log2_num_nodes = 0;
22827bd09bSSatish Balay 
23827bd09bSSatish Balay /* global program control variables */
2452f87cdaSBarry Smith static PetscInt p_init = 0;
2552f87cdaSBarry Smith static PetscInt modfl_num_nodes;
2652f87cdaSBarry Smith static PetscInt edge_not_pow_2;
27827bd09bSSatish Balay 
2852f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt) * 32];
29827bd09bSSatish Balay 
307b1ae94cSBarry Smith /***********************************comm.c*************************************/
PCTFS_comm_init(void)31d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_comm_init(void)
32d71ae5a4SJacob Faibussowitsch {
33362febeeSStefano Zampini   PetscFunctionBegin;
343ba16761SJacob Faibussowitsch   if (p_init++) PetscFunctionReturn(PETSC_SUCCESS);
35827bd09bSSatish Balay 
363ba16761SJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &PCTFS_num_nodes));
373ba16761SJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &PCTFS_my_id));
38827bd09bSSatish Balay 
39467446fbSPierre Jolivet   PetscCheck(PCTFS_num_nodes <= (INT_MAX >> 1), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Can't have more than MAX_INT/2 nodes!!!");
40827bd09bSSatish Balay 
413ba16761SJacob Faibussowitsch   PetscCall(PCTFS_ivec_zero((PetscInt *)edge_node, sizeof(PetscInt) * 32));
42827bd09bSSatish Balay 
43b1c944f5SJed Brown   PCTFS_floor_num_nodes  = 1;
44b1c944f5SJed Brown   PCTFS_i_log2_num_nodes = modfl_num_nodes = 0;
45db4deed7SKarl Rupp   while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) {
46b1c944f5SJed Brown     edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes;
47b1c944f5SJed Brown     PCTFS_floor_num_nodes <<= 1;
48b1c944f5SJed Brown     PCTFS_i_log2_num_nodes++;
49827bd09bSSatish Balay   }
50827bd09bSSatish Balay 
51b1c944f5SJed Brown   PCTFS_i_log2_num_nodes--;
52b1c944f5SJed Brown   PCTFS_floor_num_nodes >>= 1;
53b1c944f5SJed Brown   modfl_num_nodes = (PCTFS_num_nodes - PCTFS_floor_num_nodes);
54827bd09bSSatish Balay 
552fa5cd67SKarl Rupp   if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) edge_not_pow_2 = ((PCTFS_my_id | PCTFS_floor_num_nodes) - 1);
562fa5cd67SKarl Rupp   else if (PCTFS_my_id >= PCTFS_floor_num_nodes) edge_not_pow_2 = ((PCTFS_my_id ^ PCTFS_floor_num_nodes) + 1);
572fa5cd67SKarl Rupp   else edge_not_pow_2 = 0;
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59827bd09bSSatish Balay }
60827bd09bSSatish Balay 
617b1ae94cSBarry Smith /***********************************comm.c*************************************/
PCTFS_giop(PetscInt * vals,PetscInt * work,PetscInt n,PetscInt * oprs)62d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
63d71ae5a4SJacob Faibussowitsch {
643fdc5746SBarry Smith   PetscInt   mask, edge;
653fdc5746SBarry Smith   PetscInt   type, dest;
66827bd09bSSatish Balay   vfp        fp;
67827bd09bSSatish Balay   MPI_Status status;
68827bd09bSSatish Balay 
693fdc5746SBarry Smith   PetscFunctionBegin;
70827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
71c75bfeddSPierre Jolivet   PetscCheck(vals && work && oprs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop() :: vals=%p, work=%p, oprs=%p", (void *)vals, (void *)work, (void *)oprs);
72827bd09bSSatish Balay 
73827bd09bSSatish Balay   /* non-uniform should have at least two entries */
7408401ef6SPierre Jolivet   PetscCheck(!(oprs[0] == NON_UNIFORM) || !(n < 2), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop() :: non_uniform and n=0,1?");
75827bd09bSSatish Balay 
76827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
773ba16761SJacob Faibussowitsch   if (!p_init) PetscCall(PCTFS_comm_init());
78827bd09bSSatish Balay 
79827bd09bSSatish Balay   /* if there's nothing to do return */
803ba16761SJacob Faibussowitsch   if ((PCTFS_num_nodes < 2) || (!n)) PetscFunctionReturn(PETSC_SUCCESS);
8171a0148aSBarry Smith 
82827bd09bSSatish Balay   /* a negative number if items to send ==> fatal */
8363a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop() :: n=%" PetscInt_FMT "<0?", n);
84827bd09bSSatish Balay 
85827bd09bSSatish Balay   /* advance to list of n operations for custom */
862fa5cd67SKarl Rupp   if ((type = oprs[0]) == NON_UNIFORM) oprs++;
87827bd09bSSatish Balay 
88827bd09bSSatish Balay   /* major league hack */
897827d75bSBarry Smith   PetscCheck(fp = (vfp)PCTFS_ivec_fct_addr(type), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop() :: Could not retrieve function pointer!");
90827bd09bSSatish Balay 
91827bd09bSSatish Balay   /* all msgs will be of the same length */
921cb3f120SPierre Jolivet   /* if not a hypercube must collapse partial dim */
93db4deed7SKarl Rupp   if (edge_not_pow_2) {
942fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
959566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_INT, edge_not_pow_2, MSGTAG0 + PCTFS_my_id, MPI_COMM_WORLD));
962fa5cd67SKarl Rupp     } else {
979566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work, n, MPIU_INT, MPI_ANY_SOURCE, MSGTAG0 + edge_not_pow_2, MPI_COMM_WORLD, &status));
983ba16761SJacob Faibussowitsch       PetscCall((*fp)(vals, work, n, oprs));
99827bd09bSSatish Balay     }
100827bd09bSSatish Balay   }
101827bd09bSSatish Balay 
102827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
103db4deed7SKarl Rupp   if (PCTFS_my_id < PCTFS_floor_num_nodes) {
104db4deed7SKarl Rupp     for (mask = 1, edge = 0; edge < PCTFS_i_log2_num_nodes; edge++, mask <<= 1) {
105b1c944f5SJed Brown       dest = PCTFS_my_id ^ mask;
1062fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
1079566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals, n, MPIU_INT, dest, MSGTAG2 + PCTFS_my_id, MPI_COMM_WORLD));
1082fa5cd67SKarl Rupp       } else {
1099566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(work, n, MPIU_INT, MPI_ANY_SOURCE, MSGTAG2 + dest, MPI_COMM_WORLD, &status));
1103ba16761SJacob Faibussowitsch         PetscCall((*fp)(vals, work, n, oprs));
111827bd09bSSatish Balay       }
112827bd09bSSatish Balay     }
113827bd09bSSatish Balay 
114b1c944f5SJed Brown     mask = PCTFS_floor_num_nodes >> 1;
115db4deed7SKarl Rupp     for (edge = 0; edge < PCTFS_i_log2_num_nodes; edge++, mask >>= 1) {
1162fa5cd67SKarl Rupp       if (PCTFS_my_id % mask) continue;
117827bd09bSSatish Balay 
118b1c944f5SJed Brown       dest = PCTFS_my_id ^ mask;
1192fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
1209566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals, n, MPIU_INT, dest, MSGTAG4 + PCTFS_my_id, MPI_COMM_WORLD));
1212fa5cd67SKarl Rupp       } else {
1229566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(vals, n, MPIU_INT, MPI_ANY_SOURCE, MSGTAG4 + dest, MPI_COMM_WORLD, &status));
123827bd09bSSatish Balay       }
124827bd09bSSatish Balay     }
125827bd09bSSatish Balay   }
126827bd09bSSatish Balay 
127827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
128db4deed7SKarl Rupp   if (edge_not_pow_2) {
1292fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
1309566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals, n, MPIU_INT, MPI_ANY_SOURCE, MSGTAG5 + edge_not_pow_2, MPI_COMM_WORLD, &status));
1312fa5cd67SKarl Rupp     } else {
1329566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_INT, edge_not_pow_2, MSGTAG5 + PCTFS_my_id, MPI_COMM_WORLD));
133827bd09bSSatish Balay     }
134827bd09bSSatish Balay   }
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136827bd09bSSatish Balay }
137827bd09bSSatish Balay 
1387b1ae94cSBarry Smith /***********************************comm.c*************************************/
PCTFS_grop(PetscScalar * vals,PetscScalar * work,PetscInt n,PetscInt * oprs)139d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
140d71ae5a4SJacob Faibussowitsch {
1413fdc5746SBarry Smith   PetscInt   mask, edge;
1423fdc5746SBarry Smith   PetscInt   type, dest;
143827bd09bSSatish Balay   vfp        fp;
144827bd09bSSatish Balay   MPI_Status status;
145827bd09bSSatish Balay 
1463fdc5746SBarry Smith   PetscFunctionBegin;
147827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
148c75bfeddSPierre Jolivet   PetscCheck(vals && work && oprs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop() :: vals=%p, work=%p, oprs=%p", (void *)vals, (void *)work, (void *)oprs);
149827bd09bSSatish Balay 
150827bd09bSSatish Balay   /* non-uniform should have at least two entries */
15108401ef6SPierre Jolivet   PetscCheck(!(oprs[0] == NON_UNIFORM) || !(n < 2), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop() :: non_uniform and n=0,1?");
152827bd09bSSatish Balay 
153827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
1543ba16761SJacob Faibussowitsch   if (!p_init) PetscCall(PCTFS_comm_init());
155827bd09bSSatish Balay 
156827bd09bSSatish Balay   /* if there's nothing to do return */
1573ba16761SJacob Faibussowitsch   if ((PCTFS_num_nodes < 2) || (!n)) PetscFunctionReturn(PETSC_SUCCESS);
158827bd09bSSatish Balay 
159827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
16063a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gdop() :: n=%" PetscInt_FMT "<0?", n);
161827bd09bSSatish Balay 
162827bd09bSSatish Balay   /* advance to list of n operations for custom */
1632fa5cd67SKarl Rupp   if ((type = oprs[0]) == NON_UNIFORM) oprs++;
164827bd09bSSatish Balay 
16557508eceSPierre Jolivet   PetscCheck(fp = (vfp)PCTFS_rvec_fct_addr(type), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop() :: Could not retrieve function pointer!");
166827bd09bSSatish Balay 
167827bd09bSSatish Balay   /* all msgs will be of the same length */
1681cb3f120SPierre Jolivet   /* if not a hypercube must collapse partial dim */
1692fa5cd67SKarl Rupp   if (edge_not_pow_2) {
1702fa5cd67SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
1719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_SCALAR, edge_not_pow_2, MSGTAG0 + PCTFS_my_id, MPI_COMM_WORLD));
1722fa5cd67SKarl Rupp     } else {
1739566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work, n, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG0 + edge_not_pow_2, MPI_COMM_WORLD, &status));
1743ba16761SJacob Faibussowitsch       PetscCall((*fp)(vals, work, n, oprs));
175827bd09bSSatish Balay     }
176827bd09bSSatish Balay   }
177827bd09bSSatish Balay 
178827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
179db4deed7SKarl Rupp   if (PCTFS_my_id < PCTFS_floor_num_nodes) {
180db4deed7SKarl Rupp     for (mask = 1, edge = 0; edge < PCTFS_i_log2_num_nodes; edge++, mask <<= 1) {
181b1c944f5SJed Brown       dest = PCTFS_my_id ^ mask;
1822fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
1839566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals, n, MPIU_SCALAR, dest, MSGTAG2 + PCTFS_my_id, MPI_COMM_WORLD));
1842fa5cd67SKarl Rupp       } else {
1859566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(work, n, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG2 + dest, MPI_COMM_WORLD, &status));
1863ba16761SJacob Faibussowitsch         PetscCall((*fp)(vals, work, n, oprs));
187827bd09bSSatish Balay       }
188827bd09bSSatish Balay     }
189827bd09bSSatish Balay 
190b1c944f5SJed Brown     mask = PCTFS_floor_num_nodes >> 1;
191db4deed7SKarl Rupp     for (edge = 0; edge < PCTFS_i_log2_num_nodes; edge++, mask >>= 1) {
1922fa5cd67SKarl Rupp       if (PCTFS_my_id % mask) continue;
193827bd09bSSatish Balay 
194b1c944f5SJed Brown       dest = PCTFS_my_id ^ mask;
1952fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
1969566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals, n, MPIU_SCALAR, dest, MSGTAG4 + PCTFS_my_id, MPI_COMM_WORLD));
1972fa5cd67SKarl Rupp       } else {
1989566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(vals, n, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG4 + dest, MPI_COMM_WORLD, &status));
199827bd09bSSatish Balay       }
200827bd09bSSatish Balay     }
201827bd09bSSatish Balay   }
202827bd09bSSatish Balay 
203827bd09bSSatish Balay   /* if not a hypercube must expand to partial dim */
204db4deed7SKarl Rupp   if (edge_not_pow_2) {
205db4deed7SKarl Rupp     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
2069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals, n, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG5 + edge_not_pow_2, MPI_COMM_WORLD, &status));
207db4deed7SKarl Rupp     } else {
2089566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_SCALAR, edge_not_pow_2, MSGTAG5 + PCTFS_my_id, MPI_COMM_WORLD));
209827bd09bSSatish Balay     }
210827bd09bSSatish Balay   }
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212827bd09bSSatish Balay }
213827bd09bSSatish Balay 
2147b1ae94cSBarry Smith /***********************************comm.c*************************************/
PCTFS_grop_hc(PetscScalar * vals,PetscScalar * work,PetscInt n,PetscInt * oprs,PetscInt dim)215d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
216d71ae5a4SJacob Faibussowitsch {
2173fdc5746SBarry Smith   PetscInt   mask, edge;
2183fdc5746SBarry Smith   PetscInt   type, dest;
219827bd09bSSatish Balay   vfp        fp;
220827bd09bSSatish Balay   MPI_Status status;
221827bd09bSSatish Balay 
2223fdc5746SBarry Smith   PetscFunctionBegin;
223827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
224c75bfeddSPierre Jolivet   PetscCheck(vals && work && oprs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop_hc() :: vals=%p, work=%p, oprs=%p", (void *)vals, (void *)work, (void *)oprs);
225827bd09bSSatish Balay 
226827bd09bSSatish Balay   /* non-uniform should have at least two entries */
22708401ef6SPierre Jolivet   PetscCheck(!(oprs[0] == NON_UNIFORM) || !(n < 2), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop_hc() :: non_uniform and n=0,1?");
228827bd09bSSatish Balay 
229827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
2303ba16761SJacob Faibussowitsch   if (!p_init) PetscCall(PCTFS_comm_init());
231827bd09bSSatish Balay 
232827bd09bSSatish Balay   /* if there's nothing to do return */
2333ba16761SJacob Faibussowitsch   if ((PCTFS_num_nodes < 2) || (!n) || (dim <= 0)) PetscFunctionReturn(PETSC_SUCCESS);
234827bd09bSSatish Balay 
235827bd09bSSatish Balay   /* the error msg says it all!!! */
23628b400f6SJacob Faibussowitsch   PetscCheck(!modfl_num_nodes, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!");
237827bd09bSSatish Balay 
238827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
23963a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop_hc() :: n=%" PetscInt_FMT "<0?", n);
240827bd09bSSatish Balay 
241827bd09bSSatish Balay   /* can't do more dimensions then exist */
242b1c944f5SJed Brown   dim = PetscMin(dim, PCTFS_i_log2_num_nodes);
243827bd09bSSatish Balay 
244827bd09bSSatish Balay   /* advance to list of n operations for custom */
2452fa5cd67SKarl Rupp   if ((type = oprs[0]) == NON_UNIFORM) oprs++;
246827bd09bSSatish Balay 
2477827d75bSBarry Smith   PetscCheck(fp = (vfp)PCTFS_rvec_fct_addr(type), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_grop_hc() :: Could not retrieve function pointer!");
248827bd09bSSatish Balay 
249db4deed7SKarl Rupp   for (mask = 1, edge = 0; edge < dim; edge++, mask <<= 1) {
250b1c944f5SJed Brown     dest = PCTFS_my_id ^ mask;
2512fa5cd67SKarl Rupp     if (PCTFS_my_id > dest) {
2529566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_SCALAR, dest, MSGTAG2 + PCTFS_my_id, MPI_COMM_WORLD));
2532fa5cd67SKarl Rupp     } else {
2549566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work, n, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG2 + dest, MPI_COMM_WORLD, &status));
2553ba16761SJacob Faibussowitsch       PetscCall((*fp)(vals, work, n, oprs));
256827bd09bSSatish Balay     }
257827bd09bSSatish Balay   }
258827bd09bSSatish Balay 
2592fa5cd67SKarl Rupp   if (edge == dim) mask >>= 1;
260db4deed7SKarl Rupp   else {
2612fa5cd67SKarl Rupp     while (++edge < dim) mask <<= 1;
262db4deed7SKarl Rupp   }
263827bd09bSSatish Balay 
264db4deed7SKarl Rupp   for (edge = 0; edge < dim; edge++, mask >>= 1) {
2652fa5cd67SKarl Rupp     if (PCTFS_my_id % mask) continue;
266827bd09bSSatish Balay 
267b1c944f5SJed Brown     dest = PCTFS_my_id ^ mask;
2682fa5cd67SKarl Rupp     if (PCTFS_my_id < dest) {
2699566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_SCALAR, dest, MSGTAG4 + PCTFS_my_id, MPI_COMM_WORLD));
2702fa5cd67SKarl Rupp     } else {
2719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals, n, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG4 + dest, MPI_COMM_WORLD, &status));
272827bd09bSSatish Balay     }
273827bd09bSSatish Balay   }
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275827bd09bSSatish Balay }
276827bd09bSSatish Balay 
2777b1ae94cSBarry Smith /******************************************************************************/
PCTFS_ssgl_radd(PetscScalar * vals,PetscScalar * work,PetscInt level,PetscInt * segs)278d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_ssgl_radd(PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs)
279d71ae5a4SJacob Faibussowitsch {
2803fdc5746SBarry Smith   PetscInt     edge, type, dest, mask;
2813fdc5746SBarry Smith   PetscInt     stage_n;
282827bd09bSSatish Balay   MPI_Status   status;
283*b8b5be36SMartin Diehl   PetscMPIInt *maxval, iflg;
284827bd09bSSatish Balay 
2853fdc5746SBarry Smith   PetscFunctionBegin;
286827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
2873ba16761SJacob Faibussowitsch   if (!p_init) PetscCall(PCTFS_comm_init());
288827bd09bSSatish Balay 
289827bd09bSSatish Balay   /* all msgs are *NOT* the same length */
290827bd09bSSatish Balay   /* implement the mesh fan in/out exchange algorithm */
291db4deed7SKarl Rupp   for (mask = 0, edge = 0; edge < level; edge++, mask++) {
292827bd09bSSatish Balay     stage_n = (segs[level] - segs[edge]);
293db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
294827bd09bSSatish Balay       dest = edge_node[edge];
295b1c944f5SJed Brown       type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes * edge);
2962fa5cd67SKarl Rupp       if (PCTFS_my_id > dest) {
2979566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals + segs[edge], stage_n, MPIU_SCALAR, dest, type, MPI_COMM_WORLD));
2982fa5cd67SKarl Rupp       } else {
299b1c944f5SJed Brown         type = type - PCTFS_my_id + dest;
3009566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(work, stage_n, MPIU_SCALAR, MPI_ANY_SOURCE, type, MPI_COMM_WORLD, &status));
3013ba16761SJacob Faibussowitsch         PetscCall(PCTFS_rvec_add(vals + segs[edge], work, stage_n));
302827bd09bSSatish Balay       }
303827bd09bSSatish Balay     }
304827bd09bSSatish Balay     mask <<= 1;
305827bd09bSSatish Balay   }
306827bd09bSSatish Balay   mask >>= 1;
307db4deed7SKarl Rupp   for (edge = 0; edge < level; edge++) {
308827bd09bSSatish Balay     stage_n = (segs[level] - segs[level - 1 - edge]);
309db4deed7SKarl Rupp     if (stage_n && !(PCTFS_my_id & mask)) {
310827bd09bSSatish Balay       dest = edge_node[level - edge - 1];
311b1c944f5SJed Brown       type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes * edge);
312*b8b5be36SMartin Diehl       PetscCallMPI(MPI_Comm_get_attr(MPI_COMM_WORLD, MPI_TAG_UB, &maxval, &iflg));
313*b8b5be36SMartin Diehl       PetscCheck(iflg, PETSC_COMM_SELF, PETSC_ERR_LIB, "MPI error: MPI_Comm_get_attr() is not returning a MPI_TAG_UB");
31408401ef6SPierre Jolivet       PetscCheck(*maxval > type, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPI_TAG_UB for your current MPI implementation is not large enough to use PCTFS");
3152fa5cd67SKarl Rupp       if (PCTFS_my_id < dest) {
3169566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Send(vals + segs[level - 1 - edge], stage_n, MPIU_SCALAR, dest, type, MPI_COMM_WORLD));
3172fa5cd67SKarl Rupp       } else {
318b1c944f5SJed Brown         type = type - PCTFS_my_id + dest;
3199566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Recv(vals + segs[level - 1 - edge], stage_n, MPIU_SCALAR, MPI_ANY_SOURCE, type, MPI_COMM_WORLD, &status));
320827bd09bSSatish Balay       }
321827bd09bSSatish Balay     }
322827bd09bSSatish Balay     mask >>= 1;
323827bd09bSSatish Balay   }
3243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
325827bd09bSSatish Balay }
326827bd09bSSatish Balay 
3277b1ae94cSBarry Smith /***********************************comm.c*************************************/
PCTFS_giop_hc(PetscInt * vals,PetscInt * work,PetscInt n,PetscInt * oprs,PetscInt dim)328d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
329d71ae5a4SJacob Faibussowitsch {
33052f87cdaSBarry Smith   PetscInt   mask, edge;
33152f87cdaSBarry Smith   PetscInt   type, dest;
332827bd09bSSatish Balay   vfp        fp;
333827bd09bSSatish Balay   MPI_Status status;
334827bd09bSSatish Balay 
3353fdc5746SBarry Smith   PetscFunctionBegin;
336827bd09bSSatish Balay   /* ok ... should have some data, work, and operator(s) */
337c75bfeddSPierre Jolivet   PetscCheck(vals && work && oprs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop_hc() :: vals=%p, work=%p, oprs=%p", (void *)vals, (void *)work, (void *)oprs);
338827bd09bSSatish Balay 
339827bd09bSSatish Balay   /* non-uniform should have at least two entries */
34008401ef6SPierre Jolivet   PetscCheck(!(oprs[0] == NON_UNIFORM) || !(n < 2), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop_hc() :: non_uniform and n=0,1?");
341827bd09bSSatish Balay 
342827bd09bSSatish Balay   /* check to make sure comm package has been initialized */
3433ba16761SJacob Faibussowitsch   if (!p_init) PetscCall(PCTFS_comm_init());
344827bd09bSSatish Balay 
345827bd09bSSatish Balay   /* if there's nothing to do return */
3463ba16761SJacob Faibussowitsch   if ((PCTFS_num_nodes < 2) || (!n) || (dim <= 0)) PetscFunctionReturn(PETSC_SUCCESS);
347827bd09bSSatish Balay 
348827bd09bSSatish Balay   /* the error msg says it all!!! */
34928b400f6SJacob Faibussowitsch   PetscCheck(!modfl_num_nodes, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!");
350827bd09bSSatish Balay 
351827bd09bSSatish Balay   /* a negative number of items to send ==> fatal */
35263a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop_hc() :: n=%" PetscInt_FMT "<0?", n);
353827bd09bSSatish Balay 
354827bd09bSSatish Balay   /* can't do more dimensions then exist */
355b1c944f5SJed Brown   dim = PetscMin(dim, PCTFS_i_log2_num_nodes);
356827bd09bSSatish Balay 
357827bd09bSSatish Balay   /* advance to list of n operations for custom */
3582fa5cd67SKarl Rupp   if ((type = oprs[0]) == NON_UNIFORM) oprs++;
359827bd09bSSatish Balay 
3607827d75bSBarry Smith   PetscCheck(fp = (vfp)PCTFS_ivec_fct_addr(type), PETSC_COMM_SELF, PETSC_ERR_PLIB, "PCTFS_giop_hc() :: Could not retrieve function pointer!");
361827bd09bSSatish Balay 
362db4deed7SKarl Rupp   for (mask = 1, edge = 0; edge < dim; edge++, mask <<= 1) {
363b1c944f5SJed Brown     dest = PCTFS_my_id ^ mask;
3642fa5cd67SKarl Rupp     if (PCTFS_my_id > dest) {
3659566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_INT, dest, MSGTAG2 + PCTFS_my_id, MPI_COMM_WORLD));
3662fa5cd67SKarl Rupp     } else {
3679566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(work, n, MPIU_INT, MPI_ANY_SOURCE, MSGTAG2 + dest, MPI_COMM_WORLD, &status));
3683ba16761SJacob Faibussowitsch       PetscCall((*fp)(vals, work, n, oprs));
369827bd09bSSatish Balay     }
370827bd09bSSatish Balay   }
371827bd09bSSatish Balay 
3722fa5cd67SKarl Rupp   if (edge == dim) mask >>= 1;
3732fa5cd67SKarl Rupp   else {
3742fa5cd67SKarl Rupp     while (++edge < dim) mask <<= 1;
3752fa5cd67SKarl Rupp   }
376827bd09bSSatish Balay 
377db4deed7SKarl Rupp   for (edge = 0; edge < dim; edge++, mask >>= 1) {
3782fa5cd67SKarl Rupp     if (PCTFS_my_id % mask) continue;
379827bd09bSSatish Balay 
380b1c944f5SJed Brown     dest = PCTFS_my_id ^ mask;
3812fa5cd67SKarl Rupp     if (PCTFS_my_id < dest) {
3829566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(vals, n, MPIU_INT, dest, MSGTAG4 + PCTFS_my_id, MPI_COMM_WORLD));
3832fa5cd67SKarl Rupp     } else {
3849566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(vals, n, MPIU_INT, MPI_ANY_SOURCE, MSGTAG4 + dest, MPI_COMM_WORLD, &status));
385827bd09bSSatish Balay     }
386827bd09bSSatish Balay   }
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388827bd09bSSatish Balay }
389