1 2 /***********************************gs.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 6.21.97 15 ************************************gs.c**************************************/ 16 17 /***********************************gs.c*************************************** 18 File Description: 19 ----------------- 20 21 ************************************gs.c**************************************/ 22 23 #include <../src/ksp/pc/impls/tfs/tfs.h> 24 25 /* default length of number of items via tree - doubles if exceeded */ 26 #define TREE_BUF_SZ 2048; 27 #define GS_VEC_SZ 1 28 29 30 31 /***********************************gs.c*************************************** 32 Type: struct gather_scatter_id 33 ------------------------------ 34 35 ************************************gs.c**************************************/ 36 typedef struct gather_scatter_id { 37 PetscInt id; 38 PetscInt nel_min; 39 PetscInt nel_max; 40 PetscInt nel_sum; 41 PetscInt negl; 42 PetscInt gl_max; 43 PetscInt gl_min; 44 PetscInt repeats; 45 PetscInt ordered; 46 PetscInt positive; 47 PetscScalar *vals; 48 49 /* bit mask info */ 50 PetscInt *my_proc_mask; 51 PetscInt mask_sz; 52 PetscInt *ngh_buf; 53 PetscInt ngh_buf_sz; 54 PetscInt *nghs; 55 PetscInt num_nghs; 56 PetscInt max_nghs; 57 PetscInt *pw_nghs; 58 PetscInt num_pw_nghs; 59 PetscInt *tree_nghs; 60 PetscInt num_tree_nghs; 61 62 PetscInt num_loads; 63 64 /* repeats == true -> local info */ 65 PetscInt nel; /* number of unique elememts */ 66 PetscInt *elms; /* of size nel */ 67 PetscInt nel_total; 68 PetscInt *local_elms; /* of size nel_total */ 69 PetscInt *companion; /* of size nel_total */ 70 71 /* local info */ 72 PetscInt num_local_total; 73 PetscInt local_strength; 74 PetscInt num_local; 75 PetscInt *num_local_reduce; 76 PetscInt **local_reduce; 77 PetscInt num_local_gop; 78 PetscInt *num_gop_local_reduce; 79 PetscInt **gop_local_reduce; 80 81 /* pairwise info */ 82 PetscInt level; 83 PetscInt num_pairs; 84 PetscInt max_pairs; 85 PetscInt loc_node_pairs; 86 PetscInt max_node_pairs; 87 PetscInt min_node_pairs; 88 PetscInt avg_node_pairs; 89 PetscInt *pair_list; 90 PetscInt *msg_sizes; 91 PetscInt **node_list; 92 PetscInt len_pw_list; 93 PetscInt *pw_elm_list; 94 PetscScalar *pw_vals; 95 96 MPI_Request *msg_ids_in; 97 MPI_Request *msg_ids_out; 98 99 PetscScalar *out; 100 PetscScalar *in; 101 PetscInt msg_total; 102 103 /* tree - crystal accumulator info */ 104 PetscInt max_left_over; 105 PetscInt *pre; 106 PetscInt *in_num; 107 PetscInt *out_num; 108 PetscInt **in_list; 109 PetscInt **out_list; 110 111 /* new tree work*/ 112 PetscInt tree_nel; 113 PetscInt *tree_elms; 114 PetscScalar *tree_buf; 115 PetscScalar *tree_work; 116 117 PetscInt tree_map_sz; 118 PetscInt *tree_map_in; 119 PetscInt *tree_map_out; 120 121 /* current memory status */ 122 PetscInt gl_bss_min; 123 PetscInt gl_perm_min; 124 125 /* max segment size for PCTFS_gs_gop_vec() */ 126 PetscInt vec_sz; 127 128 /* hack to make paul happy */ 129 MPI_Comm PCTFS_gs_comm; 130 131 } PCTFS_gs_id; 132 133 static PCTFS_gs_id *gsi_check_args(PetscInt *elms, PetscInt nel, PetscInt level); 134 static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs); 135 static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs); 136 static PetscErrorCode set_pairwise(PCTFS_gs_id *gs); 137 static PCTFS_gs_id *gsi_new(void); 138 static PetscErrorCode set_tree(PCTFS_gs_id *gs); 139 140 /* same for all but vector flavor */ 141 static PetscErrorCode PCTFS_gs_gop_local_out(PCTFS_gs_id *gs, PetscScalar *vals); 142 /* vector flavor */ 143 static PetscErrorCode PCTFS_gs_gop_vec_local_out(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 144 145 static PetscErrorCode PCTFS_gs_gop_vec_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step); 146 static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step); 147 static PetscErrorCode PCTFS_gs_gop_vec_local_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 148 static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 149 static PetscErrorCode PCTFS_gs_gop_vec_tree_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 150 151 152 static PetscErrorCode PCTFS_gs_gop_local_plus(PCTFS_gs_id *gs, PetscScalar *vals); 153 static PetscErrorCode PCTFS_gs_gop_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals); 154 155 static PetscErrorCode PCTFS_gs_gop_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim); 156 static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim); 157 static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim); 158 159 /* global vars */ 160 /* from comm.c module */ 161 162 static PetscInt num_gs_ids = 0; 163 164 /* should make this dynamic ... later */ 165 static PetscInt msg_buf =MAX_MSG_BUF; 166 static PetscInt vec_sz =GS_VEC_SZ; 167 static PetscInt *tree_buf =NULL; 168 static PetscInt tree_buf_sz=0; 169 static PetscInt ntree =0; 170 171 /***************************************************************************/ 172 PetscErrorCode PCTFS_gs_init_vec_sz(PetscInt size) 173 { 174 PetscFunctionBegin; 175 vec_sz = size; 176 PetscFunctionReturn(0); 177 } 178 179 /******************************************************************************/ 180 PetscErrorCode PCTFS_gs_init_msg_buf_sz(PetscInt buf_size) 181 { 182 PetscFunctionBegin; 183 msg_buf = buf_size; 184 PetscFunctionReturn(0); 185 } 186 187 /******************************************************************************/ 188 PCTFS_gs_id *PCTFS_gs_init(PetscInt *elms, PetscInt nel, PetscInt level) 189 { 190 PCTFS_gs_id *gs; 191 MPI_Group PCTFS_gs_group; 192 MPI_Comm PCTFS_gs_comm; 193 PetscErrorCode ierr; 194 195 PetscFunctionBegin; 196 /* ensure that communication package has been initialized */ 197 PCTFS_comm_init(); 198 199 200 /* determines if we have enough dynamic/semi-static memory */ 201 /* checks input, allocs and sets gd_id template */ 202 gs = gsi_check_args(elms,nel,level); 203 204 /* only bit mask version up and working for the moment */ 205 /* LATER :: get int list version working for sparse pblms */ 206 ierr = gsi_via_bit_mask(gs);CHKERRABORT(PETSC_COMM_WORLD,ierr); 207 208 209 ierr = MPI_Comm_group(MPI_COMM_WORLD,&PCTFS_gs_group);CHKERRABORT(PETSC_COMM_WORLD,ierr); 210 ierr = MPI_Comm_create(MPI_COMM_WORLD,PCTFS_gs_group,&PCTFS_gs_comm);CHKERRABORT(PETSC_COMM_WORLD,ierr); 211 212 gs->PCTFS_gs_comm=PCTFS_gs_comm; 213 214 return(gs); 215 } 216 217 /******************************************************************************/ 218 static PCTFS_gs_id *gsi_new(void) 219 { 220 PetscErrorCode ierr; 221 PCTFS_gs_id *gs; 222 gs = (PCTFS_gs_id*) malloc(sizeof(PCTFS_gs_id)); 223 ierr = PetscMemzero(gs,sizeof(PCTFS_gs_id));CHKERRABORT(PETSC_COMM_WORLD,ierr); 224 return(gs); 225 } 226 227 /******************************************************************************/ 228 static PCTFS_gs_id *gsi_check_args(PetscInt *in_elms, PetscInt nel, PetscInt level) 229 { 230 PetscInt i, j, k, t2; 231 PetscInt *companion, *elms, *unique, *iptr; 232 PetscInt num_local=0, *num_to_reduce, **local_reduce; 233 PetscInt oprs[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_MIN,GL_B_AND}; 234 PetscInt vals[sizeof(oprs)/sizeof(oprs[0])-1]; 235 PetscInt work[sizeof(oprs)/sizeof(oprs[0])-1]; 236 PCTFS_gs_id *gs; 237 PetscErrorCode ierr; 238 239 240 if (!in_elms) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"elms point to nothing!!!\n"); 241 if (nel<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"can't have fewer than 0 elms!!!\n"); 242 243 if (nel==0) { ierr = PetscInfo(0,"I don't have any elements!!!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); } 244 245 /* get space for gs template */ 246 gs = gsi_new(); 247 gs->id = ++num_gs_ids; 248 249 /* hmt 6.4.99 */ 250 /* caller can set global ids that don't participate to 0 */ 251 /* PCTFS_gs_init ignores all zeros in elm list */ 252 /* negative global ids are still invalid */ 253 for (i=j=0; i<nel; i++) { 254 if (in_elms[i]!=0) j++; 255 } 256 257 k=nel; nel=j; 258 259 /* copy over in_elms list and create inverse map */ 260 elms = (PetscInt*) malloc((nel+1)*sizeof(PetscInt)); 261 companion = (PetscInt*) malloc(nel*sizeof(PetscInt)); 262 263 for (i=j=0; i<k; i++) { 264 if (in_elms[i]!=0) { elms[j] = in_elms[i]; companion[j++] = i; } 265 } 266 267 if (j!=nel) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"nel j mismatch!\n"); 268 269 /* pre-pass ... check to see if sorted */ 270 elms[nel] = INT_MAX; 271 iptr = elms; 272 unique = elms+1; 273 j =0; 274 while (*iptr!=INT_MAX) { 275 if (*iptr++>*unique++) { j=1; break; } 276 } 277 278 /* set up inverse map */ 279 if (j) { 280 ierr = PetscInfo(0,"gsi_check_args() :: elm list *not* sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); 281 ierr = PCTFS_SMI_sort((void*)elms, (void*)companion, nel, SORT_INTEGER);CHKERRABORT(PETSC_COMM_WORLD,ierr); 282 } else { ierr = PetscInfo(0,"gsi_check_args() :: elm list sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); } 283 elms[nel] = INT_MIN; 284 285 /* first pass */ 286 /* determine number of unique elements, check pd */ 287 for (i=k=0; i<nel; i+=j) { 288 t2 = elms[i]; 289 j = ++i; 290 291 /* clump 'em for now */ 292 while (elms[j]==t2) j++; 293 294 /* how many together and num local */ 295 if (j-=i) { num_local++; k+=j; } 296 } 297 298 /* how many unique elements? */ 299 gs->repeats = k; 300 gs->nel = nel-k; 301 302 303 /* number of repeats? */ 304 gs->num_local = num_local; 305 num_local += 2; 306 gs->local_reduce = local_reduce=(PetscInt**)malloc(num_local*sizeof(PetscInt*)); 307 gs->num_local_reduce = num_to_reduce=(PetscInt*) malloc(num_local*sizeof(PetscInt)); 308 309 unique = (PetscInt*) malloc((gs->nel+1)*sizeof(PetscInt)); 310 gs->elms = unique; 311 gs->nel_total = nel; 312 gs->local_elms = elms; 313 gs->companion = companion; 314 315 /* compess map as well as keep track of local ops */ 316 for (num_local=i=j=0; i<gs->nel; i++) { 317 k = j; 318 t2 = unique[i] = elms[j]; 319 companion[i] = companion[j]; 320 321 while (elms[j]==t2) j++; 322 323 if ((t2=(j-k))>1) { 324 /* number together */ 325 num_to_reduce[num_local] = t2++; 326 327 iptr = local_reduce[num_local++] = (PetscInt*)malloc(t2*sizeof(PetscInt)); 328 329 /* to use binary searching don't remap until we check intersection */ 330 *iptr++ = i; 331 332 /* note that we're skipping the first one */ 333 while (++k<j) *(iptr++) = companion[k]; 334 *iptr = -1; 335 } 336 } 337 338 /* sentinel for ngh_buf */ 339 unique[gs->nel]=INT_MAX; 340 341 /* for two partition sort hack */ 342 num_to_reduce[num_local] = 0; 343 local_reduce[num_local] = NULL; 344 num_to_reduce[++num_local] = 0; 345 local_reduce[num_local] = NULL; 346 347 /* load 'em up */ 348 /* note one extra to hold NON_UNIFORM flag!!! */ 349 vals[2] = vals[1] = vals[0] = nel; 350 if (gs->nel>0) { 351 vals[3] = unique[0]; 352 vals[4] = unique[gs->nel-1]; 353 } else { 354 vals[3] = INT_MAX; 355 vals[4] = INT_MIN; 356 } 357 vals[5] = level; 358 vals[6] = num_gs_ids; 359 360 /* GLOBAL: send 'em out */ 361 ierr = PCTFS_giop(vals,work,sizeof(oprs)/sizeof(oprs[0])-1,oprs);CHKERRABORT(PETSC_COMM_WORLD,ierr); 362 363 /* must be semi-pos def - only pairwise depends on this */ 364 /* LATER - remove this restriction */ 365 if (vals[3]<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system not semi-pos def \n"); 366 if (vals[4]==INT_MAX) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system ub too large !\n"); 367 368 gs->nel_min = vals[0]; 369 gs->nel_max = vals[1]; 370 gs->nel_sum = vals[2]; 371 gs->gl_min = vals[3]; 372 gs->gl_max = vals[4]; 373 gs->negl = vals[4]-vals[3]+1; 374 375 if (gs->negl<=0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system empty or neg :: %d\n"); 376 377 /* LATER :: add level == -1 -> program selects level */ 378 if (vals[5]<0) vals[5]=0; 379 else if (vals[5]>PCTFS_num_nodes) vals[5]=PCTFS_num_nodes; 380 gs->level = vals[5]; 381 382 return(gs); 383 } 384 385 /******************************************************************************/ 386 static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs) 387 { 388 PetscInt i, nel, *elms; 389 PetscInt t1; 390 PetscInt **reduce; 391 PetscInt *map; 392 PetscErrorCode ierr; 393 394 PetscFunctionBegin; 395 /* totally local removes ... PCTFS_ct_bits == 0 */ 396 get_ngh_buf(gs); 397 398 if (gs->level) set_pairwise(gs); 399 if (gs->max_left_over) set_tree(gs); 400 401 /* intersection local and pairwise/tree? */ 402 gs->num_local_total = gs->num_local; 403 gs->gop_local_reduce = gs->local_reduce; 404 gs->num_gop_local_reduce = gs->num_local_reduce; 405 406 map = gs->companion; 407 408 /* is there any local compression */ 409 if (!gs->num_local) { 410 gs->local_strength = NONE; 411 gs->num_local_gop = 0; 412 } else { 413 /* ok find intersection */ 414 map = gs->companion; 415 reduce = gs->local_reduce; 416 for (i=0, t1=0; i<gs->num_local; i++, reduce++) { 417 if ((PCTFS_ivec_binary_search(**reduce,gs->pw_elm_list,gs->len_pw_list)>=0) || PCTFS_ivec_binary_search(**reduce,gs->tree_map_in,gs->tree_map_sz)>=0) { 418 t1++; 419 if (gs->num_local_reduce[i]<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nobody in list?"); 420 gs->num_local_reduce[i] *= -1; 421 } 422 **reduce=map[**reduce]; 423 } 424 425 /* intersection is empty */ 426 if (!t1) { 427 gs->local_strength = FULL; 428 gs->num_local_gop = 0; 429 } else { /* intersection not empty */ 430 gs->local_strength = PARTIAL; 431 432 ierr = PCTFS_SMI_sort((void*)gs->num_local_reduce, (void*)gs->local_reduce, gs->num_local + 1, SORT_INT_PTR);CHKERRQ(ierr); 433 434 gs->num_local_gop = t1; 435 gs->num_local_total = gs->num_local; 436 gs->num_local -= t1; 437 gs->gop_local_reduce = gs->local_reduce; 438 gs->num_gop_local_reduce = gs->num_local_reduce; 439 440 for (i=0; i<t1; i++) { 441 if (gs->num_gop_local_reduce[i]>=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"they aren't negative?"); 442 gs->num_gop_local_reduce[i] *= -1; 443 gs->local_reduce++; 444 gs->num_local_reduce++; 445 } 446 gs->local_reduce++; 447 gs->num_local_reduce++; 448 } 449 } 450 451 elms = gs->pw_elm_list; 452 nel = gs->len_pw_list; 453 for (i=0; i<nel; i++) elms[i] = map[elms[i]]; 454 455 elms = gs->tree_map_in; 456 nel = gs->tree_map_sz; 457 for (i=0; i<nel; i++) elms[i] = map[elms[i]]; 458 459 /* clean up */ 460 free((void*) gs->local_elms); 461 free((void*) gs->companion); 462 free((void*) gs->elms); 463 free((void*) gs->ngh_buf); 464 gs->local_elms = gs->companion = gs->elms = gs->ngh_buf = NULL; 465 PetscFunctionReturn(0); 466 } 467 468 /******************************************************************************/ 469 static PetscErrorCode place_in_tree(PetscInt elm) 470 { 471 PetscInt *tp, n; 472 473 PetscFunctionBegin; 474 if (ntree==tree_buf_sz) { 475 if (tree_buf_sz) { 476 tp = tree_buf; 477 n = tree_buf_sz; 478 tree_buf_sz<<=1; 479 tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt)); 480 PCTFS_ivec_copy(tree_buf,tp,n); 481 free(tp); 482 } else { 483 tree_buf_sz = TREE_BUF_SZ; 484 tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt)); 485 } 486 } 487 488 tree_buf[ntree++] = elm; 489 PetscFunctionReturn(0); 490 } 491 492 /******************************************************************************/ 493 static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs) 494 { 495 PetscInt i, j, npw=0, ntree_map=0; 496 PetscInt p_mask_size, ngh_buf_size, buf_size; 497 PetscInt *p_mask, *sh_proc_mask, *pw_sh_proc_mask; 498 PetscInt *ngh_buf, *buf1, *buf2; 499 PetscInt offset, per_load, num_loads, or_ct, start, end; 500 PetscInt *ptr1, *ptr2, i_start, negl, nel, *elms; 501 PetscInt oper=GL_B_OR; 502 PetscInt *ptr3, *t_mask, level, ct1, ct2; 503 PetscErrorCode ierr; 504 505 PetscFunctionBegin; 506 /* to make life easier */ 507 nel = gs->nel; 508 elms = gs->elms; 509 level = gs->level; 510 511 /* det #bytes needed for processor bit masks and init w/mask cor. to PCTFS_my_id */ 512 p_mask = (PetscInt*) malloc(p_mask_size=PCTFS_len_bit_mask(PCTFS_num_nodes)); 513 ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr); 514 515 /* allocate space for masks and info bufs */ 516 gs->nghs = sh_proc_mask = (PetscInt*) malloc(p_mask_size); 517 gs->pw_nghs = pw_sh_proc_mask = (PetscInt*) malloc(p_mask_size); 518 gs->ngh_buf_sz = ngh_buf_size = p_mask_size*nel; 519 t_mask = (PetscInt*) malloc(p_mask_size); 520 gs->ngh_buf = ngh_buf = (PetscInt*) malloc(ngh_buf_size); 521 522 /* comm buffer size ... memory usage bounded by ~2*msg_buf */ 523 /* had thought I could exploit rendezvous threshold */ 524 525 /* default is one pass */ 526 per_load = negl = gs->negl; 527 gs->num_loads = num_loads = 1; 528 i = p_mask_size*negl; 529 530 /* possible overflow on buffer size */ 531 /* overflow hack */ 532 if (i<0) i=INT_MAX; 533 534 buf_size = PetscMin(msg_buf,i); 535 536 /* can we do it? */ 537 if (p_mask_size>buf_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"get_ngh_buf() :: buf<pms :: %d>%d\n",p_mask_size,buf_size); 538 539 /* get PCTFS_giop buf space ... make *only* one malloc */ 540 buf1 = (PetscInt*) malloc(buf_size<<1); 541 542 /* more than one gior exchange needed? */ 543 if (buf_size!=i) { 544 per_load = buf_size/p_mask_size; 545 buf_size = per_load*p_mask_size; 546 gs->num_loads = num_loads = negl/per_load + (negl%per_load>0); 547 } 548 549 550 /* convert buf sizes from #bytes to #ints - 32 bit only! */ 551 p_mask_size/=sizeof(PetscInt); ngh_buf_size/=sizeof(PetscInt); buf_size/=sizeof(PetscInt); 552 553 /* find PCTFS_giop work space */ 554 buf2 = buf1+buf_size; 555 556 /* hold #ints needed for processor masks */ 557 gs->mask_sz=p_mask_size; 558 559 /* init buffers */ 560 ierr = PCTFS_ivec_zero(sh_proc_mask,p_mask_size);CHKERRQ(ierr); 561 ierr = PCTFS_ivec_zero(pw_sh_proc_mask,p_mask_size);CHKERRQ(ierr); 562 ierr = PCTFS_ivec_zero(ngh_buf,ngh_buf_size);CHKERRQ(ierr); 563 564 /* HACK reset tree info */ 565 tree_buf = NULL; 566 tree_buf_sz = ntree = 0; 567 568 /* ok do it */ 569 for (ptr1=ngh_buf,ptr2=elms,end=gs->gl_min,or_ct=i=0; or_ct<num_loads; or_ct++) { 570 /* identity for bitwise or is 000...000 */ 571 PCTFS_ivec_zero(buf1,buf_size); 572 573 /* load msg buffer */ 574 for (start=end,end+=per_load,i_start=i; (offset=*ptr2)<end; i++, ptr2++) { 575 offset = (offset-start)*p_mask_size; 576 PCTFS_ivec_copy(buf1+offset,p_mask,p_mask_size); 577 } 578 579 /* GLOBAL: pass buffer */ 580 ierr = PCTFS_giop(buf1,buf2,buf_size,&oper);CHKERRQ(ierr); 581 582 583 /* unload buffer into ngh_buf */ 584 ptr2=(elms+i_start); 585 for (ptr3=buf1,j=start; j<end; ptr3+=p_mask_size,j++) { 586 /* I own it ... may have to pairwise it */ 587 if (j==*ptr2) { 588 /* do i share it w/anyone? */ 589 ct1 = PCTFS_ct_bits((char*)ptr3,p_mask_size*sizeof(PetscInt)); 590 /* guess not */ 591 if (ct1<2) { ptr2++; ptr1+=p_mask_size; continue; } 592 593 /* i do ... so keep info and turn off my bit */ 594 PCTFS_ivec_copy(ptr1,ptr3,p_mask_size); 595 ierr = PCTFS_ivec_xor(ptr1,p_mask,p_mask_size);CHKERRQ(ierr); 596 ierr = PCTFS_ivec_or(sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr); 597 598 /* is it to be done pairwise? */ 599 if (--ct1<=level) { 600 npw++; 601 602 /* turn on high bit to indicate pw need to process */ 603 *ptr2++ |= TOP_BIT; 604 ierr = PCTFS_ivec_or(pw_sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr); 605 ptr1 += p_mask_size; 606 continue; 607 } 608 609 /* get set for next and note that I have a tree contribution */ 610 /* could save exact elm index for tree here -> save a search */ 611 ptr2++; ptr1+=p_mask_size; ntree_map++; 612 } else { /* i don't but still might be involved in tree */ 613 614 /* shared by how many? */ 615 ct1 = PCTFS_ct_bits((char*)ptr3,p_mask_size*sizeof(PetscInt)); 616 617 /* none! */ 618 if (ct1<2) continue; 619 620 /* is it going to be done pairwise? but not by me of course!*/ 621 if (--ct1<=level) continue; 622 } 623 /* LATER we're going to have to process it NOW */ 624 /* nope ... tree it */ 625 ierr = place_in_tree(j);CHKERRQ(ierr); 626 } 627 } 628 629 free((void*)t_mask); 630 free((void*)buf1); 631 632 gs->len_pw_list = npw; 633 gs->num_nghs = PCTFS_ct_bits((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt)); 634 635 /* expand from bit mask list to int list and save ngh list */ 636 gs->nghs = (PetscInt*) malloc(gs->num_nghs * sizeof(PetscInt)); 637 PCTFS_bm_to_proc((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt),gs->nghs); 638 639 gs->num_pw_nghs = PCTFS_ct_bits((char*)pw_sh_proc_mask,p_mask_size*sizeof(PetscInt)); 640 641 oper = GL_MAX; 642 ct1 = gs->num_nghs; 643 ierr = PCTFS_giop(&ct1,&ct2,1,&oper);CHKERRQ(ierr); 644 gs->max_nghs = ct1; 645 646 gs->tree_map_sz = ntree_map; 647 gs->max_left_over=ntree; 648 649 free((void*)p_mask); 650 free((void*)sh_proc_mask); 651 PetscFunctionReturn(0); 652 } 653 654 /******************************************************************************/ 655 static PetscErrorCode set_pairwise(PCTFS_gs_id *gs) 656 { 657 PetscInt i, j; 658 PetscInt p_mask_size; 659 PetscInt *p_mask, *sh_proc_mask, *tmp_proc_mask; 660 PetscInt *ngh_buf, *buf2; 661 PetscInt offset; 662 PetscInt *msg_list, *msg_size, **msg_nodes, nprs; 663 PetscInt *pairwise_elm_list, len_pair_list=0; 664 PetscInt *iptr, t1, i_start, nel, *elms; 665 PetscInt ct; 666 PetscErrorCode ierr; 667 668 PetscFunctionBegin; 669 /* to make life easier */ 670 nel = gs->nel; 671 elms = gs->elms; 672 ngh_buf = gs->ngh_buf; 673 sh_proc_mask = gs->pw_nghs; 674 675 /* need a few temp masks */ 676 p_mask_size = PCTFS_len_bit_mask(PCTFS_num_nodes); 677 p_mask = (PetscInt*) malloc(p_mask_size); 678 tmp_proc_mask = (PetscInt*) malloc(p_mask_size); 679 680 /* set mask to my PCTFS_my_id's bit mask */ 681 ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr); 682 683 p_mask_size /= sizeof(PetscInt); 684 685 len_pair_list = gs->len_pw_list; 686 gs->pw_elm_list = pairwise_elm_list=(PetscInt*)malloc((len_pair_list+1)*sizeof(PetscInt)); 687 688 /* how many processors (nghs) do we have to exchange with? */ 689 nprs = gs->num_pairs = PCTFS_ct_bits((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt)); 690 691 692 /* allocate space for PCTFS_gs_gop() info */ 693 gs->pair_list = msg_list = (PetscInt*) malloc(sizeof(PetscInt)*nprs); 694 gs->msg_sizes = msg_size = (PetscInt*) malloc(sizeof(PetscInt)*nprs); 695 gs->node_list = msg_nodes = (PetscInt**) malloc(sizeof(PetscInt*)*(nprs+1)); 696 697 /* init msg_size list */ 698 ierr = PCTFS_ivec_zero(msg_size,nprs);CHKERRQ(ierr); 699 700 /* expand from bit mask list to int list */ 701 ierr = PCTFS_bm_to_proc((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt),msg_list);CHKERRQ(ierr); 702 703 /* keep list of elements being handled pairwise */ 704 for (i=j=0; i<nel; i++) { 705 if (elms[i] & TOP_BIT) { elms[i] ^= TOP_BIT; pairwise_elm_list[j++] = i; } 706 } 707 pairwise_elm_list[j] = -1; 708 709 gs->msg_ids_out = (MPI_Request*) malloc(sizeof(MPI_Request)*(nprs+1)); 710 gs->msg_ids_out[nprs] = MPI_REQUEST_NULL; 711 gs->msg_ids_in = (MPI_Request*) malloc(sizeof(MPI_Request)*(nprs+1)); 712 gs->msg_ids_in[nprs] = MPI_REQUEST_NULL; 713 gs->pw_vals = (PetscScalar*) malloc(sizeof(PetscScalar)*len_pair_list*vec_sz); 714 715 /* find who goes to each processor */ 716 for (i_start=i=0; i<nprs; i++) { 717 /* processor i's mask */ 718 ierr = PCTFS_set_bit_mask(p_mask,p_mask_size*sizeof(PetscInt),msg_list[i]);CHKERRQ(ierr); 719 720 /* det # going to processor i */ 721 for (ct=j=0; j<len_pair_list; j++) { 722 buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size); 723 ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr); 724 if (PCTFS_ct_bits((char*)tmp_proc_mask,p_mask_size*sizeof(PetscInt))) ct++; 725 } 726 msg_size[i] = ct; 727 i_start = PetscMax(i_start,ct); 728 729 /*space to hold nodes in message to first neighbor */ 730 msg_nodes[i] = iptr = (PetscInt*) malloc(sizeof(PetscInt)*(ct+1)); 731 732 for (j=0;j<len_pair_list;j++) { 733 buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size); 734 ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr); 735 if (PCTFS_ct_bits((char*)tmp_proc_mask,p_mask_size*sizeof(PetscInt))) *iptr++ = j; 736 } 737 *iptr = -1; 738 } 739 msg_nodes[nprs] = NULL; 740 741 j = gs->loc_node_pairs=i_start; 742 t1 = GL_MAX; 743 ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 744 gs->max_node_pairs = i_start; 745 746 i_start = j; 747 t1 = GL_MIN; 748 ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 749 gs->min_node_pairs = i_start; 750 751 i_start = j; 752 t1 = GL_ADD; 753 ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 754 gs->avg_node_pairs = i_start/PCTFS_num_nodes + 1; 755 756 i_start = nprs; 757 t1 = GL_MAX; 758 PCTFS_giop(&i_start,&offset,1,&t1); 759 gs->max_pairs = i_start; 760 761 762 /* remap pairwise in tail of gsi_via_bit_mask() */ 763 gs->msg_total = PCTFS_ivec_sum(gs->msg_sizes,nprs); 764 gs->out = (PetscScalar*) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz); 765 gs->in = (PetscScalar*) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz); 766 767 /* reset malloc pool */ 768 free((void*)p_mask); 769 free((void*)tmp_proc_mask); 770 PetscFunctionReturn(0); 771 } 772 773 /* to do pruned tree just save ngh buf copy for each one and decode here! 774 ******************************************************************************/ 775 static PetscErrorCode set_tree(PCTFS_gs_id *gs) 776 { 777 PetscInt i, j, n, nel; 778 PetscInt *iptr_in, *iptr_out, *tree_elms, *elms; 779 780 PetscFunctionBegin; 781 /* local work ptrs */ 782 elms = gs->elms; 783 nel = gs->nel; 784 785 /* how many via tree */ 786 gs->tree_nel = n = ntree; 787 gs->tree_elms = tree_elms = iptr_in = tree_buf; 788 gs->tree_buf = (PetscScalar*) malloc(sizeof(PetscScalar)*n*vec_sz); 789 gs->tree_work = (PetscScalar*) malloc(sizeof(PetscScalar)*n*vec_sz); 790 j = gs->tree_map_sz; 791 gs->tree_map_in = iptr_in = (PetscInt*) malloc(sizeof(PetscInt)*(j+1)); 792 gs->tree_map_out = iptr_out = (PetscInt*) malloc(sizeof(PetscInt)*(j+1)); 793 794 /* search the longer of the two lists */ 795 /* note ... could save this info in get_ngh_buf and save searches */ 796 if (n<=nel) { 797 /* bijective fct w/remap - search elm list */ 798 for (i=0; i<n; i++) { 799 if ((j=PCTFS_ivec_binary_search(*tree_elms++,elms,nel))>=0) {*iptr_in++ = j; *iptr_out++ = i;} 800 } 801 } else { 802 for (i=0; i<nel; i++) { 803 if ((j=PCTFS_ivec_binary_search(*elms++,tree_elms,n))>=0) {*iptr_in++ = i; *iptr_out++ = j;} 804 } 805 } 806 807 /* sentinel */ 808 *iptr_in = *iptr_out = -1; 809 PetscFunctionReturn(0); 810 } 811 812 /******************************************************************************/ 813 static PetscErrorCode PCTFS_gs_gop_local_out(PCTFS_gs_id *gs, PetscScalar *vals) 814 { 815 PetscInt *num, *map, **reduce; 816 PetscScalar tmp; 817 818 PetscFunctionBegin; 819 num = gs->num_gop_local_reduce; 820 reduce = gs->gop_local_reduce; 821 while ((map = *reduce++)) { 822 /* wall */ 823 if (*num == 2) { 824 num++; 825 vals[map[1]] = vals[map[0]]; 826 } else if (*num == 3) { /* corner shared by three elements */ 827 num++; 828 vals[map[2]] = vals[map[1]] = vals[map[0]]; 829 } else if (*num == 4) { /* corner shared by four elements */ 830 num++; 831 vals[map[3]] = vals[map[2]] = vals[map[1]] = vals[map[0]]; 832 } else { /* general case ... odd geoms ... 3D*/ 833 num++; 834 tmp = *(vals + *map++); 835 while (*map >= 0) *(vals + *map++) = tmp; 836 } 837 } 838 PetscFunctionReturn(0); 839 } 840 841 /******************************************************************************/ 842 static PetscErrorCode PCTFS_gs_gop_local_plus(PCTFS_gs_id *gs, PetscScalar *vals) 843 { 844 PetscInt *num, *map, **reduce; 845 PetscScalar tmp; 846 847 PetscFunctionBegin; 848 num = gs->num_local_reduce; 849 reduce = gs->local_reduce; 850 while ((map = *reduce)) { 851 /* wall */ 852 if (*num == 2) { 853 num++; reduce++; 854 vals[map[1]] = vals[map[0]] += vals[map[1]]; 855 } else if (*num == 3) { /* corner shared by three elements */ 856 num++; reduce++; 857 vals[map[2]]=vals[map[1]]=vals[map[0]]+=(vals[map[1]]+vals[map[2]]); 858 } else if (*num == 4) { /* corner shared by four elements */ 859 num++; reduce++; 860 vals[map[1]]=vals[map[2]]=vals[map[3]]=vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]); 861 } else { /* general case ... odd geoms ... 3D*/ 862 num++; 863 tmp = 0.0; 864 while (*map >= 0) tmp += *(vals + *map++); 865 866 map = *reduce++; 867 while (*map >= 0) *(vals + *map++) = tmp; 868 } 869 } 870 PetscFunctionReturn(0); 871 } 872 873 /******************************************************************************/ 874 static PetscErrorCode PCTFS_gs_gop_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals) 875 { 876 PetscInt *num, *map, **reduce; 877 PetscScalar *base; 878 879 PetscFunctionBegin; 880 num = gs->num_gop_local_reduce; 881 reduce = gs->gop_local_reduce; 882 while ((map = *reduce++)) { 883 /* wall */ 884 if (*num == 2) { 885 num++; 886 vals[map[0]] += vals[map[1]]; 887 } else if (*num == 3) { /* corner shared by three elements */ 888 num++; 889 vals[map[0]] += (vals[map[1]] + vals[map[2]]); 890 } else if (*num == 4) { /* corner shared by four elements */ 891 num++; 892 vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]); 893 } else { /* general case ... odd geoms ... 3D*/ 894 num++; 895 base = vals + *map++; 896 while (*map >= 0) *base += *(vals + *map++); 897 } 898 } 899 PetscFunctionReturn(0); 900 } 901 902 /******************************************************************************/ 903 PetscErrorCode PCTFS_gs_free(PCTFS_gs_id *gs) 904 { 905 PetscInt i; 906 907 PetscFunctionBegin; 908 if (gs->nghs) free((void*) gs->nghs); 909 if (gs->pw_nghs) free((void*) gs->pw_nghs); 910 911 /* tree */ 912 if (gs->max_left_over) { 913 if (gs->tree_elms) free((void*) gs->tree_elms); 914 if (gs->tree_buf) free((void*) gs->tree_buf); 915 if (gs->tree_work) free((void*) gs->tree_work); 916 if (gs->tree_map_in) free((void*) gs->tree_map_in); 917 if (gs->tree_map_out) free((void*) gs->tree_map_out); 918 } 919 920 /* pairwise info */ 921 if (gs->num_pairs) { 922 /* should be NULL already */ 923 if (gs->ngh_buf) free((void*) gs->ngh_buf); 924 if (gs->elms) free((void*) gs->elms); 925 if (gs->local_elms) free((void*) gs->local_elms); 926 if (gs->companion) free((void*) gs->companion); 927 928 /* only set if pairwise */ 929 if (gs->vals) free((void*) gs->vals); 930 if (gs->in) free((void*) gs->in); 931 if (gs->out) free((void*) gs->out); 932 if (gs->msg_ids_in) free((void*) gs->msg_ids_in); 933 if (gs->msg_ids_out) free((void*) gs->msg_ids_out); 934 if (gs->pw_vals) free((void*) gs->pw_vals); 935 if (gs->pw_elm_list) free((void*) gs->pw_elm_list); 936 if (gs->node_list) { 937 for (i=0;i<gs->num_pairs;i++) { 938 if (gs->node_list[i]) { 939 free((void*) gs->node_list[i]); 940 } 941 } 942 free((void*) gs->node_list); 943 } 944 if (gs->msg_sizes) free((void*) gs->msg_sizes); 945 if (gs->pair_list) free((void*) gs->pair_list); 946 } 947 948 /* local info */ 949 if (gs->num_local_total>=0) { 950 for (i=0;i<gs->num_local_total+1;i++) { 951 if (gs->num_gop_local_reduce[i]) free((void*) gs->gop_local_reduce[i]); 952 } 953 } 954 955 /* if intersection tree/pairwise and local isn't empty */ 956 if (gs->gop_local_reduce) free((void*) gs->gop_local_reduce); 957 if (gs->num_gop_local_reduce) free((void*) gs->num_gop_local_reduce); 958 959 free((void*) gs); 960 PetscFunctionReturn(0); 961 } 962 963 /******************************************************************************/ 964 PetscErrorCode PCTFS_gs_gop_vec(PCTFS_gs_id *gs, PetscScalar *vals, const char *op, PetscInt step) 965 { 966 PetscErrorCode ierr; 967 968 PetscFunctionBegin; 969 switch (*op) { 970 case '+': 971 PCTFS_gs_gop_vec_plus(gs,vals,step); 972 break; 973 default: 974 ierr = PetscInfo1(0,"PCTFS_gs_gop_vec() :: %c is not a valid op\n",op[0]);CHKERRQ(ierr); 975 ierr = PetscInfo(0,"PCTFS_gs_gop_vec() :: default :: plus\n");CHKERRQ(ierr); 976 PCTFS_gs_gop_vec_plus(gs,vals,step); 977 break; 978 } 979 PetscFunctionReturn(0); 980 } 981 982 /******************************************************************************/ 983 static PetscErrorCode PCTFS_gs_gop_vec_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 984 { 985 PetscFunctionBegin; 986 if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_gs_gop_vec() passed NULL gs handle!!!"); 987 988 /* local only operations!!! */ 989 if (gs->num_local) PCTFS_gs_gop_vec_local_plus(gs,vals,step); 990 991 /* if intersection tree/pairwise and local isn't empty */ 992 if (gs->num_local_gop) { 993 PCTFS_gs_gop_vec_local_in_plus(gs,vals,step); 994 995 /* pairwise */ 996 if (gs->num_pairs) PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step); 997 998 /* tree */ 999 else if (gs->max_left_over) PCTFS_gs_gop_vec_tree_plus(gs,vals,step); 1000 1001 PCTFS_gs_gop_vec_local_out(gs,vals,step); 1002 } else { /* if intersection tree/pairwise and local is empty */ 1003 /* pairwise */ 1004 if (gs->num_pairs) PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step); 1005 1006 /* tree */ 1007 else if (gs->max_left_over) PCTFS_gs_gop_vec_tree_plus(gs,vals,step); 1008 } 1009 PetscFunctionReturn(0); 1010 } 1011 1012 /******************************************************************************/ 1013 static PetscErrorCode PCTFS_gs_gop_vec_local_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1014 { 1015 PetscInt *num, *map, **reduce; 1016 PetscScalar *base; 1017 1018 PetscFunctionBegin; 1019 num = gs->num_local_reduce; 1020 reduce = gs->local_reduce; 1021 while ((map = *reduce)) { 1022 base = vals + map[0] * step; 1023 1024 /* wall */ 1025 if (*num == 2) { 1026 num++; reduce++; 1027 PCTFS_rvec_add (base,vals+map[1]*step,step); 1028 PCTFS_rvec_copy(vals+map[1]*step,base,step); 1029 } else if (*num == 3) { /* corner shared by three elements */ 1030 num++; reduce++; 1031 PCTFS_rvec_add (base,vals+map[1]*step,step); 1032 PCTFS_rvec_add (base,vals+map[2]*step,step); 1033 PCTFS_rvec_copy(vals+map[2]*step,base,step); 1034 PCTFS_rvec_copy(vals+map[1]*step,base,step); 1035 } else if (*num == 4) { /* corner shared by four elements */ 1036 num++; reduce++; 1037 PCTFS_rvec_add (base,vals+map[1]*step,step); 1038 PCTFS_rvec_add (base,vals+map[2]*step,step); 1039 PCTFS_rvec_add (base,vals+map[3]*step,step); 1040 PCTFS_rvec_copy(vals+map[3]*step,base,step); 1041 PCTFS_rvec_copy(vals+map[2]*step,base,step); 1042 PCTFS_rvec_copy(vals+map[1]*step,base,step); 1043 } else { /* general case ... odd geoms ... 3D */ 1044 num++; 1045 while (*++map >= 0) PCTFS_rvec_add (base,vals+*map*step,step); 1046 1047 map = *reduce; 1048 while (*++map >= 0) PCTFS_rvec_copy(vals+*map*step,base,step); 1049 1050 reduce++; 1051 } 1052 } 1053 PetscFunctionReturn(0); 1054 } 1055 1056 /******************************************************************************/ 1057 static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1058 { 1059 PetscInt *num, *map, **reduce; 1060 PetscScalar *base; 1061 1062 PetscFunctionBegin; 1063 num = gs->num_gop_local_reduce; 1064 reduce = gs->gop_local_reduce; 1065 while ((map = *reduce++)) { 1066 base = vals + map[0] * step; 1067 1068 /* wall */ 1069 if (*num == 2) { 1070 num++; 1071 PCTFS_rvec_add(base,vals+map[1]*step,step); 1072 } else if (*num == 3) { /* corner shared by three elements */ 1073 num++; 1074 PCTFS_rvec_add(base,vals+map[1]*step,step); 1075 PCTFS_rvec_add(base,vals+map[2]*step,step); 1076 } else if (*num == 4) { /* corner shared by four elements */ 1077 num++; 1078 PCTFS_rvec_add(base,vals+map[1]*step,step); 1079 PCTFS_rvec_add(base,vals+map[2]*step,step); 1080 PCTFS_rvec_add(base,vals+map[3]*step,step); 1081 } else { /* general case ... odd geoms ... 3D*/ 1082 num++; 1083 while (*++map >= 0) PCTFS_rvec_add(base,vals+*map*step,step); 1084 } 1085 } 1086 PetscFunctionReturn(0); 1087 } 1088 1089 /******************************************************************************/ 1090 static PetscErrorCode PCTFS_gs_gop_vec_local_out(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1091 { 1092 PetscInt *num, *map, **reduce; 1093 PetscScalar *base; 1094 1095 PetscFunctionBegin; 1096 num = gs->num_gop_local_reduce; 1097 reduce = gs->gop_local_reduce; 1098 while ((map = *reduce++)) { 1099 base = vals + map[0] * step; 1100 1101 /* wall */ 1102 if (*num == 2) { 1103 num++; 1104 PCTFS_rvec_copy(vals+map[1]*step,base,step); 1105 } else if (*num == 3) { /* corner shared by three elements */ 1106 num++; 1107 PCTFS_rvec_copy(vals+map[1]*step,base,step); 1108 PCTFS_rvec_copy(vals+map[2]*step,base,step); 1109 } else if (*num == 4) { /* corner shared by four elements */ 1110 num++; 1111 PCTFS_rvec_copy(vals+map[1]*step,base,step); 1112 PCTFS_rvec_copy(vals+map[2]*step,base,step); 1113 PCTFS_rvec_copy(vals+map[3]*step,base,step); 1114 } else { /* general case ... odd geoms ... 3D*/ 1115 num++; 1116 while (*++map >= 0) PCTFS_rvec_copy(vals+*map*step,base,step); 1117 } 1118 } 1119 PetscFunctionReturn(0); 1120 } 1121 1122 /******************************************************************************/ 1123 static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step) 1124 { 1125 PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2; 1126 PetscInt *iptr, *msg_list, *msg_size, **msg_nodes; 1127 PetscInt *pw, *list, *size, **nodes; 1128 MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out; 1129 MPI_Status status; 1130 PetscBLASInt i1 = 1,dstep; 1131 PetscErrorCode ierr; 1132 1133 PetscFunctionBegin; 1134 /* strip and load s */ 1135 msg_list = list = gs->pair_list; 1136 msg_size = size = gs->msg_sizes; 1137 msg_nodes = nodes = gs->node_list; 1138 iptr = pw = gs->pw_elm_list; 1139 dptr1 = dptr3 = gs->pw_vals; 1140 msg_ids_in = ids_in = gs->msg_ids_in; 1141 msg_ids_out = ids_out = gs->msg_ids_out; 1142 dptr2 = gs->out; 1143 in1=in2 = gs->in; 1144 1145 /* post the receives */ 1146 /* msg_nodes=nodes; */ 1147 do { 1148 /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the 1149 second one *list and do list++ afterwards */ 1150 ierr = MPI_Irecv(in1, *size *step, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr); 1151 list++;msg_ids_in++; 1152 in1 += *size++ *step; 1153 } while (*++msg_nodes); 1154 msg_nodes=nodes; 1155 1156 /* load gs values into in out gs buffers */ 1157 while (*iptr >= 0) { 1158 PCTFS_rvec_copy(dptr3,in_vals + *iptr*step,step); 1159 dptr3+=step; 1160 iptr++; 1161 } 1162 1163 /* load out buffers and post the sends */ 1164 while ((iptr = *msg_nodes++)) { 1165 dptr3 = dptr2; 1166 while (*iptr >= 0) { 1167 PCTFS_rvec_copy(dptr2,dptr1 + *iptr*step,step); 1168 dptr2+=step; 1169 iptr++; 1170 } 1171 ierr = MPI_Isend(dptr3, *msg_size *step, MPIU_SCALAR, *msg_list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr); 1172 msg_size++; msg_list++;msg_ids_out++; 1173 } 1174 1175 /* tree */ 1176 if (gs->max_left_over) PCTFS_gs_gop_vec_tree_plus(gs,in_vals,step); 1177 1178 /* process the received data */ 1179 msg_nodes=nodes; 1180 while ((iptr = *nodes++)) { 1181 PetscScalar d1 = 1.0; 1182 1183 /* Should I check the return value of MPI_Wait() or status? */ 1184 /* Can this loop be replaced by a call to MPI_Waitall()? */ 1185 ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr); 1186 ids_in++; 1187 while (*iptr >= 0) { 1188 ierr = PetscBLASIntCast(step,&dstep);CHKERRQ(ierr); 1189 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dstep,&d1,in2,&i1,dptr1 + *iptr*step,&i1)); 1190 in2+=step; 1191 iptr++; 1192 } 1193 } 1194 1195 /* replace vals */ 1196 while (*pw >= 0) { 1197 PCTFS_rvec_copy(in_vals + *pw*step,dptr1,step); 1198 dptr1+=step; 1199 pw++; 1200 } 1201 1202 /* clear isend message handles */ 1203 /* This changed for clarity though it could be the same */ 1204 1205 /* Should I check the return value of MPI_Wait() or status? */ 1206 /* Can this loop be replaced by a call to MPI_Waitall()? */ 1207 while (*msg_nodes++) { 1208 ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr); 1209 ids_out++; 1210 } 1211 PetscFunctionReturn(0); 1212 } 1213 1214 /******************************************************************************/ 1215 static PetscErrorCode PCTFS_gs_gop_vec_tree_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1216 { 1217 PetscInt size, *in, *out; 1218 PetscScalar *buf, *work; 1219 PetscInt op[] = {GL_ADD,0}; 1220 PetscBLASInt i1 = 1; 1221 PetscErrorCode ierr; 1222 PetscBLASInt dstep; 1223 1224 PetscFunctionBegin; 1225 /* copy over to local variables */ 1226 in = gs->tree_map_in; 1227 out = gs->tree_map_out; 1228 buf = gs->tree_buf; 1229 work = gs->tree_work; 1230 size = gs->tree_nel*step; 1231 1232 /* zero out collection buffer */ 1233 PCTFS_rvec_zero(buf,size); 1234 1235 1236 /* copy over my contributions */ 1237 while (*in >= 0) { 1238 ierr = PetscBLASIntCast(step,&dstep);CHKERRQ(ierr); 1239 PetscStackCallBLAS("BLAScopy",BLAScopy_(&dstep,vals + *in++ * step,&i1,buf + *out++ * step,&i1)); 1240 } 1241 1242 /* perform fan in/out on full buffer */ 1243 /* must change PCTFS_grop to handle the blas */ 1244 PCTFS_grop(buf,work,size,op); 1245 1246 /* reset */ 1247 in = gs->tree_map_in; 1248 out = gs->tree_map_out; 1249 1250 /* get the portion of the results I need */ 1251 while (*in >= 0) { 1252 ierr = PetscBLASIntCast(step,&dstep);CHKERRQ(ierr); 1253 PetscStackCallBLAS("BLAScopy",BLAScopy_(&dstep,buf + *out++ * step,&i1,vals + *in++ * step,&i1)); 1254 } 1255 PetscFunctionReturn(0); 1256 } 1257 1258 /******************************************************************************/ 1259 PetscErrorCode PCTFS_gs_gop_hc(PCTFS_gs_id *gs, PetscScalar *vals, const char *op, PetscInt dim) 1260 { 1261 PetscErrorCode ierr; 1262 1263 PetscFunctionBegin; 1264 switch (*op) { 1265 case '+': 1266 PCTFS_gs_gop_plus_hc(gs,vals,dim); 1267 break; 1268 default: 1269 ierr = PetscInfo1(0,"PCTFS_gs_gop_hc() :: %c is not a valid op\n",op[0]);CHKERRQ(ierr); 1270 ierr = PetscInfo(0,"PCTFS_gs_gop_hc() :: default :: plus\n");CHKERRQ(ierr); 1271 PCTFS_gs_gop_plus_hc(gs,vals,dim); 1272 break; 1273 } 1274 PetscFunctionReturn(0); 1275 } 1276 1277 /******************************************************************************/ 1278 static PetscErrorCode PCTFS_gs_gop_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim) 1279 { 1280 PetscFunctionBegin; 1281 /* if there's nothing to do return */ 1282 if (dim<=0) PetscFunctionReturn(0); 1283 1284 /* can't do more dimensions then exist */ 1285 dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 1286 1287 /* local only operations!!! */ 1288 if (gs->num_local) PCTFS_gs_gop_local_plus(gs,vals); 1289 1290 /* if intersection tree/pairwise and local isn't empty */ 1291 if (gs->num_local_gop) { 1292 PCTFS_gs_gop_local_in_plus(gs,vals); 1293 1294 /* pairwise will do tree inside ... */ 1295 if (gs->num_pairs) PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim); /* tree only */ 1296 else if (gs->max_left_over) PCTFS_gs_gop_tree_plus_hc(gs,vals,dim); 1297 1298 PCTFS_gs_gop_local_out(gs,vals); 1299 } else { /* if intersection tree/pairwise and local is empty */ 1300 /* pairwise will do tree inside */ 1301 if (gs->num_pairs) PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim); /* tree */ 1302 else if (gs->max_left_over) PCTFS_gs_gop_tree_plus_hc(gs,vals,dim); 1303 } 1304 PetscFunctionReturn(0); 1305 } 1306 1307 /******************************************************************************/ 1308 static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim) 1309 { 1310 PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2; 1311 PetscInt *iptr, *msg_list, *msg_size, **msg_nodes; 1312 PetscInt *pw, *list, *size, **nodes; 1313 MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out; 1314 MPI_Status status; 1315 PetscInt i, mask=1; 1316 PetscErrorCode ierr; 1317 1318 PetscFunctionBegin; 1319 for (i=1; i<dim; i++) { mask<<=1; mask++; } 1320 1321 /* strip and load s */ 1322 msg_list = list = gs->pair_list; 1323 msg_size = size = gs->msg_sizes; 1324 msg_nodes = nodes = gs->node_list; 1325 iptr = pw = gs->pw_elm_list; 1326 dptr1 = dptr3 = gs->pw_vals; 1327 msg_ids_in = ids_in = gs->msg_ids_in; 1328 msg_ids_out = ids_out = gs->msg_ids_out; 1329 dptr2 = gs->out; 1330 in1 = in2 = gs->in; 1331 1332 /* post the receives */ 1333 /* msg_nodes=nodes; */ 1334 do { 1335 /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the 1336 second one *list and do list++ afterwards */ 1337 if ((PCTFS_my_id|mask)==(*list|mask)) { 1338 ierr = MPI_Irecv(in1, *size, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr); 1339 list++; msg_ids_in++;in1 += *size++; 1340 } else { list++; size++; } 1341 } while (*++msg_nodes); 1342 1343 /* load gs values into in out gs buffers */ 1344 while (*iptr >= 0) *dptr3++ = *(in_vals + *iptr++); 1345 1346 /* load out buffers and post the sends */ 1347 msg_nodes=nodes; 1348 list = msg_list; 1349 while ((iptr = *msg_nodes++)) { 1350 if ((PCTFS_my_id|mask)==(*list|mask)) { 1351 dptr3 = dptr2; 1352 while (*iptr >= 0) *dptr2++ = *(dptr1 + *iptr++); 1353 /* CHECK PERSISTENT COMMS MODE FOR ALL THIS STUFF */ 1354 /* is msg_ids_out++ correct? */ 1355 ierr = MPI_Isend(dptr3, *msg_size, MPIU_SCALAR, *list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr); 1356 msg_size++;list++;msg_ids_out++; 1357 } else {list++; msg_size++;} 1358 } 1359 1360 /* do the tree while we're waiting */ 1361 if (gs->max_left_over) PCTFS_gs_gop_tree_plus_hc(gs,in_vals,dim); 1362 1363 /* process the received data */ 1364 msg_nodes=nodes; 1365 list = msg_list; 1366 while ((iptr = *nodes++)) { 1367 if ((PCTFS_my_id|mask)==(*list|mask)) { 1368 /* Should I check the return value of MPI_Wait() or status? */ 1369 /* Can this loop be replaced by a call to MPI_Waitall()? */ 1370 ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr); 1371 ids_in++; 1372 while (*iptr >= 0) *(dptr1 + *iptr++) += *in2++; 1373 } 1374 list++; 1375 } 1376 1377 /* replace vals */ 1378 while (*pw >= 0) *(in_vals + *pw++) = *dptr1++; 1379 1380 /* clear isend message handles */ 1381 /* This changed for clarity though it could be the same */ 1382 while (*msg_nodes++) { 1383 if ((PCTFS_my_id|mask)==(*msg_list|mask)) { 1384 /* Should I check the return value of MPI_Wait() or status? */ 1385 /* Can this loop be replaced by a call to MPI_Waitall()? */ 1386 ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr); 1387 ids_out++; 1388 } 1389 msg_list++; 1390 } 1391 PetscFunctionReturn(0); 1392 } 1393 1394 /******************************************************************************/ 1395 static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim) 1396 { 1397 PetscInt size; 1398 PetscInt *in, *out; 1399 PetscScalar *buf, *work; 1400 PetscInt op[] = {GL_ADD,0}; 1401 1402 PetscFunctionBegin; 1403 in = gs->tree_map_in; 1404 out = gs->tree_map_out; 1405 buf = gs->tree_buf; 1406 work = gs->tree_work; 1407 size = gs->tree_nel; 1408 1409 PCTFS_rvec_zero(buf,size); 1410 1411 while (*in >= 0) *(buf + *out++) = *(vals + *in++); 1412 1413 in = gs->tree_map_in; 1414 out = gs->tree_map_out; 1415 1416 PCTFS_grop_hc(buf,work,size,op,dim); 1417 1418 while (*in >= 0) *(vals + *in++) = *(buf + *out++); 1419 PetscFunctionReturn(0); 1420 } 1421 1422 1423 1424