1 2 /***********************************comm.c************************************* 3 4 Author: Henry M. Tufo III 5 6 e-mail: hmt@cs.brown.edu 7 8 snail-mail: 9 Division of Applied Mathematics 10 Brown University 11 Providence, RI 02912 12 13 Last Modification: 14 11.21.97 15 ***********************************comm.c*************************************/ 16 #include <../src/ksp/pc/impls/tfs/tfs.h> 17 18 /* global program control variables - explicitly exported */ 19 PetscMPIInt PCTFS_my_id = 0; 20 PetscMPIInt PCTFS_num_nodes = 1; 21 PetscMPIInt PCTFS_floor_num_nodes = 0; 22 PetscMPIInt PCTFS_i_log2_num_nodes = 0; 23 24 /* global program control variables */ 25 static PetscInt p_init = 0; 26 static PetscInt modfl_num_nodes; 27 static PetscInt edge_not_pow_2; 28 29 static PetscInt edge_node[sizeof(PetscInt)*32]; 30 31 /***********************************comm.c*************************************/ 32 PetscErrorCode PCTFS_comm_init(void) 33 { 34 PetscFunctionBegin; 35 if (p_init++) PetscFunctionReturn(0); 36 37 MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes); 38 MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id); 39 40 if (PCTFS_num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!"); 41 42 PCTFS_ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32); 43 44 PCTFS_floor_num_nodes = 1; 45 PCTFS_i_log2_num_nodes = modfl_num_nodes = 0; 46 while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) { 47 edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes; 48 PCTFS_floor_num_nodes <<= 1; 49 PCTFS_i_log2_num_nodes++; 50 } 51 52 PCTFS_i_log2_num_nodes--; 53 PCTFS_floor_num_nodes >>= 1; 54 modfl_num_nodes = (PCTFS_num_nodes - PCTFS_floor_num_nodes); 55 56 if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1); 57 else if (PCTFS_my_id >= PCTFS_floor_num_nodes) edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1); 58 else edge_not_pow_2 = 0; 59 PetscFunctionReturn(0); 60 } 61 62 /***********************************comm.c*************************************/ 63 PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs) 64 { 65 PetscInt mask, edge; 66 PetscInt type, dest; 67 vfp fp; 68 MPI_Status status; 69 PetscInt ierr; 70 71 PetscFunctionBegin; 72 /* ok ... should have some data, work, and operator(s) */ 73 if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 74 75 /* non-uniform should have at least two entries */ 76 if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?"); 77 78 /* check to make sure comm package has been initialized */ 79 if (!p_init) PCTFS_comm_init(); 80 81 /* if there's nothing to do return */ 82 if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0); 83 84 /* a negative number if items to send ==> fatal */ 85 if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n); 86 87 /* advance to list of n operations for custom */ 88 if ((type=oprs[0])==NON_UNIFORM) oprs++; 89 90 /* major league hack */ 91 if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: Could not retrieve function pointer!\n"); 92 93 /* all msgs will be of the same length */ 94 /* if not a hypercube must colapse partial dim */ 95 if (edge_not_pow_2) { 96 if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 97 ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 98 } else { 99 ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 100 (*fp)(vals,work,n,oprs); 101 } 102 } 103 104 /* implement the mesh fan in/out exchange algorithm */ 105 if (PCTFS_my_id<PCTFS_floor_num_nodes) { 106 for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) { 107 dest = PCTFS_my_id^mask; 108 if (PCTFS_my_id > dest) { 109 ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 110 } else { 111 ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 112 (*fp)(vals, work, n, oprs); 113 } 114 } 115 116 mask=PCTFS_floor_num_nodes>>1; 117 for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) { 118 if (PCTFS_my_id%mask) continue; 119 120 dest = PCTFS_my_id^mask; 121 if (PCTFS_my_id < dest) { 122 ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 123 } else { 124 ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 125 } 126 } 127 } 128 129 /* if not a hypercube must expand to partial dim */ 130 if (edge_not_pow_2) { 131 if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 132 ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 133 } else { 134 ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 135 } 136 } 137 PetscFunctionReturn(0); 138 } 139 140 /***********************************comm.c*************************************/ 141 PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs) 142 { 143 PetscInt mask, edge; 144 PetscInt type, dest; 145 vfp fp; 146 MPI_Status status; 147 PetscErrorCode ierr; 148 149 PetscFunctionBegin; 150 /* ok ... should have some data, work, and operator(s) */ 151 if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 152 153 /* non-uniform should have at least two entries */ 154 if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?"); 155 156 /* check to make sure comm package has been initialized */ 157 if (!p_init) PCTFS_comm_init(); 158 159 /* if there's nothing to do return */ 160 if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0); 161 162 /* a negative number of items to send ==> fatal */ 163 if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n); 164 165 /* advance to list of n operations for custom */ 166 if ((type=oprs[0])==NON_UNIFORM) oprs++; 167 168 if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: Could not retrieve function pointer!\n"); 169 170 /* all msgs will be of the same length */ 171 /* if not a hypercube must colapse partial dim */ 172 if (edge_not_pow_2) { 173 if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 174 ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 175 } else { 176 ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 177 (*fp)(vals,work,n,oprs); 178 } 179 } 180 181 /* implement the mesh fan in/out exchange algorithm */ 182 if (PCTFS_my_id<PCTFS_floor_num_nodes) { 183 for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) { 184 dest = PCTFS_my_id^mask; 185 if (PCTFS_my_id > dest) { 186 ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 187 } else { 188 ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 189 (*fp)(vals, work, n, oprs); 190 } 191 } 192 193 mask=PCTFS_floor_num_nodes>>1; 194 for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) { 195 if (PCTFS_my_id%mask) continue; 196 197 dest = PCTFS_my_id^mask; 198 if (PCTFS_my_id < dest) { 199 ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 200 } else { 201 ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 202 } 203 } 204 } 205 206 /* if not a hypercube must expand to partial dim */ 207 if (edge_not_pow_2) { 208 if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 209 ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 210 } else { 211 ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 212 } 213 } 214 PetscFunctionReturn(0); 215 } 216 217 /***********************************comm.c*************************************/ 218 PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim) 219 { 220 PetscInt mask, edge; 221 PetscInt type, dest; 222 vfp fp; 223 MPI_Status status; 224 PetscErrorCode ierr; 225 226 PetscFunctionBegin; 227 /* ok ... should have some data, work, and operator(s) */ 228 if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 229 230 /* non-uniform should have at least two entries */ 231 if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?"); 232 233 /* check to make sure comm package has been initialized */ 234 if (!p_init) PCTFS_comm_init(); 235 236 /* if there's nothing to do return */ 237 if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0); 238 239 /* the error msg says it all!!! */ 240 if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!"); 241 242 /* a negative number of items to send ==> fatal */ 243 if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n); 244 245 /* can't do more dimensions then exist */ 246 dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 247 248 /* advance to list of n operations for custom */ 249 if ((type=oprs[0])==NON_UNIFORM) oprs++; 250 251 if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: Could not retrieve function pointer!\n"); 252 253 for (mask=1,edge=0; edge<dim; edge++,mask<<=1) { 254 dest = PCTFS_my_id^mask; 255 if (PCTFS_my_id > dest) { 256 ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 257 } else { 258 ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 259 (*fp)(vals, work, n, oprs); 260 } 261 } 262 263 if (edge==dim) mask>>=1; 264 else { 265 while (++edge<dim) mask<<=1; 266 } 267 268 for (edge=0; edge<dim; edge++,mask>>=1) { 269 if (PCTFS_my_id%mask) continue; 270 271 dest = PCTFS_my_id^mask; 272 if (PCTFS_my_id < dest) { 273 ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 274 } else { 275 ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 276 } 277 } 278 PetscFunctionReturn(0); 279 } 280 281 /******************************************************************************/ 282 PetscErrorCode PCTFS_ssgl_radd(PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs) 283 { 284 PetscInt edge, type, dest, mask; 285 PetscInt stage_n; 286 MPI_Status status; 287 PetscErrorCode ierr; 288 289 PetscFunctionBegin; 290 /* check to make sure comm package has been initialized */ 291 if (!p_init) PCTFS_comm_init(); 292 293 /* all msgs are *NOT* the same length */ 294 /* implement the mesh fan in/out exchange algorithm */ 295 for (mask=0, edge=0; edge<level; edge++, mask++) { 296 stage_n = (segs[level] - segs[edge]); 297 if (stage_n && !(PCTFS_my_id & mask)) { 298 dest = edge_node[edge]; 299 type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge); 300 if (PCTFS_my_id>dest) { 301 ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRMPI(ierr); 302 } else { 303 type = type - PCTFS_my_id + dest; 304 ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 305 PCTFS_rvec_add(vals+segs[edge], work, stage_n); 306 } 307 } 308 mask <<= 1; 309 } 310 mask>>=1; 311 for (edge=0; edge<level; edge++) { 312 stage_n = (segs[level] - segs[level-1-edge]); 313 if (stage_n && !(PCTFS_my_id & mask)) { 314 dest = edge_node[level-edge-1]; 315 type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge); 316 if (PCTFS_my_id<dest) { 317 ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRMPI(ierr); 318 } else { 319 type = type - PCTFS_my_id + dest; 320 ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 321 } 322 } 323 mask >>= 1; 324 } 325 PetscFunctionReturn(0); 326 } 327 328 /***********************************comm.c*************************************/ 329 PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim) 330 { 331 PetscInt mask, edge; 332 PetscInt type, dest; 333 vfp fp; 334 MPI_Status status; 335 PetscErrorCode ierr; 336 337 PetscFunctionBegin; 338 /* ok ... should have some data, work, and operator(s) */ 339 if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 340 341 /* non-uniform should have at least two entries */ 342 if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?"); 343 344 /* check to make sure comm package has been initialized */ 345 if (!p_init) PCTFS_comm_init(); 346 347 /* if there's nothing to do return */ 348 if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0); 349 350 /* the error msg says it all!!! */ 351 if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!"); 352 353 /* a negative number of items to send ==> fatal */ 354 if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n); 355 356 /* can't do more dimensions then exist */ 357 dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 358 359 /* advance to list of n operations for custom */ 360 if ((type=oprs[0])==NON_UNIFORM) oprs++; 361 362 if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: Could not retrieve function pointer!\n"); 363 364 for (mask=1,edge=0; edge<dim; edge++,mask<<=1) { 365 dest = PCTFS_my_id^mask; 366 if (PCTFS_my_id > dest) { 367 ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 368 } else { 369 ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 370 (*fp)(vals, work, n, oprs); 371 } 372 } 373 374 if (edge==dim) mask>>=1; 375 else { 376 while (++edge<dim) mask<<=1; 377 } 378 379 for (edge=0; edge<dim; edge++,mask>>=1) { 380 if (PCTFS_my_id%mask) continue; 381 382 dest = PCTFS_my_id^mask; 383 if (PCTFS_my_id < dest) { 384 ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 385 } else { 386 ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 387 } 388 } 389 PetscFunctionReturn(0); 390 } 391