1827bd09bSSatish Balay 2827bd09bSSatish Balay /***********************************gs.c*************************************** 3827bd09bSSatish Balay 4827bd09bSSatish Balay Author: Henry M. Tufo III 5827bd09bSSatish Balay 6827bd09bSSatish Balay e-mail: hmt@cs.brown.edu 7827bd09bSSatish Balay 8827bd09bSSatish Balay snail-mail: 9827bd09bSSatish Balay Division of Applied Mathematics 10827bd09bSSatish Balay Brown University 11827bd09bSSatish Balay Providence, RI 02912 12827bd09bSSatish Balay 13827bd09bSSatish Balay Last Modification: 14827bd09bSSatish Balay 6.21.97 15827bd09bSSatish Balay ************************************gs.c**************************************/ 16827bd09bSSatish Balay 17827bd09bSSatish Balay /***********************************gs.c*************************************** 18827bd09bSSatish Balay File Description: 19827bd09bSSatish Balay ----------------- 20827bd09bSSatish Balay 21827bd09bSSatish Balay ************************************gs.c**************************************/ 22827bd09bSSatish Balay 23c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h> 2439945688SSatish Balay 25827bd09bSSatish Balay /* default length of number of items via tree - doubles if exceeded */ 26827bd09bSSatish Balay #define TREE_BUF_SZ 2048; 27827bd09bSSatish Balay #define GS_VEC_SZ 1 28827bd09bSSatish Balay 29827bd09bSSatish Balay 30827bd09bSSatish Balay 31827bd09bSSatish Balay /***********************************gs.c*************************************** 32827bd09bSSatish Balay Type: struct gather_scatter_id 33827bd09bSSatish Balay ------------------------------ 34827bd09bSSatish Balay 35827bd09bSSatish Balay ************************************gs.c**************************************/ 36827bd09bSSatish Balay typedef struct gather_scatter_id { 3752f87cdaSBarry Smith PetscInt id; 3852f87cdaSBarry Smith PetscInt nel_min; 3952f87cdaSBarry Smith PetscInt nel_max; 4052f87cdaSBarry Smith PetscInt nel_sum; 4152f87cdaSBarry Smith PetscInt negl; 4252f87cdaSBarry Smith PetscInt gl_max; 4352f87cdaSBarry Smith PetscInt gl_min; 4452f87cdaSBarry Smith PetscInt repeats; 4552f87cdaSBarry Smith PetscInt ordered; 4652f87cdaSBarry Smith PetscInt positive; 47a501084fSBarry Smith PetscScalar *vals; 48827bd09bSSatish Balay 49827bd09bSSatish Balay /* bit mask info */ 5052f87cdaSBarry Smith PetscInt *my_proc_mask; 5152f87cdaSBarry Smith PetscInt mask_sz; 5252f87cdaSBarry Smith PetscInt *ngh_buf; 5352f87cdaSBarry Smith PetscInt ngh_buf_sz; 5452f87cdaSBarry Smith PetscInt *nghs; 5552f87cdaSBarry Smith PetscInt num_nghs; 5652f87cdaSBarry Smith PetscInt max_nghs; 5752f87cdaSBarry Smith PetscInt *pw_nghs; 5852f87cdaSBarry Smith PetscInt num_pw_nghs; 5952f87cdaSBarry Smith PetscInt *tree_nghs; 6052f87cdaSBarry Smith PetscInt num_tree_nghs; 61827bd09bSSatish Balay 6252f87cdaSBarry Smith PetscInt num_loads; 63827bd09bSSatish Balay 64827bd09bSSatish Balay /* repeats == true -> local info */ 6552f87cdaSBarry Smith PetscInt nel; /* number of unique elememts */ 6652f87cdaSBarry Smith PetscInt *elms; /* of size nel */ 6752f87cdaSBarry Smith PetscInt nel_total; 6852f87cdaSBarry Smith PetscInt *local_elms; /* of size nel_total */ 6952f87cdaSBarry Smith PetscInt *companion; /* of size nel_total */ 70827bd09bSSatish Balay 71827bd09bSSatish Balay /* local info */ 7252f87cdaSBarry Smith PetscInt num_local_total; 7352f87cdaSBarry Smith PetscInt local_strength; 7452f87cdaSBarry Smith PetscInt num_local; 7552f87cdaSBarry Smith PetscInt *num_local_reduce; 7652f87cdaSBarry Smith PetscInt **local_reduce; 7752f87cdaSBarry Smith PetscInt num_local_gop; 7852f87cdaSBarry Smith PetscInt *num_gop_local_reduce; 7952f87cdaSBarry Smith PetscInt **gop_local_reduce; 80827bd09bSSatish Balay 81827bd09bSSatish Balay /* pairwise info */ 8252f87cdaSBarry Smith PetscInt level; 8352f87cdaSBarry Smith PetscInt num_pairs; 8452f87cdaSBarry Smith PetscInt max_pairs; 8552f87cdaSBarry Smith PetscInt loc_node_pairs; 8652f87cdaSBarry Smith PetscInt max_node_pairs; 8752f87cdaSBarry Smith PetscInt min_node_pairs; 8852f87cdaSBarry Smith PetscInt avg_node_pairs; 8952f87cdaSBarry Smith PetscInt *pair_list; 9052f87cdaSBarry Smith PetscInt *msg_sizes; 9152f87cdaSBarry Smith PetscInt **node_list; 9252f87cdaSBarry Smith PetscInt len_pw_list; 9352f87cdaSBarry Smith PetscInt *pw_elm_list; 94a501084fSBarry Smith PetscScalar *pw_vals; 95827bd09bSSatish Balay 96827bd09bSSatish Balay MPI_Request *msg_ids_in; 97827bd09bSSatish Balay MPI_Request *msg_ids_out; 98827bd09bSSatish Balay 99a501084fSBarry Smith PetscScalar *out; 100a501084fSBarry Smith PetscScalar *in; 10152f87cdaSBarry Smith PetscInt msg_total; 102827bd09bSSatish Balay 103827bd09bSSatish Balay /* tree - crystal accumulator info */ 10452f87cdaSBarry Smith PetscInt max_left_over; 10552f87cdaSBarry Smith PetscInt *pre; 10652f87cdaSBarry Smith PetscInt *in_num; 10752f87cdaSBarry Smith PetscInt *out_num; 10852f87cdaSBarry Smith PetscInt **in_list; 10952f87cdaSBarry Smith PetscInt **out_list; 110827bd09bSSatish Balay 111827bd09bSSatish Balay /* new tree work*/ 11252f87cdaSBarry Smith PetscInt tree_nel; 11352f87cdaSBarry Smith PetscInt *tree_elms; 114a501084fSBarry Smith PetscScalar *tree_buf; 115a501084fSBarry Smith PetscScalar *tree_work; 116827bd09bSSatish Balay 11752f87cdaSBarry Smith PetscInt tree_map_sz; 11852f87cdaSBarry Smith PetscInt *tree_map_in; 11952f87cdaSBarry Smith PetscInt *tree_map_out; 120827bd09bSSatish Balay 121827bd09bSSatish Balay /* current memory status */ 12252f87cdaSBarry Smith PetscInt gl_bss_min; 12352f87cdaSBarry Smith PetscInt gl_perm_min; 124827bd09bSSatish Balay 125ca8e9878SJed Brown /* max segment size for PCTFS_gs_gop_vec() */ 12652f87cdaSBarry Smith PetscInt vec_sz; 127827bd09bSSatish Balay 128827bd09bSSatish Balay /* hack to make paul happy */ 129ca8e9878SJed Brown MPI_Comm PCTFS_gs_comm; 130827bd09bSSatish Balay 131ca8e9878SJed Brown } PCTFS_gs_id; 132827bd09bSSatish Balay 133ca8e9878SJed Brown static PCTFS_gs_id *gsi_check_args(PetscInt *elms, PetscInt nel, PetscInt level); 134ca8e9878SJed Brown static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs); 135ca8e9878SJed Brown static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs); 136ca8e9878SJed Brown static PetscErrorCode set_pairwise(PCTFS_gs_id *gs); 137ca8e9878SJed Brown static PCTFS_gs_id *gsi_new(void); 138ca8e9878SJed Brown static PetscErrorCode set_tree(PCTFS_gs_id *gs); 139827bd09bSSatish Balay 140827bd09bSSatish Balay /* same for all but vector flavor */ 141ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_out(PCTFS_gs_id *gs, PetscScalar *vals); 142827bd09bSSatish Balay /* vector flavor */ 143ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_out(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 144827bd09bSSatish Balay 145ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step); 146ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step); 147ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 148ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 149ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_tree_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step); 150827bd09bSSatish Balay 151827bd09bSSatish Balay 152ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_plus(PCTFS_gs_id *gs, PetscScalar *vals); 153ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals); 154827bd09bSSatish Balay 155ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim); 156ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim); 157ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim); 158827bd09bSSatish Balay 159827bd09bSSatish Balay /* global vars */ 160827bd09bSSatish Balay /* from comm.c module */ 161827bd09bSSatish Balay 16252f87cdaSBarry Smith static PetscInt num_gs_ids = 0; 163827bd09bSSatish Balay 164827bd09bSSatish Balay /* should make this dynamic ... later */ 16552f87cdaSBarry Smith static PetscInt msg_buf =MAX_MSG_BUF; 16652f87cdaSBarry Smith static PetscInt vec_sz =GS_VEC_SZ; 16752f87cdaSBarry Smith static PetscInt *tree_buf =NULL; 16852f87cdaSBarry Smith static PetscInt tree_buf_sz=0; 16952f87cdaSBarry Smith static PetscInt ntree =0; 170827bd09bSSatish Balay 171f1ed62a8SBarry Smith /***************************************************************************/ 172ca8e9878SJed Brown PetscErrorCode PCTFS_gs_init_vec_sz(PetscInt size) 173827bd09bSSatish Balay { 1743fdc5746SBarry Smith PetscFunctionBegin; 175827bd09bSSatish Balay vec_sz = size; 1763fdc5746SBarry Smith PetscFunctionReturn(0); 177827bd09bSSatish Balay } 178827bd09bSSatish Balay 179f1ed62a8SBarry Smith /******************************************************************************/ 180ca8e9878SJed Brown PetscErrorCode PCTFS_gs_init_msg_buf_sz(PetscInt buf_size) 181827bd09bSSatish Balay { 1823fdc5746SBarry Smith PetscFunctionBegin; 183827bd09bSSatish Balay msg_buf = buf_size; 1843fdc5746SBarry Smith PetscFunctionReturn(0); 185827bd09bSSatish Balay } 186827bd09bSSatish Balay 187f1ed62a8SBarry Smith /******************************************************************************/ 188ca8e9878SJed Brown PCTFS_gs_id *PCTFS_gs_init(PetscInt *elms, PetscInt nel, PetscInt level) 189827bd09bSSatish Balay { 190ca8e9878SJed Brown PCTFS_gs_id *gs; 191ca8e9878SJed Brown MPI_Group PCTFS_gs_group; 192ca8e9878SJed Brown MPI_Comm PCTFS_gs_comm; 193f1ed62a8SBarry Smith PetscErrorCode ierr; 194827bd09bSSatish Balay 1953fdc5746SBarry Smith PetscFunctionBegin; 196827bd09bSSatish Balay /* ensure that communication package has been initialized */ 197b1c944f5SJed Brown PCTFS_comm_init(); 198827bd09bSSatish Balay 199827bd09bSSatish Balay 200827bd09bSSatish Balay /* determines if we have enough dynamic/semi-static memory */ 201827bd09bSSatish Balay /* checks input, allocs and sets gd_id template */ 202827bd09bSSatish Balay gs = gsi_check_args(elms,nel,level); 203827bd09bSSatish Balay 204827bd09bSSatish Balay /* only bit mask version up and working for the moment */ 205827bd09bSSatish Balay /* LATER :: get int list version working for sparse pblms */ 206f1ed62a8SBarry Smith ierr = gsi_via_bit_mask(gs);CHKERRABORT(PETSC_COMM_WORLD,ierr); 207827bd09bSSatish Balay 208827bd09bSSatish Balay 209ca8e9878SJed Brown ierr = MPI_Comm_group(MPI_COMM_WORLD,&PCTFS_gs_group);CHKERRABORT(PETSC_COMM_WORLD,ierr); 210ca8e9878SJed Brown ierr = MPI_Comm_create(MPI_COMM_WORLD,PCTFS_gs_group,&PCTFS_gs_comm);CHKERRABORT(PETSC_COMM_WORLD,ierr); 211*2fa5cd67SKarl Rupp 212ca8e9878SJed Brown gs->PCTFS_gs_comm=PCTFS_gs_comm; 213827bd09bSSatish Balay 214827bd09bSSatish Balay return(gs); 215827bd09bSSatish Balay } 216827bd09bSSatish Balay 217f1ed62a8SBarry Smith /******************************************************************************/ 218ca8e9878SJed Brown static PCTFS_gs_id *gsi_new(void) 219827bd09bSSatish Balay { 220f1ed62a8SBarry Smith PetscErrorCode ierr; 221ca8e9878SJed Brown PCTFS_gs_id *gs; 222ca8e9878SJed Brown gs = (PCTFS_gs_id*) malloc(sizeof(PCTFS_gs_id)); 223ca8e9878SJed Brown ierr = PetscMemzero(gs,sizeof(PCTFS_gs_id));CHKERRABORT(PETSC_COMM_WORLD,ierr); 224827bd09bSSatish Balay return(gs); 225827bd09bSSatish Balay } 226827bd09bSSatish Balay 227f1ed62a8SBarry Smith /******************************************************************************/ 228ca8e9878SJed Brown static PCTFS_gs_id *gsi_check_args(PetscInt *in_elms, PetscInt nel, PetscInt level) 229827bd09bSSatish Balay { 23052f87cdaSBarry Smith PetscInt i, j, k, t2; 23152f87cdaSBarry Smith PetscInt *companion, *elms, *unique, *iptr; 23252f87cdaSBarry Smith PetscInt num_local=0, *num_to_reduce, **local_reduce; 23352f87cdaSBarry Smith PetscInt oprs[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_MIN,GL_B_AND}; 23452f87cdaSBarry Smith PetscInt vals[sizeof(oprs)/sizeof(oprs[0])-1]; 23552f87cdaSBarry Smith PetscInt work[sizeof(oprs)/sizeof(oprs[0])-1]; 236ca8e9878SJed Brown PCTFS_gs_id *gs; 237d1528f56SBarry Smith PetscErrorCode ierr; 238827bd09bSSatish Balay 239827bd09bSSatish Balay 240c1235816SBarry Smith if (!in_elms) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"elms point to nothing!!!\n"); 241c1235816SBarry Smith if (nel<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"can't have fewer than 0 elms!!!\n"); 242827bd09bSSatish Balay 243db4deed7SKarl Rupp if (nel==0) { ierr = PetscInfo(0,"I don't have any elements!!!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); } 244827bd09bSSatish Balay 245827bd09bSSatish Balay /* get space for gs template */ 246827bd09bSSatish Balay gs = gsi_new(); 247827bd09bSSatish Balay gs->id = ++num_gs_ids; 248827bd09bSSatish Balay 249827bd09bSSatish Balay /* hmt 6.4.99 */ 250827bd09bSSatish Balay /* caller can set global ids that don't participate to 0 */ 251ca8e9878SJed Brown /* PCTFS_gs_init ignores all zeros in elm list */ 252827bd09bSSatish Balay /* negative global ids are still invalid */ 253*2fa5cd67SKarl Rupp for (i=j=0; i<nel; i++) { 254*2fa5cd67SKarl Rupp if (in_elms[i]!=0) j++; 255*2fa5cd67SKarl Rupp } 256827bd09bSSatish Balay 257827bd09bSSatish Balay k=nel; nel=j; 258827bd09bSSatish Balay 259827bd09bSSatish Balay /* copy over in_elms list and create inverse map */ 26052f87cdaSBarry Smith elms = (PetscInt*) malloc((nel+1)*sizeof(PetscInt)); 26152f87cdaSBarry Smith companion = (PetscInt*) malloc(nel*sizeof(PetscInt)); 2621d7d0905SBarry Smith 263db4deed7SKarl Rupp for (i=j=0; i<k; i++) { 264db4deed7SKarl Rupp if (in_elms[i]!=0) { elms[j] = in_elms[i]; companion[j++] = i; } 265827bd09bSSatish Balay } 266827bd09bSSatish Balay 267c1235816SBarry Smith if (j!=nel) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"nel j mismatch!\n"); 268827bd09bSSatish Balay 269827bd09bSSatish Balay /* pre-pass ... check to see if sorted */ 270827bd09bSSatish Balay elms[nel] = INT_MAX; 271827bd09bSSatish Balay iptr = elms; 272827bd09bSSatish Balay unique = elms+1; 273827bd09bSSatish Balay j =0; 274db4deed7SKarl Rupp while (*iptr!=INT_MAX) { 275db4deed7SKarl Rupp if (*iptr++>*unique++) { j=1; break; } 276827bd09bSSatish Balay } 277827bd09bSSatish Balay 278827bd09bSSatish Balay /* set up inverse map */ 279db4deed7SKarl Rupp if (j) { 280f1ed62a8SBarry Smith ierr = PetscInfo(0,"gsi_check_args() :: elm list *not* sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); 281ca8e9878SJed Brown ierr = PCTFS_SMI_sort((void*)elms, (void*)companion, nel, SORT_INTEGER);CHKERRABORT(PETSC_COMM_WORLD,ierr); 282*2fa5cd67SKarl Rupp } else { ierr = PetscInfo(0,"gsi_check_args() :: elm list sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); } 283827bd09bSSatish Balay elms[nel] = INT_MIN; 284827bd09bSSatish Balay 285827bd09bSSatish Balay /* first pass */ 286827bd09bSSatish Balay /* determine number of unique elements, check pd */ 287db4deed7SKarl Rupp for (i=k=0; i<nel; i+=j) { 288827bd09bSSatish Balay t2 = elms[i]; 289827bd09bSSatish Balay j = ++i; 290827bd09bSSatish Balay 291827bd09bSSatish Balay /* clump 'em for now */ 292*2fa5cd67SKarl Rupp while (elms[j]==t2) j++; 293827bd09bSSatish Balay 294827bd09bSSatish Balay /* how many together and num local */ 295db4deed7SKarl Rupp if (j-=i) { num_local++; k+=j; } 296827bd09bSSatish Balay } 297827bd09bSSatish Balay 298827bd09bSSatish Balay /* how many unique elements? */ 299827bd09bSSatish Balay gs->repeats = k; 300827bd09bSSatish Balay gs->nel = nel-k; 301827bd09bSSatish Balay 302827bd09bSSatish Balay 303827bd09bSSatish Balay /* number of repeats? */ 304827bd09bSSatish Balay gs->num_local = num_local; 305827bd09bSSatish Balay num_local += 2; 30652f87cdaSBarry Smith gs->local_reduce = local_reduce=(PetscInt**)malloc(num_local*sizeof(PetscInt*)); 30752f87cdaSBarry Smith gs->num_local_reduce = num_to_reduce=(PetscInt*) malloc(num_local*sizeof(PetscInt)); 308827bd09bSSatish Balay 30952f87cdaSBarry Smith unique = (PetscInt*) malloc((gs->nel+1)*sizeof(PetscInt)); 310827bd09bSSatish Balay gs->elms = unique; 311827bd09bSSatish Balay gs->nel_total = nel; 312827bd09bSSatish Balay gs->local_elms = elms; 313827bd09bSSatish Balay gs->companion = companion; 314827bd09bSSatish Balay 315827bd09bSSatish Balay /* compess map as well as keep track of local ops */ 316db4deed7SKarl Rupp for (num_local=i=j=0; i<gs->nel; i++) { 317827bd09bSSatish Balay k = j; 318827bd09bSSatish Balay t2 = unique[i] = elms[j]; 319827bd09bSSatish Balay companion[i] = companion[j]; 320827bd09bSSatish Balay 321*2fa5cd67SKarl Rupp while (elms[j]==t2) j++; 322827bd09bSSatish Balay 323db4deed7SKarl Rupp if ((t2=(j-k))>1) { 324827bd09bSSatish Balay /* number together */ 325827bd09bSSatish Balay num_to_reduce[num_local] = t2++; 326*2fa5cd67SKarl Rupp 32752f87cdaSBarry Smith iptr = local_reduce[num_local++] = (PetscInt*)malloc(t2*sizeof(PetscInt)); 328827bd09bSSatish Balay 329827bd09bSSatish Balay /* to use binary searching don't remap until we check intersection */ 330827bd09bSSatish Balay *iptr++ = i; 331827bd09bSSatish Balay 332827bd09bSSatish Balay /* note that we're skipping the first one */ 333*2fa5cd67SKarl Rupp while (++k<j) *(iptr++) = companion[k]; 334827bd09bSSatish Balay *iptr = -1; 335827bd09bSSatish Balay } 336827bd09bSSatish Balay } 337827bd09bSSatish Balay 338827bd09bSSatish Balay /* sentinel for ngh_buf */ 339827bd09bSSatish Balay unique[gs->nel]=INT_MAX; 340827bd09bSSatish Balay 341827bd09bSSatish Balay /* for two partition sort hack */ 342827bd09bSSatish Balay num_to_reduce[num_local] = 0; 343827bd09bSSatish Balay local_reduce[num_local] = NULL; 344827bd09bSSatish Balay num_to_reduce[++num_local] = 0; 345827bd09bSSatish Balay local_reduce[num_local] = NULL; 346827bd09bSSatish Balay 347827bd09bSSatish Balay /* load 'em up */ 348827bd09bSSatish Balay /* note one extra to hold NON_UNIFORM flag!!! */ 349827bd09bSSatish Balay vals[2] = vals[1] = vals[0] = nel; 350db4deed7SKarl Rupp if (gs->nel>0) { 3511d7d0905SBarry Smith vals[3] = unique[0]; 3521d7d0905SBarry Smith vals[4] = unique[gs->nel-1]; 353db4deed7SKarl Rupp } else { 3541d7d0905SBarry Smith vals[3] = INT_MAX; 3551d7d0905SBarry Smith vals[4] = INT_MIN; 356827bd09bSSatish Balay } 357827bd09bSSatish Balay vals[5] = level; 358827bd09bSSatish Balay vals[6] = num_gs_ids; 359827bd09bSSatish Balay 360827bd09bSSatish Balay /* GLOBAL: send 'em out */ 361b1c944f5SJed Brown ierr = PCTFS_giop(vals,work,sizeof(oprs)/sizeof(oprs[0])-1,oprs);CHKERRABORT(PETSC_COMM_WORLD,ierr); 362827bd09bSSatish Balay 363827bd09bSSatish Balay /* must be semi-pos def - only pairwise depends on this */ 364827bd09bSSatish Balay /* LATER - remove this restriction */ 365c1235816SBarry Smith if (vals[3]<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system not semi-pos def \n"); 366c1235816SBarry Smith if (vals[4]==INT_MAX) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system ub too large !\n"); 367827bd09bSSatish Balay 368827bd09bSSatish Balay gs->nel_min = vals[0]; 369827bd09bSSatish Balay gs->nel_max = vals[1]; 370827bd09bSSatish Balay gs->nel_sum = vals[2]; 371827bd09bSSatish Balay gs->gl_min = vals[3]; 372827bd09bSSatish Balay gs->gl_max = vals[4]; 373827bd09bSSatish Balay gs->negl = vals[4]-vals[3]+1; 374827bd09bSSatish Balay 375c1235816SBarry Smith if (gs->negl<=0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system empty or neg :: %d\n"); 376827bd09bSSatish Balay 377827bd09bSSatish Balay /* LATER :: add level == -1 -> program selects level */ 378*2fa5cd67SKarl Rupp if (vals[5]<0) vals[5]=0; 379*2fa5cd67SKarl Rupp else if (vals[5]>PCTFS_num_nodes) vals[5]=PCTFS_num_nodes; 380827bd09bSSatish Balay gs->level = vals[5]; 381827bd09bSSatish Balay 382827bd09bSSatish Balay return(gs); 383827bd09bSSatish Balay } 384827bd09bSSatish Balay 385f1ed62a8SBarry Smith /******************************************************************************/ 386ca8e9878SJed Brown static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs) 387827bd09bSSatish Balay { 38852f87cdaSBarry Smith PetscInt i, nel, *elms; 38952f87cdaSBarry Smith PetscInt t1; 39052f87cdaSBarry Smith PetscInt **reduce; 39152f87cdaSBarry Smith PetscInt *map; 392f1ed62a8SBarry Smith PetscErrorCode ierr; 393827bd09bSSatish Balay 394f1ed62a8SBarry Smith PetscFunctionBegin; 395ca8e9878SJed Brown /* totally local removes ... PCTFS_ct_bits == 0 */ 396827bd09bSSatish Balay get_ngh_buf(gs); 397827bd09bSSatish Balay 39894dd86cdSBarry Smith if (gs->level) set_pairwise(gs); 39994dd86cdSBarry Smith if (gs->max_left_over) set_tree(gs); 400827bd09bSSatish Balay 401827bd09bSSatish Balay /* intersection local and pairwise/tree? */ 402827bd09bSSatish Balay gs->num_local_total = gs->num_local; 403827bd09bSSatish Balay gs->gop_local_reduce = gs->local_reduce; 404827bd09bSSatish Balay gs->num_gop_local_reduce = gs->num_local_reduce; 405827bd09bSSatish Balay 406827bd09bSSatish Balay map = gs->companion; 407827bd09bSSatish Balay 408827bd09bSSatish Balay /* is there any local compression */ 409d890fc11SSatish Balay if (!gs->num_local) { 410827bd09bSSatish Balay gs->local_strength = NONE; 411827bd09bSSatish Balay gs->num_local_gop = 0; 412d890fc11SSatish Balay } else { 413827bd09bSSatish Balay /* ok find intersection */ 414827bd09bSSatish Balay map = gs->companion; 415827bd09bSSatish Balay reduce = gs->local_reduce; 416827bd09bSSatish Balay for (i=0, t1=0; i<gs->num_local; i++, reduce++) 417827bd09bSSatish Balay { 418ca8e9878SJed Brown if ((PCTFS_ivec_binary_search(**reduce,gs->pw_elm_list,gs->len_pw_list)>=0) 419827bd09bSSatish Balay || 420db4deed7SKarl Rupp PCTFS_ivec_binary_search(**reduce,gs->tree_map_in,gs->tree_map_sz)>=0) { 421827bd09bSSatish Balay t1++; 422e32f2f54SBarry Smith if (gs->num_local_reduce[i]<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nobody in list?"); 423827bd09bSSatish Balay gs->num_local_reduce[i] *= -1; 424827bd09bSSatish Balay } 425827bd09bSSatish Balay **reduce=map[**reduce]; 426827bd09bSSatish Balay } 427827bd09bSSatish Balay 428827bd09bSSatish Balay /* intersection is empty */ 429db4deed7SKarl Rupp if (!t1) { 430827bd09bSSatish Balay gs->local_strength = FULL; 431827bd09bSSatish Balay gs->num_local_gop = 0; 432db4deed7SKarl Rupp } else { /* intersection not empty */ 433827bd09bSSatish Balay gs->local_strength = PARTIAL; 434*2fa5cd67SKarl Rupp 435ca8e9878SJed Brown ierr = PCTFS_SMI_sort((void*)gs->num_local_reduce, (void*)gs->local_reduce, gs->num_local + 1, SORT_INT_PTR);CHKERRQ(ierr); 436827bd09bSSatish Balay 437827bd09bSSatish Balay gs->num_local_gop = t1; 438827bd09bSSatish Balay gs->num_local_total = gs->num_local; 439827bd09bSSatish Balay gs->num_local -= t1; 440827bd09bSSatish Balay gs->gop_local_reduce = gs->local_reduce; 441827bd09bSSatish Balay gs->num_gop_local_reduce = gs->num_local_reduce; 442827bd09bSSatish Balay 443*2fa5cd67SKarl Rupp for (i=0; i<t1; i++) { 444e32f2f54SBarry Smith if (gs->num_gop_local_reduce[i]>=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"they aren't negative?"); 445827bd09bSSatish Balay gs->num_gop_local_reduce[i] *= -1; 446827bd09bSSatish Balay gs->local_reduce++; 447827bd09bSSatish Balay gs->num_local_reduce++; 448827bd09bSSatish Balay } 449827bd09bSSatish Balay gs->local_reduce++; 450827bd09bSSatish Balay gs->num_local_reduce++; 451827bd09bSSatish Balay } 452827bd09bSSatish Balay } 453827bd09bSSatish Balay 454827bd09bSSatish Balay elms = gs->pw_elm_list; 455827bd09bSSatish Balay nel = gs->len_pw_list; 456*2fa5cd67SKarl Rupp for (i=0; i<nel; i++) elms[i] = map[elms[i]]; 457827bd09bSSatish Balay 458827bd09bSSatish Balay elms = gs->tree_map_in; 459827bd09bSSatish Balay nel = gs->tree_map_sz; 460*2fa5cd67SKarl Rupp for (i=0; i<nel; i++) elms[i] = map[elms[i]]; 461827bd09bSSatish Balay 462827bd09bSSatish Balay /* clean up */ 463a501084fSBarry Smith free((void*) gs->local_elms); 464a501084fSBarry Smith free((void*) gs->companion); 465a501084fSBarry Smith free((void*) gs->elms); 466a501084fSBarry Smith free((void*) gs->ngh_buf); 467827bd09bSSatish Balay gs->local_elms = gs->companion = gs->elms = gs->ngh_buf = NULL; 4683fdc5746SBarry Smith PetscFunctionReturn(0); 469827bd09bSSatish Balay } 470827bd09bSSatish Balay 471f1ed62a8SBarry Smith /******************************************************************************/ 47252f87cdaSBarry Smith static PetscErrorCode place_in_tree(PetscInt elm) 473827bd09bSSatish Balay { 47452f87cdaSBarry Smith PetscInt *tp, n; 475827bd09bSSatish Balay 4763fdc5746SBarry Smith PetscFunctionBegin; 477*2fa5cd67SKarl Rupp if (ntree==tree_buf_sz) { 478db4deed7SKarl Rupp if (tree_buf_sz) { 479827bd09bSSatish Balay tp = tree_buf; 480827bd09bSSatish Balay n = tree_buf_sz; 481827bd09bSSatish Balay tree_buf_sz<<=1; 48252f87cdaSBarry Smith tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt)); 483ca8e9878SJed Brown PCTFS_ivec_copy(tree_buf,tp,n); 484a501084fSBarry Smith free(tp); 485db4deed7SKarl Rupp } else { 486827bd09bSSatish Balay tree_buf_sz = TREE_BUF_SZ; 48752f87cdaSBarry Smith tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt)); 488827bd09bSSatish Balay } 489827bd09bSSatish Balay } 490827bd09bSSatish Balay 491827bd09bSSatish Balay tree_buf[ntree++] = elm; 4923fdc5746SBarry Smith PetscFunctionReturn(0); 493827bd09bSSatish Balay } 494827bd09bSSatish Balay 495f1ed62a8SBarry Smith /******************************************************************************/ 496ca8e9878SJed Brown static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs) 497827bd09bSSatish Balay { 49852f87cdaSBarry Smith PetscInt i, j, npw=0, ntree_map=0; 49952f87cdaSBarry Smith PetscInt p_mask_size, ngh_buf_size, buf_size; 50052f87cdaSBarry Smith PetscInt *p_mask, *sh_proc_mask, *pw_sh_proc_mask; 50152f87cdaSBarry Smith PetscInt *ngh_buf, *buf1, *buf2; 50252f87cdaSBarry Smith PetscInt offset, per_load, num_loads, or_ct, start, end; 50352f87cdaSBarry Smith PetscInt *ptr1, *ptr2, i_start, negl, nel, *elms; 50452f87cdaSBarry Smith PetscInt oper=GL_B_OR; 50552f87cdaSBarry Smith PetscInt *ptr3, *t_mask, level, ct1, ct2; 506f1ed62a8SBarry Smith PetscErrorCode ierr; 507827bd09bSSatish Balay 5083fdc5746SBarry Smith PetscFunctionBegin; 509827bd09bSSatish Balay /* to make life easier */ 510827bd09bSSatish Balay nel = gs->nel; 511827bd09bSSatish Balay elms = gs->elms; 512827bd09bSSatish Balay level = gs->level; 513827bd09bSSatish Balay 514b1c944f5SJed Brown /* det #bytes needed for processor bit masks and init w/mask cor. to PCTFS_my_id */ 515ca8e9878SJed Brown p_mask = (PetscInt*) malloc(p_mask_size=PCTFS_len_bit_mask(PCTFS_num_nodes)); 516ca8e9878SJed Brown ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr); 517827bd09bSSatish Balay 518827bd09bSSatish Balay /* allocate space for masks and info bufs */ 51952f87cdaSBarry Smith gs->nghs = sh_proc_mask = (PetscInt*) malloc(p_mask_size); 52052f87cdaSBarry Smith gs->pw_nghs = pw_sh_proc_mask = (PetscInt*) malloc(p_mask_size); 521827bd09bSSatish Balay gs->ngh_buf_sz = ngh_buf_size = p_mask_size*nel; 52252f87cdaSBarry Smith t_mask = (PetscInt*) malloc(p_mask_size); 52352f87cdaSBarry Smith gs->ngh_buf = ngh_buf = (PetscInt*) malloc(ngh_buf_size); 524827bd09bSSatish Balay 525827bd09bSSatish Balay /* comm buffer size ... memory usage bounded by ~2*msg_buf */ 526827bd09bSSatish Balay /* had thought I could exploit rendezvous threshold */ 527827bd09bSSatish Balay 528827bd09bSSatish Balay /* default is one pass */ 529827bd09bSSatish Balay per_load = negl = gs->negl; 530827bd09bSSatish Balay gs->num_loads = num_loads = 1; 531827bd09bSSatish Balay i = p_mask_size*negl; 532827bd09bSSatish Balay 533827bd09bSSatish Balay /* possible overflow on buffer size */ 534827bd09bSSatish Balay /* overflow hack */ 535*2fa5cd67SKarl Rupp if (i<0) i=INT_MAX; 536827bd09bSSatish Balay 53739945688SSatish Balay buf_size = PetscMin(msg_buf,i); 538827bd09bSSatish Balay 539827bd09bSSatish Balay /* can we do it? */ 540e32f2f54SBarry Smith 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); 541827bd09bSSatish Balay 542b1c944f5SJed Brown /* get PCTFS_giop buf space ... make *only* one malloc */ 54352f87cdaSBarry Smith buf1 = (PetscInt*) malloc(buf_size<<1); 544827bd09bSSatish Balay 545827bd09bSSatish Balay /* more than one gior exchange needed? */ 546db4deed7SKarl Rupp if (buf_size!=i) { 547827bd09bSSatish Balay per_load = buf_size/p_mask_size; 548827bd09bSSatish Balay buf_size = per_load*p_mask_size; 549827bd09bSSatish Balay gs->num_loads = num_loads = negl/per_load + (negl%per_load>0); 550827bd09bSSatish Balay } 551827bd09bSSatish Balay 552827bd09bSSatish Balay 553827bd09bSSatish Balay /* convert buf sizes from #bytes to #ints - 32 bit only! */ 554a501084fSBarry Smith p_mask_size/=sizeof(PetscInt); ngh_buf_size/=sizeof(PetscInt); buf_size/=sizeof(PetscInt); 555827bd09bSSatish Balay 556b1c944f5SJed Brown /* find PCTFS_giop work space */ 557827bd09bSSatish Balay buf2 = buf1+buf_size; 558827bd09bSSatish Balay 559827bd09bSSatish Balay /* hold #ints needed for processor masks */ 560827bd09bSSatish Balay gs->mask_sz=p_mask_size; 561827bd09bSSatish Balay 562827bd09bSSatish Balay /* init buffers */ 563ca8e9878SJed Brown ierr = PCTFS_ivec_zero(sh_proc_mask,p_mask_size);CHKERRQ(ierr); 564ca8e9878SJed Brown ierr = PCTFS_ivec_zero(pw_sh_proc_mask,p_mask_size);CHKERRQ(ierr); 565ca8e9878SJed Brown ierr = PCTFS_ivec_zero(ngh_buf,ngh_buf_size);CHKERRQ(ierr); 566827bd09bSSatish Balay 567827bd09bSSatish Balay /* HACK reset tree info */ 568827bd09bSSatish Balay tree_buf = NULL; 569827bd09bSSatish Balay tree_buf_sz = ntree = 0; 570827bd09bSSatish Balay 571827bd09bSSatish Balay /* ok do it */ 572db4deed7SKarl Rupp for (ptr1=ngh_buf,ptr2=elms,end=gs->gl_min,or_ct=i=0; or_ct<num_loads; or_ct++) { 573827bd09bSSatish Balay /* identity for bitwise or is 000...000 */ 574ca8e9878SJed Brown PCTFS_ivec_zero(buf1,buf_size); 575827bd09bSSatish Balay 576827bd09bSSatish Balay /* load msg buffer */ 577db4deed7SKarl Rupp for (start=end,end+=per_load,i_start=i; (offset=*ptr2)<end; i++, ptr2++) { 578827bd09bSSatish Balay offset = (offset-start)*p_mask_size; 579ca8e9878SJed Brown PCTFS_ivec_copy(buf1+offset,p_mask,p_mask_size); 580827bd09bSSatish Balay } 581827bd09bSSatish Balay 582827bd09bSSatish Balay /* GLOBAL: pass buffer */ 583b1c944f5SJed Brown ierr = PCTFS_giop(buf1,buf2,buf_size,&oper);CHKERRQ(ierr); 584827bd09bSSatish Balay 585827bd09bSSatish Balay 586827bd09bSSatish Balay /* unload buffer into ngh_buf */ 587827bd09bSSatish Balay ptr2=(elms+i_start); 588db4deed7SKarl Rupp for (ptr3=buf1,j=start; j<end; ptr3+=p_mask_size,j++) { 589827bd09bSSatish Balay /* I own it ... may have to pairwise it */ 590db4deed7SKarl Rupp if (j==*ptr2) { 591827bd09bSSatish Balay /* do i share it w/anyone? */ 592ca8e9878SJed Brown ct1 = PCTFS_ct_bits((char*)ptr3,p_mask_size*sizeof(PetscInt)); 593827bd09bSSatish Balay /* guess not */ 594db4deed7SKarl Rupp if (ct1<2) { ptr2++; ptr1+=p_mask_size; continue; } 595827bd09bSSatish Balay 596827bd09bSSatish Balay /* i do ... so keep info and turn off my bit */ 597ca8e9878SJed Brown PCTFS_ivec_copy(ptr1,ptr3,p_mask_size); 598ca8e9878SJed Brown ierr = PCTFS_ivec_xor(ptr1,p_mask,p_mask_size);CHKERRQ(ierr); 599ca8e9878SJed Brown ierr = PCTFS_ivec_or(sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr); 600827bd09bSSatish Balay 601827bd09bSSatish Balay /* is it to be done pairwise? */ 602db4deed7SKarl Rupp if (--ct1<=level) { 603827bd09bSSatish Balay npw++; 604827bd09bSSatish Balay 605827bd09bSSatish Balay /* turn on high bit to indicate pw need to process */ 606827bd09bSSatish Balay *ptr2++ |= TOP_BIT; 607ca8e9878SJed Brown ierr = PCTFS_ivec_or(pw_sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr); 608827bd09bSSatish Balay ptr1 += p_mask_size; 609827bd09bSSatish Balay continue; 610827bd09bSSatish Balay } 611827bd09bSSatish Balay 612827bd09bSSatish Balay /* get set for next and note that I have a tree contribution */ 613827bd09bSSatish Balay /* could save exact elm index for tree here -> save a search */ 614827bd09bSSatish Balay ptr2++; ptr1+=p_mask_size; ntree_map++; 615db4deed7SKarl Rupp } else { /* i don't but still might be involved in tree */ 616827bd09bSSatish Balay 617827bd09bSSatish Balay /* shared by how many? */ 618ca8e9878SJed Brown ct1 = PCTFS_ct_bits((char*)ptr3,p_mask_size*sizeof(PetscInt)); 619827bd09bSSatish Balay 620827bd09bSSatish Balay /* none! */ 621f1ed62a8SBarry Smith if (ct1<2) continue; 622827bd09bSSatish Balay 623827bd09bSSatish Balay /* is it going to be done pairwise? but not by me of course!*/ 624f1ed62a8SBarry Smith if (--ct1<=level) continue; 625827bd09bSSatish Balay } 626827bd09bSSatish Balay /* LATER we're going to have to process it NOW */ 627827bd09bSSatish Balay /* nope ... tree it */ 628f1ed62a8SBarry Smith ierr = place_in_tree(j);CHKERRQ(ierr); 629827bd09bSSatish Balay } 630827bd09bSSatish Balay } 631827bd09bSSatish Balay 632a501084fSBarry Smith free((void*)t_mask); 633a501084fSBarry Smith free((void*)buf1); 634827bd09bSSatish Balay 635827bd09bSSatish Balay gs->len_pw_list = npw; 636ca8e9878SJed Brown gs->num_nghs = PCTFS_ct_bits((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt)); 637827bd09bSSatish Balay 638827bd09bSSatish Balay /* expand from bit mask list to int list and save ngh list */ 63952f87cdaSBarry Smith gs->nghs = (PetscInt*) malloc(gs->num_nghs * sizeof(PetscInt)); 640ca8e9878SJed Brown PCTFS_bm_to_proc((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt),gs->nghs); 641827bd09bSSatish Balay 642ca8e9878SJed Brown gs->num_pw_nghs = PCTFS_ct_bits((char*)pw_sh_proc_mask,p_mask_size*sizeof(PetscInt)); 643827bd09bSSatish Balay 644827bd09bSSatish Balay oper = GL_MAX; 645827bd09bSSatish Balay ct1 = gs->num_nghs; 646b1c944f5SJed Brown ierr = PCTFS_giop(&ct1,&ct2,1,&oper);CHKERRQ(ierr); 647827bd09bSSatish Balay gs->max_nghs = ct1; 648827bd09bSSatish Balay 649827bd09bSSatish Balay gs->tree_map_sz = ntree_map; 650827bd09bSSatish Balay gs->max_left_over=ntree; 651827bd09bSSatish Balay 652a501084fSBarry Smith free((void*)p_mask); 653a501084fSBarry Smith free((void*)sh_proc_mask); 6543fdc5746SBarry Smith PetscFunctionReturn(0); 655827bd09bSSatish Balay } 656827bd09bSSatish Balay 657f1ed62a8SBarry Smith /******************************************************************************/ 658ca8e9878SJed Brown static PetscErrorCode set_pairwise(PCTFS_gs_id *gs) 659827bd09bSSatish Balay { 66052f87cdaSBarry Smith PetscInt i, j; 66152f87cdaSBarry Smith PetscInt p_mask_size; 66252f87cdaSBarry Smith PetscInt *p_mask, *sh_proc_mask, *tmp_proc_mask; 66352f87cdaSBarry Smith PetscInt *ngh_buf, *buf2; 66452f87cdaSBarry Smith PetscInt offset; 66552f87cdaSBarry Smith PetscInt *msg_list, *msg_size, **msg_nodes, nprs; 66652f87cdaSBarry Smith PetscInt *pairwise_elm_list, len_pair_list=0; 66752f87cdaSBarry Smith PetscInt *iptr, t1, i_start, nel, *elms; 66852f87cdaSBarry Smith PetscInt ct; 669f1ed62a8SBarry Smith PetscErrorCode ierr; 670827bd09bSSatish Balay 6713fdc5746SBarry Smith PetscFunctionBegin; 672827bd09bSSatish Balay /* to make life easier */ 673827bd09bSSatish Balay nel = gs->nel; 674827bd09bSSatish Balay elms = gs->elms; 675827bd09bSSatish Balay ngh_buf = gs->ngh_buf; 676827bd09bSSatish Balay sh_proc_mask = gs->pw_nghs; 677827bd09bSSatish Balay 678827bd09bSSatish Balay /* need a few temp masks */ 679ca8e9878SJed Brown p_mask_size = PCTFS_len_bit_mask(PCTFS_num_nodes); 68052f87cdaSBarry Smith p_mask = (PetscInt*) malloc(p_mask_size); 68152f87cdaSBarry Smith tmp_proc_mask = (PetscInt*) malloc(p_mask_size); 682827bd09bSSatish Balay 683b1c944f5SJed Brown /* set mask to my PCTFS_my_id's bit mask */ 684ca8e9878SJed Brown ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr); 685827bd09bSSatish Balay 686a501084fSBarry Smith p_mask_size /= sizeof(PetscInt); 687827bd09bSSatish Balay 688827bd09bSSatish Balay len_pair_list = gs->len_pw_list; 68952f87cdaSBarry Smith gs->pw_elm_list = pairwise_elm_list=(PetscInt*)malloc((len_pair_list+1)*sizeof(PetscInt)); 690827bd09bSSatish Balay 691827bd09bSSatish Balay /* how many processors (nghs) do we have to exchange with? */ 692ca8e9878SJed Brown nprs = gs->num_pairs = PCTFS_ct_bits((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt)); 693827bd09bSSatish Balay 694827bd09bSSatish Balay 695ca8e9878SJed Brown /* allocate space for PCTFS_gs_gop() info */ 69652f87cdaSBarry Smith gs->pair_list = msg_list = (PetscInt*) malloc(sizeof(PetscInt)*nprs); 69752f87cdaSBarry Smith gs->msg_sizes = msg_size = (PetscInt*) malloc(sizeof(PetscInt)*nprs); 69852f87cdaSBarry Smith gs->node_list = msg_nodes = (PetscInt**) malloc(sizeof(PetscInt*)*(nprs+1)); 699827bd09bSSatish Balay 700827bd09bSSatish Balay /* init msg_size list */ 701ca8e9878SJed Brown ierr = PCTFS_ivec_zero(msg_size,nprs);CHKERRQ(ierr); 702827bd09bSSatish Balay 703827bd09bSSatish Balay /* expand from bit mask list to int list */ 704ca8e9878SJed Brown ierr = PCTFS_bm_to_proc((char*)sh_proc_mask,p_mask_size*sizeof(PetscInt),msg_list);CHKERRQ(ierr); 705827bd09bSSatish Balay 706827bd09bSSatish Balay /* keep list of elements being handled pairwise */ 707db4deed7SKarl Rupp for (i=j=0; i<nel; i++) { 708db4deed7SKarl Rupp if (elms[i] & TOP_BIT) { elms[i] ^= TOP_BIT; pairwise_elm_list[j++] = i; } 709827bd09bSSatish Balay } 710827bd09bSSatish Balay pairwise_elm_list[j] = -1; 711827bd09bSSatish Balay 712a501084fSBarry Smith gs->msg_ids_out = (MPI_Request*) malloc(sizeof(MPI_Request)*(nprs+1)); 713827bd09bSSatish Balay gs->msg_ids_out[nprs] = MPI_REQUEST_NULL; 714a501084fSBarry Smith gs->msg_ids_in = (MPI_Request*) malloc(sizeof(MPI_Request)*(nprs+1)); 715827bd09bSSatish Balay gs->msg_ids_in[nprs] = MPI_REQUEST_NULL; 716a501084fSBarry Smith gs->pw_vals = (PetscScalar*) malloc(sizeof(PetscScalar)*len_pair_list*vec_sz); 717827bd09bSSatish Balay 718827bd09bSSatish Balay /* find who goes to each processor */ 719db4deed7SKarl Rupp for (i_start=i=0; i<nprs; i++) { 720827bd09bSSatish Balay /* processor i's mask */ 721ca8e9878SJed Brown ierr = PCTFS_set_bit_mask(p_mask,p_mask_size*sizeof(PetscInt),msg_list[i]);CHKERRQ(ierr); 722827bd09bSSatish Balay 723827bd09bSSatish Balay /* det # going to processor i */ 724db4deed7SKarl Rupp for (ct=j=0; j<len_pair_list; j++) { 725827bd09bSSatish Balay buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size); 726ca8e9878SJed Brown ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr); 727*2fa5cd67SKarl Rupp if (PCTFS_ct_bits((char*)tmp_proc_mask,p_mask_size*sizeof(PetscInt))) ct++; 728827bd09bSSatish Balay } 729827bd09bSSatish Balay msg_size[i] = ct; 73039945688SSatish Balay i_start = PetscMax(i_start,ct); 731827bd09bSSatish Balay 732827bd09bSSatish Balay /*space to hold nodes in message to first neighbor */ 73352f87cdaSBarry Smith msg_nodes[i] = iptr = (PetscInt*) malloc(sizeof(PetscInt)*(ct+1)); 734827bd09bSSatish Balay 735db4deed7SKarl Rupp for (j=0;j<len_pair_list;j++) { 736827bd09bSSatish Balay buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size); 737ca8e9878SJed Brown ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr); 738*2fa5cd67SKarl Rupp if (PCTFS_ct_bits((char*)tmp_proc_mask,p_mask_size*sizeof(PetscInt))) *iptr++ = j; 739827bd09bSSatish Balay } 740827bd09bSSatish Balay *iptr = -1; 741827bd09bSSatish Balay } 742827bd09bSSatish Balay msg_nodes[nprs] = NULL; 743827bd09bSSatish Balay 744827bd09bSSatish Balay j = gs->loc_node_pairs=i_start; 745827bd09bSSatish Balay t1 = GL_MAX; 746b1c944f5SJed Brown ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 747827bd09bSSatish Balay gs->max_node_pairs = i_start; 748827bd09bSSatish Balay 749827bd09bSSatish Balay i_start = j; 750827bd09bSSatish Balay t1 = GL_MIN; 751b1c944f5SJed Brown ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 752827bd09bSSatish Balay gs->min_node_pairs = i_start; 753827bd09bSSatish Balay 754827bd09bSSatish Balay i_start = j; 755827bd09bSSatish Balay t1 = GL_ADD; 756b1c944f5SJed Brown ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 757b1c944f5SJed Brown gs->avg_node_pairs = i_start/PCTFS_num_nodes + 1; 758827bd09bSSatish Balay 759827bd09bSSatish Balay i_start = nprs; 760827bd09bSSatish Balay t1 = GL_MAX; 761b1c944f5SJed Brown PCTFS_giop(&i_start,&offset,1,&t1); 762827bd09bSSatish Balay gs->max_pairs = i_start; 763827bd09bSSatish Balay 764827bd09bSSatish Balay 765827bd09bSSatish Balay /* remap pairwise in tail of gsi_via_bit_mask() */ 766ca8e9878SJed Brown gs->msg_total = PCTFS_ivec_sum(gs->msg_sizes,nprs); 767a501084fSBarry Smith gs->out = (PetscScalar*) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz); 768a501084fSBarry Smith gs->in = (PetscScalar*) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz); 769827bd09bSSatish Balay 770827bd09bSSatish Balay /* reset malloc pool */ 771a501084fSBarry Smith free((void*)p_mask); 772a501084fSBarry Smith free((void*)tmp_proc_mask); 7733fdc5746SBarry Smith PetscFunctionReturn(0); 774827bd09bSSatish Balay } 775827bd09bSSatish Balay 776f1ed62a8SBarry Smith /* to do pruned tree just save ngh buf copy for each one and decode here! 777827bd09bSSatish Balay ******************************************************************************/ 778ca8e9878SJed Brown static PetscErrorCode set_tree(PCTFS_gs_id *gs) 779827bd09bSSatish Balay { 78052f87cdaSBarry Smith PetscInt i, j, n, nel; 78152f87cdaSBarry Smith PetscInt *iptr_in, *iptr_out, *tree_elms, *elms; 782827bd09bSSatish Balay 7833fdc5746SBarry Smith PetscFunctionBegin; 784827bd09bSSatish Balay /* local work ptrs */ 785827bd09bSSatish Balay elms = gs->elms; 786827bd09bSSatish Balay nel = gs->nel; 787827bd09bSSatish Balay 788827bd09bSSatish Balay /* how many via tree */ 789827bd09bSSatish Balay gs->tree_nel = n = ntree; 790827bd09bSSatish Balay gs->tree_elms = tree_elms = iptr_in = tree_buf; 791a501084fSBarry Smith gs->tree_buf = (PetscScalar*) malloc(sizeof(PetscScalar)*n*vec_sz); 792a501084fSBarry Smith gs->tree_work = (PetscScalar*) malloc(sizeof(PetscScalar)*n*vec_sz); 793827bd09bSSatish Balay j = gs->tree_map_sz; 79452f87cdaSBarry Smith gs->tree_map_in = iptr_in = (PetscInt*) malloc(sizeof(PetscInt)*(j+1)); 79552f87cdaSBarry Smith gs->tree_map_out = iptr_out = (PetscInt*) malloc(sizeof(PetscInt)*(j+1)); 796827bd09bSSatish Balay 797827bd09bSSatish Balay /* search the longer of the two lists */ 798827bd09bSSatish Balay /* note ... could save this info in get_ngh_buf and save searches */ 799db4deed7SKarl Rupp if (n<=nel) { 800827bd09bSSatish Balay /* bijective fct w/remap - search elm list */ 801db4deed7SKarl Rupp for (i=0; i<n; i++) { 802db4deed7SKarl Rupp if ((j=PCTFS_ivec_binary_search(*tree_elms++,elms,nel))>=0) {*iptr_in++ = j; *iptr_out++ = i;} 803827bd09bSSatish Balay } 804db4deed7SKarl Rupp } else { 805db4deed7SKarl Rupp for (i=0; i<nel; i++) { 806db4deed7SKarl Rupp if ((j=PCTFS_ivec_binary_search(*elms++,tree_elms,n))>=0) {*iptr_in++ = i; *iptr_out++ = j;} 807827bd09bSSatish Balay } 808827bd09bSSatish Balay } 809827bd09bSSatish Balay 810827bd09bSSatish Balay /* sentinel */ 811827bd09bSSatish Balay *iptr_in = *iptr_out = -1; 8123fdc5746SBarry Smith PetscFunctionReturn(0); 813827bd09bSSatish Balay } 814827bd09bSSatish Balay 815f1ed62a8SBarry Smith /******************************************************************************/ 816ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_out(PCTFS_gs_id *gs, PetscScalar *vals) 817827bd09bSSatish Balay { 81852f87cdaSBarry Smith PetscInt *num, *map, **reduce; 819a501084fSBarry Smith PetscScalar tmp; 820827bd09bSSatish Balay 8213fdc5746SBarry Smith PetscFunctionBegin; 822827bd09bSSatish Balay num = gs->num_gop_local_reduce; 823827bd09bSSatish Balay reduce = gs->gop_local_reduce; 824db4deed7SKarl Rupp while ((map = *reduce++)) { 825827bd09bSSatish Balay /* wall */ 826db4deed7SKarl Rupp if (*num == 2) { 827827bd09bSSatish Balay num++; 828827bd09bSSatish Balay vals[map[1]] = vals[map[0]]; 829db4deed7SKarl Rupp } else if (*num == 3) { /* corner shared by three elements */ 830827bd09bSSatish Balay num++; 831827bd09bSSatish Balay vals[map[2]] = vals[map[1]] = vals[map[0]]; 832db4deed7SKarl Rupp } else if (*num == 4) { /* corner shared by four elements */ 833827bd09bSSatish Balay num++; 834827bd09bSSatish Balay vals[map[3]] = vals[map[2]] = vals[map[1]] = vals[map[0]]; 835db4deed7SKarl Rupp } else { /* general case ... odd geoms ... 3D*/ 836827bd09bSSatish Balay num++; 837827bd09bSSatish Balay tmp = *(vals + *map++); 838*2fa5cd67SKarl Rupp while (*map >= 0) *(vals + *map++) = tmp; 839827bd09bSSatish Balay } 840827bd09bSSatish Balay } 8413fdc5746SBarry Smith PetscFunctionReturn(0); 842827bd09bSSatish Balay } 843827bd09bSSatish Balay 8447b1ae94cSBarry Smith /******************************************************************************/ 845ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_plus(PCTFS_gs_id *gs, PetscScalar *vals) 846827bd09bSSatish Balay { 84752f87cdaSBarry Smith PetscInt *num, *map, **reduce; 848a501084fSBarry Smith PetscScalar tmp; 849827bd09bSSatish Balay 8503fdc5746SBarry Smith PetscFunctionBegin; 851827bd09bSSatish Balay num = gs->num_local_reduce; 852827bd09bSSatish Balay reduce = gs->local_reduce; 853db4deed7SKarl Rupp while ((map = *reduce)) { 854827bd09bSSatish Balay /* wall */ 855db4deed7SKarl Rupp if (*num == 2) { 856827bd09bSSatish Balay num++; reduce++; 857827bd09bSSatish Balay vals[map[1]] = vals[map[0]] += vals[map[1]]; 858db4deed7SKarl Rupp } else if (*num == 3) { /* corner shared by three elements */ 859827bd09bSSatish Balay num++; reduce++; 860827bd09bSSatish Balay vals[map[2]]=vals[map[1]]=vals[map[0]]+=(vals[map[1]]+vals[map[2]]); 861db4deed7SKarl Rupp } else if (*num == 4) { /* corner shared by four elements */ 862827bd09bSSatish Balay num++; reduce++; 863*2fa5cd67SKarl Rupp vals[map[1]]=vals[map[2]]=vals[map[3]]=vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]); 864db4deed7SKarl Rupp } else { /* general case ... odd geoms ... 3D*/ 865827bd09bSSatish Balay num++; 866827bd09bSSatish Balay tmp = 0.0; 867*2fa5cd67SKarl Rupp while (*map >= 0) tmp += *(vals + *map++); 868827bd09bSSatish Balay 869827bd09bSSatish Balay map = *reduce++; 870*2fa5cd67SKarl Rupp while (*map >= 0) *(vals + *map++) = tmp; 871827bd09bSSatish Balay } 872827bd09bSSatish Balay } 8733fdc5746SBarry Smith PetscFunctionReturn(0); 874827bd09bSSatish Balay } 875827bd09bSSatish Balay 8767b1ae94cSBarry Smith /******************************************************************************/ 877ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals) 878827bd09bSSatish Balay { 87952f87cdaSBarry Smith PetscInt *num, *map, **reduce; 880a501084fSBarry Smith PetscScalar *base; 881827bd09bSSatish Balay 8823fdc5746SBarry Smith PetscFunctionBegin; 883827bd09bSSatish Balay num = gs->num_gop_local_reduce; 884827bd09bSSatish Balay reduce = gs->gop_local_reduce; 885db4deed7SKarl Rupp while ((map = *reduce++)) { 886827bd09bSSatish Balay /* wall */ 887db4deed7SKarl Rupp if (*num == 2) { 888827bd09bSSatish Balay num++; 889827bd09bSSatish Balay vals[map[0]] += vals[map[1]]; 890db4deed7SKarl Rupp } else if (*num == 3) { /* corner shared by three elements */ 891827bd09bSSatish Balay num++; 892827bd09bSSatish Balay vals[map[0]] += (vals[map[1]] + vals[map[2]]); 893db4deed7SKarl Rupp } else if (*num == 4) { /* corner shared by four elements */ 894827bd09bSSatish Balay num++; 895827bd09bSSatish Balay vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]); 896db4deed7SKarl Rupp } else { /* general case ... odd geoms ... 3D*/ 897827bd09bSSatish Balay num++; 898827bd09bSSatish Balay base = vals + *map++; 899*2fa5cd67SKarl Rupp while (*map >= 0) *base += *(vals + *map++); 900827bd09bSSatish Balay } 901827bd09bSSatish Balay } 9023fdc5746SBarry Smith PetscFunctionReturn(0); 903827bd09bSSatish Balay } 904827bd09bSSatish Balay 9057b1ae94cSBarry Smith /******************************************************************************/ 906ca8e9878SJed Brown PetscErrorCode PCTFS_gs_free(PCTFS_gs_id *gs) 907827bd09bSSatish Balay { 90852f87cdaSBarry Smith PetscInt i; 909827bd09bSSatish Balay 9103fdc5746SBarry Smith PetscFunctionBegin; 911*2fa5cd67SKarl Rupp if (gs->nghs) free((void*) gs->nghs); 912*2fa5cd67SKarl Rupp if (gs->pw_nghs) free((void*) gs->pw_nghs); 913827bd09bSSatish Balay 914827bd09bSSatish Balay /* tree */ 915*2fa5cd67SKarl Rupp if (gs->max_left_over) { 916*2fa5cd67SKarl Rupp if (gs->tree_elms) free((void*) gs->tree_elms); 917*2fa5cd67SKarl Rupp if (gs->tree_buf) free((void*) gs->tree_buf); 918*2fa5cd67SKarl Rupp if (gs->tree_work) free((void*) gs->tree_work); 919*2fa5cd67SKarl Rupp if (gs->tree_map_in) free((void*) gs->tree_map_in); 920*2fa5cd67SKarl Rupp if (gs->tree_map_out) free((void*) gs->tree_map_out); 921827bd09bSSatish Balay } 922827bd09bSSatish Balay 923827bd09bSSatish Balay /* pairwise info */ 924*2fa5cd67SKarl Rupp if (gs->num_pairs) { 925827bd09bSSatish Balay /* should be NULL already */ 926*2fa5cd67SKarl Rupp if (gs->ngh_buf) free((void*) gs->ngh_buf); 927*2fa5cd67SKarl Rupp if (gs->elms) free((void*) gs->elms); 928*2fa5cd67SKarl Rupp if (gs->local_elms) free((void*) gs->local_elms); 929*2fa5cd67SKarl Rupp if (gs->companion) free((void*) gs->companion); 930827bd09bSSatish Balay 931827bd09bSSatish Balay /* only set if pairwise */ 932*2fa5cd67SKarl Rupp if (gs->vals) free((void*) gs->vals); 933*2fa5cd67SKarl Rupp if (gs->in) free((void*) gs->in); 934*2fa5cd67SKarl Rupp if (gs->out) free((void*) gs->out); 935*2fa5cd67SKarl Rupp if (gs->msg_ids_in) free((void*) gs->msg_ids_in); 936*2fa5cd67SKarl Rupp if (gs->msg_ids_out) free((void*) gs->msg_ids_out); 937*2fa5cd67SKarl Rupp if (gs->pw_vals) free((void*) gs->pw_vals); 938*2fa5cd67SKarl Rupp if (gs->pw_elm_list) free((void*) gs->pw_elm_list); 939db4deed7SKarl Rupp if (gs->node_list) { 940db4deed7SKarl Rupp for (i=0;i<gs->num_pairs;i++) { 941db4deed7SKarl Rupp if (gs->node_list[i]) { 942db4deed7SKarl Rupp free((void*) gs->node_list[i]); 943db4deed7SKarl Rupp } 944db4deed7SKarl Rupp } 945a501084fSBarry Smith free((void*) gs->node_list); 946827bd09bSSatish Balay } 947*2fa5cd67SKarl Rupp if (gs->msg_sizes) free((void*) gs->msg_sizes); 948*2fa5cd67SKarl Rupp if (gs->pair_list) free((void*) gs->pair_list); 949827bd09bSSatish Balay } 950827bd09bSSatish Balay 951827bd09bSSatish Balay /* local info */ 952db4deed7SKarl Rupp if (gs->num_local_total>=0) { 953827bd09bSSatish Balay /* for (i=0;i<gs->num_local_total;i++) */ 954db4deed7SKarl Rupp for (i=0;i<gs->num_local_total+1;i++) { 955*2fa5cd67SKarl Rupp if (gs->num_gop_local_reduce[i]) free((void*) gs->gop_local_reduce[i]); 956827bd09bSSatish Balay } 957827bd09bSSatish Balay } 958827bd09bSSatish Balay 959827bd09bSSatish Balay /* if intersection tree/pairwise and local isn't empty */ 960*2fa5cd67SKarl Rupp if (gs->gop_local_reduce) free((void*) gs->gop_local_reduce); 961*2fa5cd67SKarl Rupp if (gs->num_gop_local_reduce) free((void*) gs->num_gop_local_reduce); 962827bd09bSSatish Balay 963a501084fSBarry Smith free((void*) gs); 9643fdc5746SBarry Smith PetscFunctionReturn(0); 965827bd09bSSatish Balay } 966827bd09bSSatish Balay 9677b1ae94cSBarry Smith /******************************************************************************/ 968ca8e9878SJed Brown PetscErrorCode PCTFS_gs_gop_vec(PCTFS_gs_id *gs, PetscScalar *vals, const char *op, PetscInt step) 969827bd09bSSatish Balay { 970d1528f56SBarry Smith PetscErrorCode ierr; 971d1528f56SBarry Smith 9723fdc5746SBarry Smith PetscFunctionBegin; 973827bd09bSSatish Balay switch (*op) { 974827bd09bSSatish Balay case '+': 975ca8e9878SJed Brown PCTFS_gs_gop_vec_plus(gs,vals,step); 976827bd09bSSatish Balay break; 977827bd09bSSatish Balay default: 978ca8e9878SJed Brown ierr = PetscInfo1(0,"PCTFS_gs_gop_vec() :: %c is not a valid op",op[0]);CHKERRQ(ierr); 979ca8e9878SJed Brown ierr = PetscInfo(0,"PCTFS_gs_gop_vec() :: default :: plus");CHKERRQ(ierr); 980ca8e9878SJed Brown PCTFS_gs_gop_vec_plus(gs,vals,step); 981827bd09bSSatish Balay break; 982827bd09bSSatish Balay } 9833fdc5746SBarry Smith PetscFunctionReturn(0); 984827bd09bSSatish Balay } 985827bd09bSSatish Balay 9867b1ae94cSBarry Smith /******************************************************************************/ 987ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 988827bd09bSSatish Balay { 9893fdc5746SBarry Smith PetscFunctionBegin; 990ca8e9878SJed Brown if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_gs_gop_vec() passed NULL gs handle!!!"); 991827bd09bSSatish Balay 992827bd09bSSatish Balay /* local only operations!!! */ 993*2fa5cd67SKarl Rupp if (gs->num_local) PCTFS_gs_gop_vec_local_plus(gs,vals,step); 994827bd09bSSatish Balay 995827bd09bSSatish Balay /* if intersection tree/pairwise and local isn't empty */ 996*2fa5cd67SKarl Rupp if (gs->num_local_gop) { 997ca8e9878SJed Brown PCTFS_gs_gop_vec_local_in_plus(gs,vals,step); 998827bd09bSSatish Balay 999827bd09bSSatish Balay /* pairwise */ 1000*2fa5cd67SKarl Rupp if (gs->num_pairs) PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step); 1001827bd09bSSatish Balay 1002827bd09bSSatish Balay /* tree */ 1003*2fa5cd67SKarl Rupp else if (gs->max_left_over) PCTFS_gs_gop_vec_tree_plus(gs,vals,step); 1004827bd09bSSatish Balay 1005ca8e9878SJed Brown PCTFS_gs_gop_vec_local_out(gs,vals,step); 1006db4deed7SKarl Rupp } else { /* if intersection tree/pairwise and local is empty */ 1007827bd09bSSatish Balay /* pairwise */ 1008*2fa5cd67SKarl Rupp if (gs->num_pairs) PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step); 1009827bd09bSSatish Balay 1010827bd09bSSatish Balay /* tree */ 1011*2fa5cd67SKarl Rupp else if (gs->max_left_over) PCTFS_gs_gop_vec_tree_plus(gs,vals,step); 1012827bd09bSSatish Balay } 10133fdc5746SBarry Smith PetscFunctionReturn(0); 1014827bd09bSSatish Balay } 1015827bd09bSSatish Balay 10167b1ae94cSBarry Smith /******************************************************************************/ 1017ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1018827bd09bSSatish Balay { 101952f87cdaSBarry Smith PetscInt *num, *map, **reduce; 1020a501084fSBarry Smith PetscScalar *base; 1021827bd09bSSatish Balay 10223fdc5746SBarry Smith PetscFunctionBegin; 1023827bd09bSSatish Balay num = gs->num_local_reduce; 1024827bd09bSSatish Balay reduce = gs->local_reduce; 1025db4deed7SKarl Rupp while ((map = *reduce)) { 1026827bd09bSSatish Balay base = vals + map[0] * step; 1027827bd09bSSatish Balay 1028827bd09bSSatish Balay /* wall */ 1029db4deed7SKarl Rupp if (*num == 2) { 1030827bd09bSSatish Balay num++; reduce++; 1031ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[1]*step,step); 1032ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1033db4deed7SKarl Rupp } else if (*num == 3) { /* corner shared by three elements */ 1034827bd09bSSatish Balay num++; reduce++; 1035ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[1]*step,step); 1036ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[2]*step,step); 1037ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1038ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1039db4deed7SKarl Rupp } else if (*num == 4) { /* corner shared by four elements */ 1040827bd09bSSatish Balay num++; reduce++; 1041ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[1]*step,step); 1042ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[2]*step,step); 1043ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[3]*step,step); 1044ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[3]*step,base,step); 1045ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1046ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1047db4deed7SKarl Rupp } else { /* general case ... odd geoms ... 3D */ 1048827bd09bSSatish Balay num++; 1049*2fa5cd67SKarl Rupp while (*++map >= 0) PCTFS_rvec_add (base,vals+*map*step,step); 1050827bd09bSSatish Balay 1051827bd09bSSatish Balay map = *reduce; 1052*2fa5cd67SKarl Rupp while (*++map >= 0) PCTFS_rvec_copy(vals+*map*step,base,step); 1053827bd09bSSatish Balay 1054827bd09bSSatish Balay reduce++; 1055827bd09bSSatish Balay } 1056827bd09bSSatish Balay } 10573fdc5746SBarry Smith PetscFunctionReturn(0); 1058827bd09bSSatish Balay } 1059827bd09bSSatish Balay 10607b1ae94cSBarry Smith /******************************************************************************/ 1061ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1062827bd09bSSatish Balay { 106352f87cdaSBarry Smith PetscInt *num, *map, **reduce; 1064a501084fSBarry Smith PetscScalar *base; 1065db4deed7SKarl Rupp 10663fdc5746SBarry Smith PetscFunctionBegin; 1067827bd09bSSatish Balay num = gs->num_gop_local_reduce; 1068827bd09bSSatish Balay reduce = gs->gop_local_reduce; 1069db4deed7SKarl Rupp while ((map = *reduce++)) { 1070827bd09bSSatish Balay base = vals + map[0] * step; 1071827bd09bSSatish Balay 1072827bd09bSSatish Balay /* wall */ 1073db4deed7SKarl Rupp if (*num == 2) { 1074827bd09bSSatish Balay num++; 1075ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[1]*step,step); 1076db4deed7SKarl Rupp } else if (*num == 3) { /* corner shared by three elements */ 1077827bd09bSSatish Balay num++; 1078ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[1]*step,step); 1079ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[2]*step,step); 1080db4deed7SKarl Rupp } else if (*num == 4) { /* corner shared by four elements */ 1081827bd09bSSatish Balay num++; 1082ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[1]*step,step); 1083ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[2]*step,step); 1084ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[3]*step,step); 1085db4deed7SKarl Rupp } else { /* general case ... odd geoms ... 3D*/ 1086827bd09bSSatish Balay num++; 1087*2fa5cd67SKarl Rupp while (*++map >= 0) PCTFS_rvec_add(base,vals+*map*step,step); 1088827bd09bSSatish Balay } 1089827bd09bSSatish Balay } 10903fdc5746SBarry Smith PetscFunctionReturn(0); 1091827bd09bSSatish Balay } 1092827bd09bSSatish Balay 10937b1ae94cSBarry Smith /******************************************************************************/ 1094ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_out(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1095827bd09bSSatish Balay { 109652f87cdaSBarry Smith PetscInt *num, *map, **reduce; 1097a501084fSBarry Smith PetscScalar *base; 1098827bd09bSSatish Balay 10993fdc5746SBarry Smith PetscFunctionBegin; 1100827bd09bSSatish Balay num = gs->num_gop_local_reduce; 1101827bd09bSSatish Balay reduce = gs->gop_local_reduce; 1102db4deed7SKarl Rupp while ((map = *reduce++)) { 1103827bd09bSSatish Balay base = vals + map[0] * step; 1104827bd09bSSatish Balay 1105827bd09bSSatish Balay /* wall */ 1106db4deed7SKarl Rupp if (*num == 2) { 1107827bd09bSSatish Balay num++; 1108ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1109db4deed7SKarl Rupp } else if (*num == 3) { /* corner shared by three elements */ 1110827bd09bSSatish Balay num++; 1111ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1112ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1113db4deed7SKarl Rupp } else if (*num == 4) { /* corner shared by four elements */ 1114827bd09bSSatish Balay num++; 1115ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1116ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1117ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[3]*step,base,step); 1118db4deed7SKarl Rupp } else { /* general case ... odd geoms ... 3D*/ 1119827bd09bSSatish Balay num++; 1120*2fa5cd67SKarl Rupp while (*++map >= 0) PCTFS_rvec_copy(vals+*map*step,base,step); 1121827bd09bSSatish Balay } 1122827bd09bSSatish Balay } 11233fdc5746SBarry Smith PetscFunctionReturn(0); 1124827bd09bSSatish Balay } 1125827bd09bSSatish Balay 11267b1ae94cSBarry Smith /******************************************************************************/ 1127ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step) 1128827bd09bSSatish Balay { 1129a501084fSBarry Smith PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2; 113052f87cdaSBarry Smith PetscInt *iptr, *msg_list, *msg_size, **msg_nodes; 113152f87cdaSBarry Smith PetscInt *pw, *list, *size, **nodes; 1132827bd09bSSatish Balay MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out; 1133827bd09bSSatish Balay MPI_Status status; 11340805154bSBarry Smith PetscBLASInt i1 = 1,dstep; 11353fdc5746SBarry Smith PetscErrorCode ierr; 1136827bd09bSSatish Balay 11373fdc5746SBarry Smith PetscFunctionBegin; 1138a501084fSBarry Smith /* strip and load s */ 1139827bd09bSSatish Balay msg_list = list = gs->pair_list; 1140827bd09bSSatish Balay msg_size = size = gs->msg_sizes; 1141827bd09bSSatish Balay msg_nodes = nodes = gs->node_list; 1142827bd09bSSatish Balay iptr = pw = gs->pw_elm_list; 1143827bd09bSSatish Balay dptr1 = dptr3 = gs->pw_vals; 1144827bd09bSSatish Balay msg_ids_in = ids_in = gs->msg_ids_in; 1145827bd09bSSatish Balay msg_ids_out = ids_out = gs->msg_ids_out; 1146827bd09bSSatish Balay dptr2 = gs->out; 1147827bd09bSSatish Balay in1=in2 = gs->in; 1148827bd09bSSatish Balay 1149827bd09bSSatish Balay /* post the receives */ 1150827bd09bSSatish Balay /* msg_nodes=nodes; */ 1151db4deed7SKarl Rupp do { 1152827bd09bSSatish Balay /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the 1153827bd09bSSatish Balay second one *list and do list++ afterwards */ 1154ca8e9878SJed Brown ierr = MPI_Irecv(in1, *size *step, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr); 11559182e22cSBarry Smith list++;msg_ids_in++; 1156827bd09bSSatish Balay in1 += *size++ *step; 1157*2fa5cd67SKarl Rupp } while (*++msg_nodes); 1158827bd09bSSatish Balay msg_nodes=nodes; 1159827bd09bSSatish Balay 1160827bd09bSSatish Balay /* load gs values into in out gs buffers */ 1161db4deed7SKarl Rupp while (*iptr >= 0) { 1162ca8e9878SJed Brown PCTFS_rvec_copy(dptr3,in_vals + *iptr*step,step); 1163827bd09bSSatish Balay dptr3+=step; 1164827bd09bSSatish Balay iptr++; 1165827bd09bSSatish Balay } 1166827bd09bSSatish Balay 1167827bd09bSSatish Balay /* load out buffers and post the sends */ 1168db4deed7SKarl Rupp while ((iptr = *msg_nodes++)) { 1169827bd09bSSatish Balay dptr3 = dptr2; 1170db4deed7SKarl Rupp while (*iptr >= 0) { 1171ca8e9878SJed Brown PCTFS_rvec_copy(dptr2,dptr1 + *iptr*step,step); 1172827bd09bSSatish Balay dptr2+=step; 1173827bd09bSSatish Balay iptr++; 1174827bd09bSSatish Balay } 1175ca8e9878SJed Brown ierr = MPI_Isend(dptr3, *msg_size *step, MPIU_SCALAR, *msg_list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr); 11769182e22cSBarry Smith msg_size++; msg_list++;msg_ids_out++; 1177827bd09bSSatish Balay } 1178827bd09bSSatish Balay 1179827bd09bSSatish Balay /* tree */ 1180*2fa5cd67SKarl Rupp if (gs->max_left_over) PCTFS_gs_gop_vec_tree_plus(gs,in_vals,step); 1181827bd09bSSatish Balay 1182827bd09bSSatish Balay /* process the received data */ 1183827bd09bSSatish Balay msg_nodes=nodes; 1184a501084fSBarry Smith while ((iptr = *nodes++)) { 1185a501084fSBarry Smith PetscScalar d1 = 1.0; 1186db4deed7SKarl Rupp 1187827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1188827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 11899182e22cSBarry Smith ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr); 11909182e22cSBarry Smith ids_in++; 1191a501084fSBarry Smith while (*iptr >= 0) { 1192c5df96a5SBarry Smith ierr = PetscBLASIntCast(step,&dstep);CHKERRQ(ierr); 11934a0de3f6SBarry Smith BLASaxpy_(&dstep,&d1,in2,&i1,dptr1 + *iptr*step,&i1); 1194827bd09bSSatish Balay in2+=step; 1195827bd09bSSatish Balay iptr++; 1196827bd09bSSatish Balay } 1197827bd09bSSatish Balay } 1198827bd09bSSatish Balay 1199827bd09bSSatish Balay /* replace vals */ 1200db4deed7SKarl Rupp while (*pw >= 0) { 1201ca8e9878SJed Brown PCTFS_rvec_copy(in_vals + *pw*step,dptr1,step); 1202827bd09bSSatish Balay dptr1+=step; 1203827bd09bSSatish Balay pw++; 1204827bd09bSSatish Balay } 1205827bd09bSSatish Balay 1206827bd09bSSatish Balay /* clear isend message handles */ 1207827bd09bSSatish Balay /* This changed for clarity though it could be the same */ 1208db4deed7SKarl Rupp 1209827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1210827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 1211*2fa5cd67SKarl Rupp while (*msg_nodes++) { 1212*2fa5cd67SKarl Rupp ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr); 1213*2fa5cd67SKarl Rupp ids_out++; 1214*2fa5cd67SKarl Rupp } 12153fdc5746SBarry Smith PetscFunctionReturn(0); 1216827bd09bSSatish Balay } 1217827bd09bSSatish Balay 12187b1ae94cSBarry Smith /******************************************************************************/ 1219ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_tree_plus(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1220827bd09bSSatish Balay { 122152f87cdaSBarry Smith PetscInt size, *in, *out; 1222a501084fSBarry Smith PetscScalar *buf, *work; 122352f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 1224a501084fSBarry Smith PetscBLASInt i1 = 1; 1225c5df96a5SBarry Smith PetscErrorCode ierr; 1226c5df96a5SBarry Smith PetscBLASInt dstep; 1227827bd09bSSatish Balay 12283fdc5746SBarry Smith PetscFunctionBegin; 1229827bd09bSSatish Balay /* copy over to local variables */ 1230827bd09bSSatish Balay in = gs->tree_map_in; 1231827bd09bSSatish Balay out = gs->tree_map_out; 1232827bd09bSSatish Balay buf = gs->tree_buf; 1233827bd09bSSatish Balay work = gs->tree_work; 1234827bd09bSSatish Balay size = gs->tree_nel*step; 1235827bd09bSSatish Balay 1236827bd09bSSatish Balay /* zero out collection buffer */ 1237ca8e9878SJed Brown PCTFS_rvec_zero(buf,size); 1238827bd09bSSatish Balay 1239827bd09bSSatish Balay 1240827bd09bSSatish Balay /* copy over my contributions */ 1241db4deed7SKarl Rupp while (*in >= 0) { 1242c5df96a5SBarry Smith ierr = PetscBLASIntCast(step,&dstep);CHKERRQ(ierr); 12436e4f4d19SBarry Smith BLAScopy_(&dstep,vals + *in++ * step,&i1,buf + *out++ * step,&i1); 1244827bd09bSSatish Balay } 1245827bd09bSSatish Balay 1246827bd09bSSatish Balay /* perform fan in/out on full buffer */ 1247b1c944f5SJed Brown /* must change PCTFS_grop to handle the blas */ 1248b1c944f5SJed Brown PCTFS_grop(buf,work,size,op); 1249827bd09bSSatish Balay 1250827bd09bSSatish Balay /* reset */ 1251827bd09bSSatish Balay in = gs->tree_map_in; 1252827bd09bSSatish Balay out = gs->tree_map_out; 1253827bd09bSSatish Balay 1254827bd09bSSatish Balay /* get the portion of the results I need */ 1255db4deed7SKarl Rupp while (*in >= 0) { 1256c5df96a5SBarry Smith ierr = PetscBLASIntCast(step,&dstep);CHKERRQ(ierr); 12576e4f4d19SBarry Smith BLAScopy_(&dstep,buf + *out++ * step,&i1,vals + *in++ * step,&i1); 1258827bd09bSSatish Balay } 12593fdc5746SBarry Smith PetscFunctionReturn(0); 1260827bd09bSSatish Balay } 1261827bd09bSSatish Balay 12627b1ae94cSBarry Smith /******************************************************************************/ 1263ca8e9878SJed Brown PetscErrorCode PCTFS_gs_gop_hc(PCTFS_gs_id *gs, PetscScalar *vals, const char *op, PetscInt dim) 1264827bd09bSSatish Balay { 1265d1528f56SBarry Smith PetscErrorCode ierr; 1266d1528f56SBarry Smith 12673fdc5746SBarry Smith PetscFunctionBegin; 1268827bd09bSSatish Balay switch (*op) { 1269827bd09bSSatish Balay case '+': 1270ca8e9878SJed Brown PCTFS_gs_gop_plus_hc(gs,vals,dim); 1271827bd09bSSatish Balay break; 1272827bd09bSSatish Balay default: 1273ca8e9878SJed Brown ierr = PetscInfo1(0,"PCTFS_gs_gop_hc() :: %c is not a valid op",op[0]);CHKERRQ(ierr); 1274ca8e9878SJed Brown ierr = PetscInfo(0,"PCTFS_gs_gop_hc() :: default :: plus\n");CHKERRQ(ierr); 1275ca8e9878SJed Brown PCTFS_gs_gop_plus_hc(gs,vals,dim); 1276827bd09bSSatish Balay break; 1277827bd09bSSatish Balay } 12783fdc5746SBarry Smith PetscFunctionReturn(0); 1279827bd09bSSatish Balay } 1280827bd09bSSatish Balay 12817b1ae94cSBarry Smith /******************************************************************************/ 1282ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim) 1283827bd09bSSatish Balay { 12843fdc5746SBarry Smith PetscFunctionBegin; 1285827bd09bSSatish Balay /* if there's nothing to do return */ 1286*2fa5cd67SKarl Rupp if (dim<=0) PetscFunctionReturn(0); 1287827bd09bSSatish Balay 1288827bd09bSSatish Balay /* can't do more dimensions then exist */ 1289b1c944f5SJed Brown dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 1290827bd09bSSatish Balay 1291827bd09bSSatish Balay /* local only operations!!! */ 1292*2fa5cd67SKarl Rupp if (gs->num_local) PCTFS_gs_gop_local_plus(gs,vals); 1293827bd09bSSatish Balay 1294827bd09bSSatish Balay /* if intersection tree/pairwise and local isn't empty */ 1295db4deed7SKarl Rupp if (gs->num_local_gop) { 1296ca8e9878SJed Brown PCTFS_gs_gop_local_in_plus(gs,vals); 1297827bd09bSSatish Balay 1298827bd09bSSatish Balay /* pairwise will do tree inside ... */ 1299*2fa5cd67SKarl Rupp if (gs->num_pairs) PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim); /* tree only */ 1300*2fa5cd67SKarl Rupp else if (gs->max_left_over) PCTFS_gs_gop_tree_plus_hc(gs,vals,dim); 1301827bd09bSSatish Balay 1302ca8e9878SJed Brown PCTFS_gs_gop_local_out(gs,vals); 1303db4deed7SKarl Rupp } else { /* if intersection tree/pairwise and local is empty */ 1304827bd09bSSatish Balay /* pairwise will do tree inside */ 1305*2fa5cd67SKarl Rupp if (gs->num_pairs) PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim); /* tree */ 1306*2fa5cd67SKarl Rupp else if (gs->max_left_over) PCTFS_gs_gop_tree_plus_hc(gs,vals,dim); 1307827bd09bSSatish Balay } 13083fdc5746SBarry Smith PetscFunctionReturn(0); 1309827bd09bSSatish Balay } 1310827bd09bSSatish Balay 13117b1ae94cSBarry Smith /******************************************************************************/ 1312ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc(PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim) 1313827bd09bSSatish Balay { 1314a501084fSBarry Smith PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2; 131552f87cdaSBarry Smith PetscInt *iptr, *msg_list, *msg_size, **msg_nodes; 131652f87cdaSBarry Smith PetscInt *pw, *list, *size, **nodes; 1317827bd09bSSatish Balay MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out; 1318827bd09bSSatish Balay MPI_Status status; 131952f87cdaSBarry Smith PetscInt i, mask=1; 13203fdc5746SBarry Smith PetscErrorCode ierr; 1321827bd09bSSatish Balay 13223fdc5746SBarry Smith PetscFunctionBegin; 1323db4deed7SKarl Rupp for (i=1; i<dim; i++) { mask<<=1; mask++; } 1324827bd09bSSatish Balay 1325a501084fSBarry Smith /* strip and load s */ 1326827bd09bSSatish Balay msg_list = list = gs->pair_list; 1327827bd09bSSatish Balay msg_size = size = gs->msg_sizes; 1328827bd09bSSatish Balay msg_nodes = nodes = gs->node_list; 1329827bd09bSSatish Balay iptr = pw = gs->pw_elm_list; 1330827bd09bSSatish Balay dptr1 = dptr3 = gs->pw_vals; 1331827bd09bSSatish Balay msg_ids_in = ids_in = gs->msg_ids_in; 1332827bd09bSSatish Balay msg_ids_out = ids_out = gs->msg_ids_out; 1333827bd09bSSatish Balay dptr2 = gs->out; 1334827bd09bSSatish Balay in1 = in2 = gs->in; 1335827bd09bSSatish Balay 1336827bd09bSSatish Balay /* post the receives */ 1337827bd09bSSatish Balay /* msg_nodes=nodes; */ 1338db4deed7SKarl Rupp do { 1339827bd09bSSatish Balay /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the 1340827bd09bSSatish Balay second one *list and do list++ afterwards */ 1341db4deed7SKarl Rupp if ((PCTFS_my_id|mask)==(*list|mask)) { 1342ca8e9878SJed Brown ierr = MPI_Irecv(in1, *size, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr); 13439182e22cSBarry Smith list++; msg_ids_in++;in1 += *size++; 1344db4deed7SKarl Rupp } else { list++; size++; } 1345*2fa5cd67SKarl Rupp } while (*++msg_nodes); 1346827bd09bSSatish Balay 1347827bd09bSSatish Balay /* load gs values into in out gs buffers */ 1348*2fa5cd67SKarl Rupp while (*iptr >= 0) *dptr3++ = *(in_vals + *iptr++); 1349827bd09bSSatish Balay 1350827bd09bSSatish Balay /* load out buffers and post the sends */ 1351827bd09bSSatish Balay msg_nodes=nodes; 1352827bd09bSSatish Balay list = msg_list; 1353db4deed7SKarl Rupp while ((iptr = *msg_nodes++)) { 1354db4deed7SKarl Rupp if ((PCTFS_my_id|mask)==(*list|mask)) { 1355827bd09bSSatish Balay dptr3 = dptr2; 1356*2fa5cd67SKarl Rupp while (*iptr >= 0) *dptr2++ = *(dptr1 + *iptr++); 1357827bd09bSSatish Balay /* CHECK PERSISTENT COMMS MODE FOR ALL THIS STUFF */ 1358827bd09bSSatish Balay /* is msg_ids_out++ correct? */ 1359ca8e9878SJed Brown ierr = MPI_Isend(dptr3, *msg_size, MPIU_SCALAR, *list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr); 13609182e22cSBarry Smith msg_size++;list++;msg_ids_out++; 1361db4deed7SKarl Rupp } else {list++; msg_size++;} 1362827bd09bSSatish Balay } 1363827bd09bSSatish Balay 1364827bd09bSSatish Balay /* do the tree while we're waiting */ 1365*2fa5cd67SKarl Rupp if (gs->max_left_over) PCTFS_gs_gop_tree_plus_hc(gs,in_vals,dim); 1366827bd09bSSatish Balay 1367827bd09bSSatish Balay /* process the received data */ 1368827bd09bSSatish Balay msg_nodes=nodes; 1369827bd09bSSatish Balay list = msg_list; 1370db4deed7SKarl Rupp while ((iptr = *nodes++)) { 1371db4deed7SKarl Rupp if ((PCTFS_my_id|mask)==(*list|mask)) { 1372827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1373827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 13749182e22cSBarry Smith ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr); 13759182e22cSBarry Smith ids_in++; 1376*2fa5cd67SKarl Rupp while (*iptr >= 0) *(dptr1 + *iptr++) += *in2++; 1377827bd09bSSatish Balay } 1378827bd09bSSatish Balay list++; 1379827bd09bSSatish Balay } 1380827bd09bSSatish Balay 1381827bd09bSSatish Balay /* replace vals */ 1382*2fa5cd67SKarl Rupp while (*pw >= 0) *(in_vals + *pw++) = *dptr1++; 1383827bd09bSSatish Balay 1384827bd09bSSatish Balay /* clear isend message handles */ 1385827bd09bSSatish Balay /* This changed for clarity though it could be the same */ 1386db4deed7SKarl Rupp while (*msg_nodes++) { 1387db4deed7SKarl Rupp if ((PCTFS_my_id|mask)==(*msg_list|mask)) { 1388827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1389827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 13909182e22cSBarry Smith ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr); 13919182e22cSBarry Smith ids_out++; 1392827bd09bSSatish Balay } 1393827bd09bSSatish Balay msg_list++; 1394827bd09bSSatish Balay } 13953fdc5746SBarry Smith PetscFunctionReturn(0); 1396827bd09bSSatish Balay } 1397827bd09bSSatish Balay 13987b1ae94cSBarry Smith /******************************************************************************/ 1399ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim) 1400827bd09bSSatish Balay { 140152f87cdaSBarry Smith PetscInt size; 140252f87cdaSBarry Smith PetscInt *in, *out; 1403a501084fSBarry Smith PetscScalar *buf, *work; 140452f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 1405827bd09bSSatish Balay 14063fdc5746SBarry Smith PetscFunctionBegin; 1407827bd09bSSatish Balay in = gs->tree_map_in; 1408827bd09bSSatish Balay out = gs->tree_map_out; 1409827bd09bSSatish Balay buf = gs->tree_buf; 1410827bd09bSSatish Balay work = gs->tree_work; 1411827bd09bSSatish Balay size = gs->tree_nel; 1412827bd09bSSatish Balay 1413ca8e9878SJed Brown PCTFS_rvec_zero(buf,size); 1414827bd09bSSatish Balay 1415*2fa5cd67SKarl Rupp while (*in >= 0) *(buf + *out++) = *(vals + *in++); 1416827bd09bSSatish Balay 1417827bd09bSSatish Balay in = gs->tree_map_in; 1418827bd09bSSatish Balay out = gs->tree_map_out; 1419827bd09bSSatish Balay 1420b1c944f5SJed Brown PCTFS_grop_hc(buf,work,size,op,dim); 1421827bd09bSSatish Balay 1422*2fa5cd67SKarl Rupp while (*in >= 0) *(vals + *in++) = *(buf + *out++); 14233fdc5746SBarry Smith PetscFunctionReturn(0); 1424827bd09bSSatish Balay } 1425827bd09bSSatish Balay 1426827bd09bSSatish Balay 1427827bd09bSSatish Balay 1428