1827bd09bSSatish Balay 2827bd09bSSatish Balay /***********************************comm.c************************************* 3827bd09bSSatish Balay SPARSE GATHER-SCATTER PACKAGE: bss_malloc bss_malloc ivec error comm gs queue 4827bd09bSSatish Balay 5827bd09bSSatish Balay Author: Henry M. Tufo III 6827bd09bSSatish Balay 7827bd09bSSatish Balay e-mail: hmt@cs.brown.edu 8827bd09bSSatish Balay 9827bd09bSSatish Balay snail-mail: 10827bd09bSSatish Balay Division of Applied Mathematics 11827bd09bSSatish Balay Brown University 12827bd09bSSatish Balay Providence, RI 02912 13827bd09bSSatish Balay 14827bd09bSSatish Balay Last Modification: 15827bd09bSSatish Balay 11.21.97 16827bd09bSSatish Balay ***********************************comm.c*************************************/ 17827bd09bSSatish Balay 18827bd09bSSatish Balay /***********************************comm.c************************************* 19827bd09bSSatish Balay File Description: 20827bd09bSSatish Balay ----------------- 21827bd09bSSatish Balay 22827bd09bSSatish Balay ***********************************comm.c*************************************/ 23d890fc11SSatish Balay #include "petsc.h" 24827bd09bSSatish Balay 25827bd09bSSatish Balay #include "const.h" 26827bd09bSSatish Balay #include "types.h" 27827bd09bSSatish Balay #include "error.h" 28827bd09bSSatish Balay #include "ivec.h" 29827bd09bSSatish Balay #include "bss_malloc.h" 30827bd09bSSatish Balay #include "comm.h" 31827bd09bSSatish Balay 32827bd09bSSatish Balay 33827bd09bSSatish Balay /* global program control variables - explicitly exported */ 34827bd09bSSatish Balay int my_id = 0; 35827bd09bSSatish Balay int num_nodes = 1; 36827bd09bSSatish Balay int floor_num_nodes = 0; 37827bd09bSSatish Balay int i_log2_num_nodes = 0; 38827bd09bSSatish Balay 39827bd09bSSatish Balay /* global program control variables */ 40827bd09bSSatish Balay static int p_init = 0; 41827bd09bSSatish Balay static int modfl_num_nodes; 42827bd09bSSatish Balay static int edge_not_pow_2; 43827bd09bSSatish Balay 44827bd09bSSatish Balay static unsigned int edge_node[INT_LEN*32]; 45827bd09bSSatish Balay 46827bd09bSSatish Balay /***********************************comm.c************************************* 47827bd09bSSatish Balay Function: giop() 48827bd09bSSatish Balay 49827bd09bSSatish Balay Input : 50827bd09bSSatish Balay Output: 51827bd09bSSatish Balay Return: 52827bd09bSSatish Balay Description: 53827bd09bSSatish Balay ***********************************comm.c*************************************/ 54827bd09bSSatish Balay void 55827bd09bSSatish Balay comm_init (void) 56827bd09bSSatish Balay { 57827bd09bSSatish Balay #ifdef DEBUG 58827bd09bSSatish Balay error_msg_warning("c_init() :: start\n"); 59827bd09bSSatish Balay #endif 60827bd09bSSatish Balay 61827bd09bSSatish Balay if (p_init++) return; 62827bd09bSSatish Balay 63827bd09bSSatish Balay #if PETSC_SIZEOF_INT != 4 64827bd09bSSatish Balay error_msg_warning("Int != 4 Bytes!"); 65827bd09bSSatish Balay #endif 66827bd09bSSatish Balay 67827bd09bSSatish Balay 68827bd09bSSatish Balay MPI_Comm_size(MPI_COMM_WORLD,&num_nodes); 69827bd09bSSatish Balay MPI_Comm_rank(MPI_COMM_WORLD,&my_id); 70827bd09bSSatish Balay 71827bd09bSSatish Balay if (num_nodes> (INT_MAX >> 1)) 72827bd09bSSatish Balay {error_msg_fatal("Can't have more then MAX_INT/2 nodes!!!");} 73827bd09bSSatish Balay 74827bd09bSSatish Balay ivec_zero((int*)edge_node,INT_LEN*32); 75827bd09bSSatish Balay 76827bd09bSSatish Balay floor_num_nodes = 1; 77827bd09bSSatish Balay i_log2_num_nodes = modfl_num_nodes = 0; 78827bd09bSSatish Balay while (floor_num_nodes <= num_nodes) 79827bd09bSSatish Balay { 80827bd09bSSatish Balay edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes; 81827bd09bSSatish Balay floor_num_nodes <<= 1; 82827bd09bSSatish Balay i_log2_num_nodes++; 83827bd09bSSatish Balay } 84827bd09bSSatish Balay 85827bd09bSSatish Balay i_log2_num_nodes--; 86827bd09bSSatish Balay floor_num_nodes >>= 1; 87827bd09bSSatish Balay modfl_num_nodes = (num_nodes - floor_num_nodes); 88827bd09bSSatish Balay 89827bd09bSSatish Balay if ((my_id > 0) && (my_id <= modfl_num_nodes)) 90827bd09bSSatish Balay {edge_not_pow_2=((my_id|floor_num_nodes)-1);} 91827bd09bSSatish Balay else if (my_id >= floor_num_nodes) 92827bd09bSSatish Balay {edge_not_pow_2=((my_id^floor_num_nodes)+1); 93827bd09bSSatish Balay } 94827bd09bSSatish Balay else 95827bd09bSSatish Balay {edge_not_pow_2 = 0;} 96827bd09bSSatish Balay 97827bd09bSSatish Balay #ifdef DEBUG 98*77431f27SBarry Smith error_msg_warning("c_init() done :: my_id=%d, num_nodes=%D",my_id,num_nodes); 99827bd09bSSatish Balay #endif 100827bd09bSSatish Balay } 101827bd09bSSatish Balay 102827bd09bSSatish Balay 103827bd09bSSatish Balay 104827bd09bSSatish Balay /***********************************comm.c************************************* 105827bd09bSSatish Balay Function: giop() 106827bd09bSSatish Balay 107827bd09bSSatish Balay Input : 108827bd09bSSatish Balay Output: 109827bd09bSSatish Balay Return: 110827bd09bSSatish Balay Description: fan-in/out version 111827bd09bSSatish Balay ***********************************comm.c*************************************/ 112827bd09bSSatish Balay void 113827bd09bSSatish Balay giop(int *vals, int *work, int n, int *oprs) 114827bd09bSSatish Balay { 115827bd09bSSatish Balay register int mask, edge; 116827bd09bSSatish Balay int type, dest; 117827bd09bSSatish Balay vfp fp; 118827bd09bSSatish Balay MPI_Status status; 119827bd09bSSatish Balay 120827bd09bSSatish Balay 121827bd09bSSatish Balay #ifdef SAFE 122827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 123827bd09bSSatish Balay if (!vals||!work||!oprs) 124*77431f27SBarry Smith {error_msg_fatal("giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);} 125827bd09bSSatish Balay 126827bd09bSSatish Balay /* non-uniform should have at least two entries */ 127827bd09bSSatish Balay if ((oprs[0] == NON_UNIFORM)&&(n<2)) 128827bd09bSSatish Balay {error_msg_fatal("giop() :: non_uniform and n=0,1?");} 129827bd09bSSatish Balay 130827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 131827bd09bSSatish Balay if (!p_init) 132827bd09bSSatish Balay {comm_init();} 133827bd09bSSatish Balay #endif 134827bd09bSSatish Balay 135827bd09bSSatish Balay /* if there's nothing to do return */ 136827bd09bSSatish Balay if ((num_nodes<2)||(!n)) 137827bd09bSSatish Balay { 138827bd09bSSatish Balay #ifdef DEBUG 139*77431f27SBarry Smith error_msg_warning("giop() :: n=%D num_nodes=%d",n,num_nodes); 140827bd09bSSatish Balay #endif 141827bd09bSSatish Balay return; 142827bd09bSSatish Balay } 143827bd09bSSatish Balay 144827bd09bSSatish Balay /* a negative number if items to send ==> fatal */ 145827bd09bSSatish Balay if (n<0) 146*77431f27SBarry Smith {error_msg_fatal("giop() :: n=%D<0?",n);} 147827bd09bSSatish Balay 148827bd09bSSatish Balay /* advance to list of n operations for custom */ 149827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 150827bd09bSSatish Balay {oprs++;} 151827bd09bSSatish Balay 152827bd09bSSatish Balay /* major league hack */ 153d890fc11SSatish Balay if (!(fp = (vfp) ivec_fct_addr(type))) { 154827bd09bSSatish Balay error_msg_warning("giop() :: hope you passed in a rbfp!\n"); 155827bd09bSSatish Balay fp = (vfp) oprs; 156827bd09bSSatish Balay } 157827bd09bSSatish Balay 158827bd09bSSatish Balay /* all msgs will be of the same length */ 159827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 160827bd09bSSatish Balay if (edge_not_pow_2) 161827bd09bSSatish Balay { 162827bd09bSSatish Balay if (my_id >= floor_num_nodes) 163827bd09bSSatish Balay {MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);} 164827bd09bSSatish Balay else 165827bd09bSSatish Balay { 166827bd09bSSatish Balay MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, 167827bd09bSSatish Balay MPI_COMM_WORLD,&status); 168827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 169827bd09bSSatish Balay } 170827bd09bSSatish Balay } 171827bd09bSSatish Balay 172827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 173827bd09bSSatish Balay if (my_id<floor_num_nodes) 174827bd09bSSatish Balay { 175827bd09bSSatish Balay for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1) 176827bd09bSSatish Balay { 177827bd09bSSatish Balay dest = my_id^mask; 178827bd09bSSatish Balay if (my_id > dest) 179827bd09bSSatish Balay {MPI_Send(vals,n,MPI_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);} 180827bd09bSSatish Balay else 181827bd09bSSatish Balay { 182827bd09bSSatish Balay MPI_Recv(work,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG2+dest, 183827bd09bSSatish Balay MPI_COMM_WORLD, &status); 184827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 185827bd09bSSatish Balay } 186827bd09bSSatish Balay } 187827bd09bSSatish Balay 188827bd09bSSatish Balay mask=floor_num_nodes>>1; 189827bd09bSSatish Balay for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1) 190827bd09bSSatish Balay { 191827bd09bSSatish Balay if (my_id%mask) 192827bd09bSSatish Balay {continue;} 193827bd09bSSatish Balay 194827bd09bSSatish Balay dest = my_id^mask; 195827bd09bSSatish Balay if (my_id < dest) 196827bd09bSSatish Balay {MPI_Send(vals,n,MPI_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);} 197827bd09bSSatish Balay else 198827bd09bSSatish Balay { 199827bd09bSSatish Balay MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG4+dest, 200827bd09bSSatish Balay MPI_COMM_WORLD, &status); 201827bd09bSSatish Balay } 202827bd09bSSatish Balay } 203827bd09bSSatish Balay } 204827bd09bSSatish Balay 205827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 206827bd09bSSatish Balay if (edge_not_pow_2) 207827bd09bSSatish Balay { 208827bd09bSSatish Balay if (my_id >= floor_num_nodes) 209827bd09bSSatish Balay { 210827bd09bSSatish Balay MPI_Recv(vals,n,MPI_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, 211827bd09bSSatish Balay MPI_COMM_WORLD,&status); 212827bd09bSSatish Balay } 213827bd09bSSatish Balay else 214827bd09bSSatish Balay {MPI_Send(vals,n,MPI_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);} 215827bd09bSSatish Balay } 216827bd09bSSatish Balay } 217827bd09bSSatish Balay 218827bd09bSSatish Balay /***********************************comm.c************************************* 219827bd09bSSatish Balay Function: grop() 220827bd09bSSatish Balay 221827bd09bSSatish Balay Input : 222827bd09bSSatish Balay Output: 223827bd09bSSatish Balay Return: 224827bd09bSSatish Balay Description: fan-in/out version 225827bd09bSSatish Balay ***********************************comm.c*************************************/ 226827bd09bSSatish Balay void 227827bd09bSSatish Balay grop(REAL *vals, REAL *work, int n, int *oprs) 228827bd09bSSatish Balay { 229827bd09bSSatish Balay register int mask, edge; 230827bd09bSSatish Balay int type, dest; 231827bd09bSSatish Balay vfp fp; 232827bd09bSSatish Balay MPI_Status status; 233827bd09bSSatish Balay 234827bd09bSSatish Balay 235827bd09bSSatish Balay #ifdef SAFE 236827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 237827bd09bSSatish Balay if (!vals||!work||!oprs) 238*77431f27SBarry Smith {error_msg_fatal("grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);} 239827bd09bSSatish Balay 240827bd09bSSatish Balay /* non-uniform should have at least two entries */ 241827bd09bSSatish Balay if ((oprs[0] == NON_UNIFORM)&&(n<2)) 242827bd09bSSatish Balay {error_msg_fatal("grop() :: non_uniform and n=0,1?");} 243827bd09bSSatish Balay 244827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 245827bd09bSSatish Balay if (!p_init) 246827bd09bSSatish Balay {comm_init();} 247827bd09bSSatish Balay #endif 248827bd09bSSatish Balay 249827bd09bSSatish Balay /* if there's nothing to do return */ 250827bd09bSSatish Balay if ((num_nodes<2)||(!n)) 251827bd09bSSatish Balay {return;} 252827bd09bSSatish Balay 253827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 254827bd09bSSatish Balay if (n<0) 255*77431f27SBarry Smith {error_msg_fatal("gdop() :: n=%D<0?",n);} 256827bd09bSSatish Balay 257827bd09bSSatish Balay /* advance to list of n operations for custom */ 258827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 259827bd09bSSatish Balay {oprs++;} 260827bd09bSSatish Balay 261d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))) { 262827bd09bSSatish Balay error_msg_warning("grop() :: hope you passed in a rbfp!\n"); 263827bd09bSSatish Balay fp = (vfp) oprs; 264827bd09bSSatish Balay } 265827bd09bSSatish Balay 266827bd09bSSatish Balay /* all msgs will be of the same length */ 267827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 268827bd09bSSatish Balay if (edge_not_pow_2) 269827bd09bSSatish Balay { 270827bd09bSSatish Balay if (my_id >= floor_num_nodes) 271827bd09bSSatish Balay {MPI_Send(vals,n,REAL_TYPE,edge_not_pow_2,MSGTAG0+my_id, 272827bd09bSSatish Balay MPI_COMM_WORLD);} 273827bd09bSSatish Balay else 274827bd09bSSatish Balay { 275827bd09bSSatish Balay MPI_Recv(work,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, 276827bd09bSSatish Balay MPI_COMM_WORLD,&status); 277827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 278827bd09bSSatish Balay } 279827bd09bSSatish Balay } 280827bd09bSSatish Balay 281827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 282827bd09bSSatish Balay if (my_id<floor_num_nodes) 283827bd09bSSatish Balay { 284827bd09bSSatish Balay for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1) 285827bd09bSSatish Balay { 286827bd09bSSatish Balay dest = my_id^mask; 287827bd09bSSatish Balay if (my_id > dest) 288827bd09bSSatish Balay {MPI_Send(vals,n,REAL_TYPE,dest,MSGTAG2+my_id,MPI_COMM_WORLD);} 289827bd09bSSatish Balay else 290827bd09bSSatish Balay { 291827bd09bSSatish Balay MPI_Recv(work,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG2+dest, 292827bd09bSSatish Balay MPI_COMM_WORLD, &status); 293827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 294827bd09bSSatish Balay } 295827bd09bSSatish Balay } 296827bd09bSSatish Balay 297827bd09bSSatish Balay mask=floor_num_nodes>>1; 298827bd09bSSatish Balay for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1) 299827bd09bSSatish Balay { 300827bd09bSSatish Balay if (my_id%mask) 301827bd09bSSatish Balay {continue;} 302827bd09bSSatish Balay 303827bd09bSSatish Balay dest = my_id^mask; 304827bd09bSSatish Balay if (my_id < dest) 305827bd09bSSatish Balay {MPI_Send(vals,n,REAL_TYPE,dest,MSGTAG4+my_id,MPI_COMM_WORLD);} 306827bd09bSSatish Balay else 307827bd09bSSatish Balay { 308827bd09bSSatish Balay MPI_Recv(vals,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG4+dest, 309827bd09bSSatish Balay MPI_COMM_WORLD, &status); 310827bd09bSSatish Balay } 311827bd09bSSatish Balay } 312827bd09bSSatish Balay } 313827bd09bSSatish Balay 314827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 315827bd09bSSatish Balay if (edge_not_pow_2) 316827bd09bSSatish Balay { 317827bd09bSSatish Balay if (my_id >= floor_num_nodes) 318827bd09bSSatish Balay { 319827bd09bSSatish Balay MPI_Recv(vals,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, 320827bd09bSSatish Balay MPI_COMM_WORLD,&status); 321827bd09bSSatish Balay } 322827bd09bSSatish Balay else 323827bd09bSSatish Balay {MPI_Send(vals,n,REAL_TYPE,edge_not_pow_2,MSGTAG5+my_id, 324827bd09bSSatish Balay MPI_COMM_WORLD);} 325827bd09bSSatish Balay } 326827bd09bSSatish Balay } 327827bd09bSSatish Balay 328827bd09bSSatish Balay 329827bd09bSSatish Balay /***********************************comm.c************************************* 330827bd09bSSatish Balay Function: grop() 331827bd09bSSatish Balay 332827bd09bSSatish Balay Input : 333827bd09bSSatish Balay Output: 334827bd09bSSatish Balay Return: 335827bd09bSSatish Balay Description: fan-in/out version 336827bd09bSSatish Balay 337827bd09bSSatish Balay note good only for num_nodes=2^k!!! 338827bd09bSSatish Balay 339827bd09bSSatish Balay ***********************************comm.c*************************************/ 340827bd09bSSatish Balay void 341827bd09bSSatish Balay grop_hc(REAL *vals, REAL *work, int n, int *oprs, int dim) 342827bd09bSSatish Balay { 343827bd09bSSatish Balay register int mask, edge; 344827bd09bSSatish Balay int type, dest; 345827bd09bSSatish Balay vfp fp; 346827bd09bSSatish Balay MPI_Status status; 347827bd09bSSatish Balay 348827bd09bSSatish Balay #ifdef SAFE 349827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 350827bd09bSSatish Balay if (!vals||!work||!oprs) 351*77431f27SBarry Smith {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);} 352827bd09bSSatish Balay 353827bd09bSSatish Balay /* non-uniform should have at least two entries */ 354827bd09bSSatish Balay if ((oprs[0] == NON_UNIFORM)&&(n<2)) 355827bd09bSSatish Balay {error_msg_fatal("grop_hc() :: non_uniform and n=0,1?");} 356827bd09bSSatish Balay 357827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 358827bd09bSSatish Balay if (!p_init) 359827bd09bSSatish Balay {comm_init();} 360827bd09bSSatish Balay #endif 361827bd09bSSatish Balay 362827bd09bSSatish Balay /* if there's nothing to do return */ 363827bd09bSSatish Balay if ((num_nodes<2)||(!n)||(dim<=0)) 364827bd09bSSatish Balay {return;} 365827bd09bSSatish Balay 366827bd09bSSatish Balay /* the error msg says it all!!! */ 367827bd09bSSatish Balay if (modfl_num_nodes) 368827bd09bSSatish Balay {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");} 369827bd09bSSatish Balay 370827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 371827bd09bSSatish Balay if (n<0) 372*77431f27SBarry Smith {error_msg_fatal("grop_hc() :: n=%D<0?",n);} 373827bd09bSSatish Balay 374827bd09bSSatish Balay /* can't do more dimensions then exist */ 375827bd09bSSatish Balay dim = MIN(dim,i_log2_num_nodes); 376827bd09bSSatish Balay 377827bd09bSSatish Balay /* advance to list of n operations for custom */ 378827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 379827bd09bSSatish Balay {oprs++;} 380827bd09bSSatish Balay 381d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))) { 382827bd09bSSatish Balay error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n"); 383827bd09bSSatish Balay fp = (vfp) oprs; 384827bd09bSSatish Balay } 385827bd09bSSatish Balay 386827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 387827bd09bSSatish Balay { 388827bd09bSSatish Balay dest = my_id^mask; 389827bd09bSSatish Balay if (my_id > dest) 390827bd09bSSatish Balay {MPI_Send(vals,n,REAL_TYPE,dest,MSGTAG2+my_id,MPI_COMM_WORLD);} 391827bd09bSSatish Balay else 392827bd09bSSatish Balay { 393827bd09bSSatish Balay MPI_Recv(work,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, 394827bd09bSSatish Balay &status); 395827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 396827bd09bSSatish Balay } 397827bd09bSSatish Balay } 398827bd09bSSatish Balay 399827bd09bSSatish Balay if (edge==dim) 400827bd09bSSatish Balay {mask>>=1;} 401827bd09bSSatish Balay else 402827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 403827bd09bSSatish Balay 404827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 405827bd09bSSatish Balay { 406827bd09bSSatish Balay if (my_id%mask) 407827bd09bSSatish Balay {continue;} 408827bd09bSSatish Balay 409827bd09bSSatish Balay dest = my_id^mask; 410827bd09bSSatish Balay if (my_id < dest) 411827bd09bSSatish Balay {MPI_Send(vals,n,REAL_TYPE,dest,MSGTAG4+my_id,MPI_COMM_WORLD);} 412827bd09bSSatish Balay else 413827bd09bSSatish Balay { 414827bd09bSSatish Balay MPI_Recv(vals,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, 415827bd09bSSatish Balay &status); 416827bd09bSSatish Balay } 417827bd09bSSatish Balay } 418827bd09bSSatish Balay } 419827bd09bSSatish Balay 420827bd09bSSatish Balay 421827bd09bSSatish Balay /***********************************comm.c************************************* 422827bd09bSSatish Balay Function: gop() 423827bd09bSSatish Balay 424827bd09bSSatish Balay Input : 425827bd09bSSatish Balay Output: 426827bd09bSSatish Balay Return: 427827bd09bSSatish Balay Description: fan-in/out version 428827bd09bSSatish Balay ***********************************comm.c*************************************/ 429c700b191SBarry Smith void gfop(void *vals, void *work, int n, vbfp fp, DATA_TYPE dt, int comm_type) 430827bd09bSSatish Balay { 431827bd09bSSatish Balay register int mask, edge; 432827bd09bSSatish Balay int dest; 433827bd09bSSatish Balay MPI_Status status; 434827bd09bSSatish Balay MPI_Op op; 435827bd09bSSatish Balay 436827bd09bSSatish Balay #ifdef SAFE 437827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 438827bd09bSSatish Balay if (!p_init) 439827bd09bSSatish Balay {comm_init();} 440827bd09bSSatish Balay 441827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 442827bd09bSSatish Balay if (!vals||!work||!fp) 443*77431f27SBarry Smith {error_msg_fatal("gop() :: v=%D, w=%D, f=%D",vals,work,fp);} 444827bd09bSSatish Balay #endif 445827bd09bSSatish Balay 446827bd09bSSatish Balay /* if there's nothing to do return */ 447827bd09bSSatish Balay if ((num_nodes<2)||(!n)) 448827bd09bSSatish Balay {return;} 449827bd09bSSatish Balay 450827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 451827bd09bSSatish Balay if (n<0) 452*77431f27SBarry Smith {error_msg_fatal("gop() :: n=%D<0?",n);} 453827bd09bSSatish Balay 454827bd09bSSatish Balay if (comm_type==MPI) 455827bd09bSSatish Balay { 456827bd09bSSatish Balay MPI_Op_create(fp,TRUE,&op); 457827bd09bSSatish Balay MPI_Allreduce (vals, work, n, dt, op, MPI_COMM_WORLD); 458827bd09bSSatish Balay MPI_Op_free(&op); 459827bd09bSSatish Balay return; 460827bd09bSSatish Balay } 461827bd09bSSatish Balay 462827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 463827bd09bSSatish Balay if (edge_not_pow_2) 464827bd09bSSatish Balay { 465827bd09bSSatish Balay if (my_id >= floor_num_nodes) 466827bd09bSSatish Balay {MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG0+my_id, 467827bd09bSSatish Balay MPI_COMM_WORLD);} 468827bd09bSSatish Balay else 469827bd09bSSatish Balay { 470827bd09bSSatish Balay MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, 471827bd09bSSatish Balay MPI_COMM_WORLD,&status); 472827bd09bSSatish Balay (*fp)(vals,work,&n,&dt); 473827bd09bSSatish Balay } 474827bd09bSSatish Balay } 475827bd09bSSatish Balay 476827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 477827bd09bSSatish Balay if (my_id<floor_num_nodes) 478827bd09bSSatish Balay { 479827bd09bSSatish Balay for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1) 480827bd09bSSatish Balay { 481827bd09bSSatish Balay dest = my_id^mask; 482827bd09bSSatish Balay if (my_id > dest) 483827bd09bSSatish Balay {MPI_Send(vals,n,dt,dest,MSGTAG2+my_id,MPI_COMM_WORLD);} 484827bd09bSSatish Balay else 485827bd09bSSatish Balay { 486827bd09bSSatish Balay MPI_Recv(work,n,dt,MPI_ANY_SOURCE,MSGTAG2+dest, 487827bd09bSSatish Balay MPI_COMM_WORLD, &status); 488827bd09bSSatish Balay (*fp)(vals, work, &n, &dt); 489827bd09bSSatish Balay } 490827bd09bSSatish Balay } 491827bd09bSSatish Balay 492827bd09bSSatish Balay mask=floor_num_nodes>>1; 493827bd09bSSatish Balay for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1) 494827bd09bSSatish Balay { 495827bd09bSSatish Balay if (my_id%mask) 496827bd09bSSatish Balay {continue;} 497827bd09bSSatish Balay 498827bd09bSSatish Balay dest = my_id^mask; 499827bd09bSSatish Balay if (my_id < dest) 500827bd09bSSatish Balay {MPI_Send(vals,n,dt,dest,MSGTAG4+my_id,MPI_COMM_WORLD);} 501827bd09bSSatish Balay else 502827bd09bSSatish Balay { 503827bd09bSSatish Balay MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG4+dest, 504827bd09bSSatish Balay MPI_COMM_WORLD, &status); 505827bd09bSSatish Balay } 506827bd09bSSatish Balay } 507827bd09bSSatish Balay } 508827bd09bSSatish Balay 509827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 510827bd09bSSatish Balay if (edge_not_pow_2) 511827bd09bSSatish Balay { 512827bd09bSSatish Balay if (my_id >= floor_num_nodes) 513827bd09bSSatish Balay { 514827bd09bSSatish Balay MPI_Recv(vals,n,dt,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, 515827bd09bSSatish Balay MPI_COMM_WORLD,&status); 516827bd09bSSatish Balay } 517827bd09bSSatish Balay else 518827bd09bSSatish Balay {MPI_Send(vals,n,dt,edge_not_pow_2,MSGTAG5+my_id, 519827bd09bSSatish Balay MPI_COMM_WORLD);} 520827bd09bSSatish Balay } 521827bd09bSSatish Balay } 522827bd09bSSatish Balay 523827bd09bSSatish Balay 524827bd09bSSatish Balay 525827bd09bSSatish Balay 526827bd09bSSatish Balay 527827bd09bSSatish Balay 528827bd09bSSatish Balay /****************************************************************************** 529827bd09bSSatish Balay Function: giop() 530827bd09bSSatish Balay 531827bd09bSSatish Balay Input : 532827bd09bSSatish Balay Output: 533827bd09bSSatish Balay Return: 534827bd09bSSatish Balay Description: 535827bd09bSSatish Balay 536827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii 537827bd09bSSatish Balay 538827bd09bSSatish Balay ******************************************************************************/ 539827bd09bSSatish Balay void 540827bd09bSSatish Balay ssgl_radd(register REAL *vals, register REAL *work, register int level, 541827bd09bSSatish Balay register int *segs) 542827bd09bSSatish Balay { 543827bd09bSSatish Balay register int edge, type, dest, mask; 544827bd09bSSatish Balay register int stage_n; 545827bd09bSSatish Balay MPI_Status status; 546827bd09bSSatish Balay 547827bd09bSSatish Balay #ifdef DEBUG 548827bd09bSSatish Balay if (level > i_log2_num_nodes) 549827bd09bSSatish Balay {error_msg_fatal("sgl_add() :: level > log_2(P)!!!");} 550827bd09bSSatish Balay #endif 551827bd09bSSatish Balay 552827bd09bSSatish Balay #ifdef SAFE 553827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 554827bd09bSSatish Balay if (!p_init) 555827bd09bSSatish Balay {comm_init();} 556827bd09bSSatish Balay #endif 557827bd09bSSatish Balay 558827bd09bSSatish Balay 559827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 560827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 561827bd09bSSatish Balay for (mask=0, edge=0; edge<level; edge++, mask++) 562827bd09bSSatish Balay { 563827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 564827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 565827bd09bSSatish Balay { 566827bd09bSSatish Balay dest = edge_node[edge]; 567827bd09bSSatish Balay type = MSGTAG3 + my_id + (num_nodes*edge); 568827bd09bSSatish Balay if (my_id>dest) 569827bd09bSSatish Balay /* {csend(type, vals+segs[edge],stage_n*REAL_LEN,dest,0);} */ 570827bd09bSSatish Balay {MPI_Send(vals+segs[edge],stage_n,REAL_TYPE,dest,type, 571827bd09bSSatish Balay MPI_COMM_WORLD);} 572827bd09bSSatish Balay else 573827bd09bSSatish Balay { 574827bd09bSSatish Balay type = type - my_id + dest; 575827bd09bSSatish Balay /* crecv(type,work,stage_n*REAL_LEN); */ 576827bd09bSSatish Balay MPI_Recv(work,stage_n,REAL_TYPE,MPI_ANY_SOURCE,type, 577827bd09bSSatish Balay MPI_COMM_WORLD,&status); 578827bd09bSSatish Balay rvec_add(vals+segs[edge], work, stage_n); 579827bd09bSSatish Balay /* daxpy(vals+segs[edge], work, stage_n); */ 580827bd09bSSatish Balay } 581827bd09bSSatish Balay } 582827bd09bSSatish Balay mask <<= 1; 583827bd09bSSatish Balay } 584827bd09bSSatish Balay mask>>=1; 585827bd09bSSatish Balay for (edge=0; edge<level; edge++) 586827bd09bSSatish Balay { 587827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 588827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 589827bd09bSSatish Balay { 590827bd09bSSatish Balay dest = edge_node[level-edge-1]; 591827bd09bSSatish Balay type = MSGTAG6 + my_id + (num_nodes*edge); 592827bd09bSSatish Balay if (my_id<dest) 593827bd09bSSatish Balay /* {csend(type,vals+segs[level-1-edge],stage_n*REAL_LEN,dest,0);} */ 594827bd09bSSatish Balay {MPI_Send(vals+segs[level-1-edge],stage_n,REAL_TYPE,dest,type, 595827bd09bSSatish Balay MPI_COMM_WORLD);} 596827bd09bSSatish Balay else 597827bd09bSSatish Balay { 598827bd09bSSatish Balay type = type - my_id + dest; 599827bd09bSSatish Balay /* crecv(type,vals+segs[level-1-edge],stage_n*REAL_LEN); */ 600827bd09bSSatish Balay MPI_Recv(vals+segs[level-1-edge],stage_n,REAL_TYPE, 601827bd09bSSatish Balay MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status); 602827bd09bSSatish Balay } 603827bd09bSSatish Balay } 604827bd09bSSatish Balay mask >>= 1; 605827bd09bSSatish Balay } 606827bd09bSSatish Balay } 607827bd09bSSatish Balay 608827bd09bSSatish Balay 609827bd09bSSatish Balay 610827bd09bSSatish Balay /***********************************comm.c************************************* 611827bd09bSSatish Balay Function: grop_hc_vvl() 612827bd09bSSatish Balay 613827bd09bSSatish Balay Input : 614827bd09bSSatish Balay Output: 615827bd09bSSatish Balay Return: 616827bd09bSSatish Balay Description: fan-in/out version 617827bd09bSSatish Balay 618827bd09bSSatish Balay note good only for num_nodes=2^k!!! 619827bd09bSSatish Balay 620827bd09bSSatish Balay ***********************************comm.c*************************************/ 621827bd09bSSatish Balay void 622827bd09bSSatish Balay grop_hc_vvl(REAL *vals, REAL *work, int *segs, int *oprs, int dim) 623827bd09bSSatish Balay { 624827bd09bSSatish Balay register int mask, edge, n; 625827bd09bSSatish Balay int type, dest; 626827bd09bSSatish Balay vfp fp; 627827bd09bSSatish Balay MPI_Status status; 628827bd09bSSatish Balay 629827bd09bSSatish Balay error_msg_fatal("grop_hc_vvl() :: is not working!\n"); 630827bd09bSSatish Balay 631827bd09bSSatish Balay #ifdef SAFE 632827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 633827bd09bSSatish Balay if (!vals||!work||!oprs||!segs) 634*77431f27SBarry Smith {error_msg_fatal("grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);} 635827bd09bSSatish Balay 636827bd09bSSatish Balay /* non-uniform should have at least two entries */ 637827bd09bSSatish Balay #if defined(not_used) 638827bd09bSSatish Balay if ((oprs[0] == NON_UNIFORM)&&(n<2)) 639827bd09bSSatish Balay {error_msg_fatal("grop_hc() :: non_uniform and n=0,1?");} 640827bd09bSSatish Balay #endif 641827bd09bSSatish Balay 642827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 643827bd09bSSatish Balay if (!p_init) 644827bd09bSSatish Balay {comm_init();} 645827bd09bSSatish Balay #endif 646827bd09bSSatish Balay 647827bd09bSSatish Balay /* if there's nothing to do return */ 648827bd09bSSatish Balay if ((num_nodes<2)||(dim<=0)) 649827bd09bSSatish Balay {return;} 650827bd09bSSatish Balay 651827bd09bSSatish Balay /* the error msg says it all!!! */ 652827bd09bSSatish Balay if (modfl_num_nodes) 653827bd09bSSatish Balay {error_msg_fatal("grop_hc() :: num_nodes not a power of 2!?!");} 654827bd09bSSatish Balay 655827bd09bSSatish Balay /* can't do more dimensions then exist */ 656827bd09bSSatish Balay dim = MIN(dim,i_log2_num_nodes); 657827bd09bSSatish Balay 658827bd09bSSatish Balay /* advance to list of n operations for custom */ 659827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 660827bd09bSSatish Balay {oprs++;} 661827bd09bSSatish Balay 662d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))){ 663827bd09bSSatish Balay error_msg_warning("grop_hc() :: hope you passed in a rbfp!\n"); 664827bd09bSSatish Balay fp = (vfp) oprs; 665827bd09bSSatish Balay } 666827bd09bSSatish Balay 667827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 668827bd09bSSatish Balay { 669827bd09bSSatish Balay n = segs[dim]-segs[edge]; 670827bd09bSSatish Balay dest = my_id^mask; 671827bd09bSSatish Balay if (my_id > dest) 672827bd09bSSatish Balay {MPI_Send(vals+segs[edge],n,REAL_TYPE,dest,MSGTAG2+my_id,MPI_COMM_WORLD);} 673827bd09bSSatish Balay else 674827bd09bSSatish Balay { 675827bd09bSSatish Balay MPI_Recv(work,n,REAL_TYPE,MPI_ANY_SOURCE,MSGTAG2+dest, 676827bd09bSSatish Balay MPI_COMM_WORLD, &status); 677827bd09bSSatish Balay (*fp)(vals+segs[edge], work, n, oprs); 678827bd09bSSatish Balay } 679827bd09bSSatish Balay } 680827bd09bSSatish Balay 681827bd09bSSatish Balay if (edge==dim) 682827bd09bSSatish Balay {mask>>=1;} 683827bd09bSSatish Balay else 684827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 685827bd09bSSatish Balay 686827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 687827bd09bSSatish Balay { 688827bd09bSSatish Balay if (my_id%mask) 689827bd09bSSatish Balay {continue;} 690827bd09bSSatish Balay 691827bd09bSSatish Balay n = (segs[dim]-segs[dim-1-edge]); 692827bd09bSSatish Balay 693827bd09bSSatish Balay dest = my_id^mask; 694827bd09bSSatish Balay if (my_id < dest) 695827bd09bSSatish Balay {MPI_Send(vals+segs[dim-1-edge],n,REAL_TYPE,dest,MSGTAG4+my_id, 696827bd09bSSatish Balay MPI_COMM_WORLD);} 697827bd09bSSatish Balay else 698827bd09bSSatish Balay { 699827bd09bSSatish Balay MPI_Recv(vals+segs[dim-1-edge],n,REAL_TYPE,MPI_ANY_SOURCE, 700827bd09bSSatish Balay MSGTAG4+dest,MPI_COMM_WORLD, &status); 701827bd09bSSatish Balay } 702827bd09bSSatish Balay } 703827bd09bSSatish Balay } 704827bd09bSSatish Balay 705827bd09bSSatish Balay /****************************************************************************** 706827bd09bSSatish Balay Function: giop() 707827bd09bSSatish Balay 708827bd09bSSatish Balay Input : 709827bd09bSSatish Balay Output: 710827bd09bSSatish Balay Return: 711827bd09bSSatish Balay Description: 712827bd09bSSatish Balay 713827bd09bSSatish Balay ii+1 entries in seg :: 0 .. ii 714827bd09bSSatish Balay 715827bd09bSSatish Balay ******************************************************************************/ 716c700b191SBarry Smith void new_ssgl_radd(register REAL *vals, register REAL *work, register int level,register int *segs) 717827bd09bSSatish Balay { 718827bd09bSSatish Balay register int edge, type, dest, mask; 719827bd09bSSatish Balay register int stage_n; 720827bd09bSSatish Balay MPI_Status status; 721827bd09bSSatish Balay 722827bd09bSSatish Balay 723827bd09bSSatish Balay #ifdef DEBUG 724827bd09bSSatish Balay if (level > i_log2_num_nodes) 725827bd09bSSatish Balay {error_msg_fatal("sgl_add() :: level > log_2(P)!!!");} 726827bd09bSSatish Balay #endif 727827bd09bSSatish Balay 728827bd09bSSatish Balay #ifdef SAFE 729827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 730827bd09bSSatish Balay if (!p_init) 731827bd09bSSatish Balay {comm_init();} 732827bd09bSSatish Balay #endif 733827bd09bSSatish Balay 734827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 735827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 736827bd09bSSatish Balay for (mask=0, edge=0; edge<level; edge++, mask++) 737827bd09bSSatish Balay { 738827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 739827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 740827bd09bSSatish Balay { 741827bd09bSSatish Balay dest = edge_node[edge]; 742827bd09bSSatish Balay type = MSGTAG3 + my_id + (num_nodes*edge); 743827bd09bSSatish Balay if (my_id>dest) 744827bd09bSSatish Balay /* {csend(type, vals+segs[edge],stage_n*REAL_LEN,dest,0);} */ 745827bd09bSSatish Balay {MPI_Send(vals+segs[edge],stage_n,REAL_TYPE,dest,type, 746c700b191SBarry Smith MPI_COMM_WORLD);} 747827bd09bSSatish Balay else 748827bd09bSSatish Balay { 749827bd09bSSatish Balay type = type - my_id + dest; 750827bd09bSSatish Balay /* crecv(type,work,stage_n*REAL_LEN); */ 751827bd09bSSatish Balay MPI_Recv(work,stage_n,REAL_TYPE,MPI_ANY_SOURCE,type, 752c700b191SBarry Smith MPI_COMM_WORLD,&status); 753827bd09bSSatish Balay rvec_add(vals+segs[edge], work, stage_n); 754827bd09bSSatish Balay /* daxpy(vals+segs[edge], work, stage_n); */ 755827bd09bSSatish Balay } 756827bd09bSSatish Balay } 757827bd09bSSatish Balay mask <<= 1; 758827bd09bSSatish Balay } 759827bd09bSSatish Balay mask>>=1; 760827bd09bSSatish Balay for (edge=0; edge<level; edge++) 761827bd09bSSatish Balay { 762827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 763827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 764827bd09bSSatish Balay { 765827bd09bSSatish Balay dest = edge_node[level-edge-1]; 766827bd09bSSatish Balay type = MSGTAG6 + my_id + (num_nodes*edge); 767827bd09bSSatish Balay if (my_id<dest) 768827bd09bSSatish Balay /* {csend(type,vals+segs[level-1-edge],stage_n*REAL_LEN,dest,0);} */ 769827bd09bSSatish Balay {MPI_Send(vals+segs[level-1-edge],stage_n,REAL_TYPE,dest,type, 770c700b191SBarry Smith MPI_COMM_WORLD);} 771827bd09bSSatish Balay else 772827bd09bSSatish Balay { 773827bd09bSSatish Balay type = type - my_id + dest; 774827bd09bSSatish Balay /* crecv(type,vals+segs[level-1-edge],stage_n*REAL_LEN); */ 775827bd09bSSatish Balay MPI_Recv(vals+segs[level-1-edge],stage_n,REAL_TYPE, 776c700b191SBarry Smith MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status); 777827bd09bSSatish Balay } 778827bd09bSSatish Balay } 779827bd09bSSatish Balay mask >>= 1; 780827bd09bSSatish Balay } 781827bd09bSSatish Balay } 782827bd09bSSatish Balay 783827bd09bSSatish Balay 784827bd09bSSatish Balay 785827bd09bSSatish Balay /***********************************comm.c************************************* 786827bd09bSSatish Balay Function: giop() 787827bd09bSSatish Balay 788827bd09bSSatish Balay Input : 789827bd09bSSatish Balay Output: 790827bd09bSSatish Balay Return: 791827bd09bSSatish Balay Description: fan-in/out version 792827bd09bSSatish Balay 793827bd09bSSatish Balay note good only for num_nodes=2^k!!! 794827bd09bSSatish Balay 795827bd09bSSatish Balay ***********************************comm.c*************************************/ 796c700b191SBarry Smith void giop_hc(int *vals, int *work, int n, int *oprs, int dim) 797827bd09bSSatish Balay { 798827bd09bSSatish Balay register int mask, edge; 799827bd09bSSatish Balay int type, dest; 800827bd09bSSatish Balay vfp fp; 801827bd09bSSatish Balay MPI_Status status; 802827bd09bSSatish Balay 803827bd09bSSatish Balay #ifdef SAFE 804827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 805827bd09bSSatish Balay if (!vals||!work||!oprs) 806*77431f27SBarry Smith {error_msg_fatal("giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);} 807827bd09bSSatish Balay 808827bd09bSSatish Balay /* non-uniform should have at least two entries */ 809827bd09bSSatish Balay if ((oprs[0] == NON_UNIFORM)&&(n<2)) 810827bd09bSSatish Balay {error_msg_fatal("giop_hc() :: non_uniform and n=0,1?");} 811827bd09bSSatish Balay 812827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 813827bd09bSSatish Balay if (!p_init) 814827bd09bSSatish Balay {comm_init();} 815827bd09bSSatish Balay #endif 816827bd09bSSatish Balay 817827bd09bSSatish Balay /* if there's nothing to do return */ 818827bd09bSSatish Balay if ((num_nodes<2)||(!n)||(dim<=0)) 819827bd09bSSatish Balay {return;} 820827bd09bSSatish Balay 821827bd09bSSatish Balay /* the error msg says it all!!! */ 822827bd09bSSatish Balay if (modfl_num_nodes) 823827bd09bSSatish Balay {error_msg_fatal("giop_hc() :: num_nodes not a power of 2!?!");} 824827bd09bSSatish Balay 825827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 826827bd09bSSatish Balay if (n<0) 827*77431f27SBarry Smith {error_msg_fatal("giop_hc() :: n=%D<0?",n);} 828827bd09bSSatish Balay 829827bd09bSSatish Balay /* can't do more dimensions then exist */ 830827bd09bSSatish Balay dim = MIN(dim,i_log2_num_nodes); 831827bd09bSSatish Balay 832827bd09bSSatish Balay /* advance to list of n operations for custom */ 833827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 834827bd09bSSatish Balay {oprs++;} 835827bd09bSSatish Balay 836d890fc11SSatish Balay if (!(fp = (vfp) ivec_fct_addr(type))){ 837827bd09bSSatish Balay error_msg_warning("giop_hc() :: hope you passed in a rbfp!\n"); 838827bd09bSSatish Balay fp = (vfp) oprs; 839827bd09bSSatish Balay } 840827bd09bSSatish Balay 841827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 842827bd09bSSatish Balay { 843827bd09bSSatish Balay dest = my_id^mask; 844827bd09bSSatish Balay if (my_id > dest) 845827bd09bSSatish Balay {MPI_Send(vals,n,INT_TYPE,dest,MSGTAG2+my_id,MPI_COMM_WORLD);} 846827bd09bSSatish Balay else 847827bd09bSSatish Balay { 848827bd09bSSatish Balay MPI_Recv(work,n,INT_TYPE,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, 849827bd09bSSatish Balay &status); 850827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 851827bd09bSSatish Balay } 852827bd09bSSatish Balay } 853827bd09bSSatish Balay 854827bd09bSSatish Balay if (edge==dim) 855827bd09bSSatish Balay {mask>>=1;} 856827bd09bSSatish Balay else 857827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 858827bd09bSSatish Balay 859827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 860827bd09bSSatish Balay { 861827bd09bSSatish Balay if (my_id%mask) 862827bd09bSSatish Balay {continue;} 863827bd09bSSatish Balay 864827bd09bSSatish Balay dest = my_id^mask; 865827bd09bSSatish Balay if (my_id < dest) 866827bd09bSSatish Balay {MPI_Send(vals,n,INT_TYPE,dest,MSGTAG4+my_id,MPI_COMM_WORLD);} 867827bd09bSSatish Balay else 868827bd09bSSatish Balay { 869827bd09bSSatish Balay MPI_Recv(vals,n,INT_TYPE,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, 870827bd09bSSatish Balay &status); 871827bd09bSSatish Balay } 872827bd09bSSatish Balay } 873827bd09bSSatish Balay } 874