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); 211ca8e9878SJed Brown gs->PCTFS_gs_comm=PCTFS_gs_comm; 212827bd09bSSatish Balay 213827bd09bSSatish Balay return(gs); 214827bd09bSSatish Balay } 215827bd09bSSatish Balay 216f1ed62a8SBarry Smith /******************************************************************************/ 217ca8e9878SJed Brown static PCTFS_gs_id *gsi_new(void) 218827bd09bSSatish Balay { 219f1ed62a8SBarry Smith PetscErrorCode ierr; 220ca8e9878SJed Brown PCTFS_gs_id *gs; 221ca8e9878SJed Brown gs = (PCTFS_gs_id *) malloc(sizeof(PCTFS_gs_id)); 222ca8e9878SJed Brown ierr = PetscMemzero(gs,sizeof(PCTFS_gs_id));CHKERRABORT(PETSC_COMM_WORLD,ierr); 223827bd09bSSatish Balay return(gs); 224827bd09bSSatish Balay } 225827bd09bSSatish Balay 226f1ed62a8SBarry Smith /******************************************************************************/ 227ca8e9878SJed Brown static PCTFS_gs_id * gsi_check_args(PetscInt *in_elms, PetscInt nel, PetscInt level) 228827bd09bSSatish Balay { 22952f87cdaSBarry Smith PetscInt i, j, k, t2; 23052f87cdaSBarry Smith PetscInt *companion, *elms, *unique, *iptr; 23152f87cdaSBarry Smith PetscInt num_local=0, *num_to_reduce, **local_reduce; 23252f87cdaSBarry Smith PetscInt oprs[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_MIN,GL_B_AND}; 23352f87cdaSBarry Smith PetscInt vals[sizeof(oprs)/sizeof(oprs[0])-1]; 23452f87cdaSBarry Smith PetscInt work[sizeof(oprs)/sizeof(oprs[0])-1]; 235ca8e9878SJed Brown PCTFS_gs_id *gs; 236d1528f56SBarry Smith PetscErrorCode ierr; 237827bd09bSSatish Balay 238827bd09bSSatish Balay 239c1235816SBarry Smith if (!in_elms) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"elms point to nothing!!!\n"); 240c1235816SBarry Smith if (nel<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"can't have fewer than 0 elms!!!\n"); 241827bd09bSSatish Balay 242827bd09bSSatish Balay if (nel==0) 243f1ed62a8SBarry Smith {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 */ 253827bd09bSSatish Balay for (i=j=0;i<nel;i++) 254827bd09bSSatish Balay {if (in_elms[i]!=0) {j++;}} 255827bd09bSSatish Balay 256827bd09bSSatish Balay k=nel; nel=j; 257827bd09bSSatish Balay 258827bd09bSSatish Balay /* copy over in_elms list and create inverse map */ 25952f87cdaSBarry Smith elms = (PetscInt*) malloc((nel+1)*sizeof(PetscInt)); 26052f87cdaSBarry Smith companion = (PetscInt*) malloc(nel*sizeof(PetscInt)); 2611d7d0905SBarry Smith 262827bd09bSSatish Balay for (i=j=0;i<k;i++) 263827bd09bSSatish Balay { 264827bd09bSSatish Balay if (in_elms[i]!=0) 265827bd09bSSatish Balay {elms[j] = in_elms[i]; companion[j++] = i;} 266827bd09bSSatish Balay } 267827bd09bSSatish Balay 268c1235816SBarry Smith if (j!=nel) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"nel j mismatch!\n"); 269827bd09bSSatish Balay 270827bd09bSSatish Balay /* pre-pass ... check to see if sorted */ 271827bd09bSSatish Balay elms[nel] = INT_MAX; 272827bd09bSSatish Balay iptr = elms; 273827bd09bSSatish Balay unique = elms+1; 274827bd09bSSatish Balay j=0; 275827bd09bSSatish Balay while (*iptr!=INT_MAX) 276827bd09bSSatish Balay { 277827bd09bSSatish Balay if (*iptr++>*unique++) 278827bd09bSSatish Balay {j=1; break;} 279827bd09bSSatish Balay } 280827bd09bSSatish Balay 281827bd09bSSatish Balay /* set up inverse map */ 282827bd09bSSatish Balay if (j) 283827bd09bSSatish Balay { 284f1ed62a8SBarry Smith ierr = PetscInfo(0,"gsi_check_args() :: elm list *not* sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr); 285ca8e9878SJed Brown ierr = PCTFS_SMI_sort((void*)elms, (void*)companion, nel, SORT_INTEGER);CHKERRABORT(PETSC_COMM_WORLD,ierr); 286827bd09bSSatish Balay } 287827bd09bSSatish Balay else 288f1ed62a8SBarry Smith {ierr = PetscInfo(0,"gsi_check_args() :: elm list sorted!\n");CHKERRABORT(PETSC_COMM_WORLD,ierr);} 289827bd09bSSatish Balay elms[nel] = INT_MIN; 290827bd09bSSatish Balay 291827bd09bSSatish Balay /* first pass */ 292827bd09bSSatish Balay /* determine number of unique elements, check pd */ 293827bd09bSSatish Balay for (i=k=0;i<nel;i+=j) 294827bd09bSSatish Balay { 295827bd09bSSatish Balay t2 = elms[i]; 296827bd09bSSatish Balay j=++i; 297827bd09bSSatish Balay 298827bd09bSSatish Balay /* clump 'em for now */ 299827bd09bSSatish Balay while (elms[j]==t2) {j++;} 300827bd09bSSatish Balay 301827bd09bSSatish Balay /* how many together and num local */ 302827bd09bSSatish Balay if (j-=i) 303827bd09bSSatish Balay {num_local++; k+=j;} 304827bd09bSSatish Balay } 305827bd09bSSatish Balay 306827bd09bSSatish Balay /* how many unique elements? */ 307827bd09bSSatish Balay gs->repeats=k; 308827bd09bSSatish Balay gs->nel = nel-k; 309827bd09bSSatish Balay 310827bd09bSSatish Balay 311827bd09bSSatish Balay /* number of repeats? */ 312827bd09bSSatish Balay gs->num_local = num_local; 313827bd09bSSatish Balay num_local+=2; 31452f87cdaSBarry Smith gs->local_reduce=local_reduce=(PetscInt **)malloc(num_local*sizeof(PetscInt*)); 31552f87cdaSBarry Smith gs->num_local_reduce=num_to_reduce=(PetscInt*) malloc(num_local*sizeof(PetscInt)); 316827bd09bSSatish Balay 31752f87cdaSBarry Smith unique = (PetscInt*) malloc((gs->nel+1)*sizeof(PetscInt)); 318827bd09bSSatish Balay gs->elms = unique; 319827bd09bSSatish Balay gs->nel_total = nel; 320827bd09bSSatish Balay gs->local_elms = elms; 321827bd09bSSatish Balay gs->companion = companion; 322827bd09bSSatish Balay 323827bd09bSSatish Balay /* compess map as well as keep track of local ops */ 324827bd09bSSatish Balay for (num_local=i=j=0;i<gs->nel;i++) 325827bd09bSSatish Balay { 326827bd09bSSatish Balay k=j; 327827bd09bSSatish Balay t2 = unique[i] = elms[j]; 328827bd09bSSatish Balay companion[i] = companion[j]; 329827bd09bSSatish Balay 330827bd09bSSatish Balay while (elms[j]==t2) {j++;} 331827bd09bSSatish Balay 332827bd09bSSatish Balay if ((t2=(j-k))>1) 333827bd09bSSatish Balay { 334827bd09bSSatish Balay /* number together */ 335827bd09bSSatish Balay num_to_reduce[num_local] = t2++; 33652f87cdaSBarry Smith iptr = local_reduce[num_local++] = (PetscInt*)malloc(t2*sizeof(PetscInt)); 337827bd09bSSatish Balay 338827bd09bSSatish Balay /* to use binary searching don't remap until we check intersection */ 339827bd09bSSatish Balay *iptr++ = i; 340827bd09bSSatish Balay 341827bd09bSSatish Balay /* note that we're skipping the first one */ 342827bd09bSSatish Balay while (++k<j) 343827bd09bSSatish Balay {*(iptr++) = companion[k];} 344827bd09bSSatish Balay *iptr = -1; 345827bd09bSSatish Balay } 346827bd09bSSatish Balay } 347827bd09bSSatish Balay 348827bd09bSSatish Balay /* sentinel for ngh_buf */ 349827bd09bSSatish Balay unique[gs->nel]=INT_MAX; 350827bd09bSSatish Balay 351827bd09bSSatish Balay /* for two partition sort hack */ 352827bd09bSSatish Balay num_to_reduce[num_local] = 0; 353827bd09bSSatish Balay local_reduce[num_local] = NULL; 354827bd09bSSatish Balay num_to_reduce[++num_local] = 0; 355827bd09bSSatish Balay local_reduce[num_local] = NULL; 356827bd09bSSatish Balay 357827bd09bSSatish Balay /* load 'em up */ 358827bd09bSSatish Balay /* note one extra to hold NON_UNIFORM flag!!! */ 359827bd09bSSatish Balay vals[2] = vals[1] = vals[0] = nel; 360827bd09bSSatish Balay if (gs->nel>0) 361827bd09bSSatish Balay { 3621d7d0905SBarry Smith vals[3] = unique[0]; 3631d7d0905SBarry Smith vals[4] = unique[gs->nel-1]; 364827bd09bSSatish Balay } 365827bd09bSSatish Balay else 366827bd09bSSatish Balay { 3671d7d0905SBarry Smith vals[3] = INT_MAX; 3681d7d0905SBarry Smith vals[4] = INT_MIN; 369827bd09bSSatish Balay } 370827bd09bSSatish Balay vals[5] = level; 371827bd09bSSatish Balay vals[6] = num_gs_ids; 372827bd09bSSatish Balay 373827bd09bSSatish Balay /* GLOBAL: send 'em out */ 374b1c944f5SJed Brown ierr = PCTFS_giop(vals,work,sizeof(oprs)/sizeof(oprs[0])-1,oprs);CHKERRABORT(PETSC_COMM_WORLD,ierr); 375827bd09bSSatish Balay 376827bd09bSSatish Balay /* must be semi-pos def - only pairwise depends on this */ 377827bd09bSSatish Balay /* LATER - remove this restriction */ 378c1235816SBarry Smith if (vals[3]<0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system not semi-pos def \n"); 379c1235816SBarry Smith if (vals[4]==INT_MAX) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system ub too large !\n"); 380827bd09bSSatish Balay 381827bd09bSSatish Balay gs->nel_min = vals[0]; 382827bd09bSSatish Balay gs->nel_max = vals[1]; 383827bd09bSSatish Balay gs->nel_sum = vals[2]; 384827bd09bSSatish Balay gs->gl_min = vals[3]; 385827bd09bSSatish Balay gs->gl_max = vals[4]; 386827bd09bSSatish Balay gs->negl = vals[4]-vals[3]+1; 387827bd09bSSatish Balay 388c1235816SBarry Smith if (gs->negl<=0) SETERRABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"gsi_check_args() :: system empty or neg :: %d\n"); 389827bd09bSSatish Balay 390827bd09bSSatish Balay /* LATER :: add level == -1 -> program selects level */ 391827bd09bSSatish Balay if (vals[5]<0) 392827bd09bSSatish Balay {vals[5]=0;} 393b1c944f5SJed Brown else if (vals[5]>PCTFS_num_nodes) 394b1c944f5SJed Brown {vals[5]=PCTFS_num_nodes;} 395827bd09bSSatish Balay gs->level = vals[5]; 396827bd09bSSatish Balay 397827bd09bSSatish Balay return(gs); 398827bd09bSSatish Balay } 399827bd09bSSatish Balay 400f1ed62a8SBarry Smith /******************************************************************************/ 401ca8e9878SJed Brown static PetscErrorCode gsi_via_bit_mask(PCTFS_gs_id *gs) 402827bd09bSSatish Balay { 40352f87cdaSBarry Smith PetscInt i, nel, *elms; 40452f87cdaSBarry Smith PetscInt t1; 40552f87cdaSBarry Smith PetscInt **reduce; 40652f87cdaSBarry Smith PetscInt *map; 407f1ed62a8SBarry Smith PetscErrorCode ierr; 408827bd09bSSatish Balay 409f1ed62a8SBarry Smith PetscFunctionBegin; 410ca8e9878SJed Brown /* totally local removes ... PCTFS_ct_bits == 0 */ 411827bd09bSSatish Balay get_ngh_buf(gs); 412827bd09bSSatish Balay 413*94dd86cdSBarry Smith if (gs->level) set_pairwise(gs); 414*94dd86cdSBarry Smith if (gs->max_left_over) set_tree(gs); 415827bd09bSSatish Balay 416827bd09bSSatish Balay /* intersection local and pairwise/tree? */ 417827bd09bSSatish Balay gs->num_local_total = gs->num_local; 418827bd09bSSatish Balay gs->gop_local_reduce = gs->local_reduce; 419827bd09bSSatish Balay gs->num_gop_local_reduce = gs->num_local_reduce; 420827bd09bSSatish Balay 421827bd09bSSatish Balay map = gs->companion; 422827bd09bSSatish Balay 423827bd09bSSatish Balay /* is there any local compression */ 424d890fc11SSatish Balay if (!gs->num_local) { 425827bd09bSSatish Balay gs->local_strength = NONE; 426827bd09bSSatish Balay gs->num_local_gop = 0; 427d890fc11SSatish Balay } else { 428827bd09bSSatish Balay /* ok find intersection */ 429827bd09bSSatish Balay map = gs->companion; 430827bd09bSSatish Balay reduce = gs->local_reduce; 431827bd09bSSatish Balay for (i=0, t1=0; i<gs->num_local; i++, reduce++) 432827bd09bSSatish Balay { 433ca8e9878SJed Brown if ((PCTFS_ivec_binary_search(**reduce,gs->pw_elm_list,gs->len_pw_list)>=0) 434827bd09bSSatish Balay || 435ca8e9878SJed Brown PCTFS_ivec_binary_search(**reduce,gs->tree_map_in,gs->tree_map_sz)>=0) 436827bd09bSSatish Balay { 437827bd09bSSatish Balay t1++; 438e32f2f54SBarry Smith if (gs->num_local_reduce[i]<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nobody in list?"); 439827bd09bSSatish Balay gs->num_local_reduce[i] *= -1; 440827bd09bSSatish Balay } 441827bd09bSSatish Balay **reduce=map[**reduce]; 442827bd09bSSatish Balay } 443827bd09bSSatish Balay 444827bd09bSSatish Balay /* intersection is empty */ 445827bd09bSSatish Balay if (!t1) 446827bd09bSSatish Balay { 447827bd09bSSatish Balay gs->local_strength = FULL; 448827bd09bSSatish Balay gs->num_local_gop = 0; 449827bd09bSSatish Balay } 450827bd09bSSatish Balay /* intersection not empty */ 451827bd09bSSatish Balay else 452827bd09bSSatish Balay { 453827bd09bSSatish Balay gs->local_strength = PARTIAL; 454ca8e9878SJed Brown ierr = PCTFS_SMI_sort((void*)gs->num_local_reduce, (void*)gs->local_reduce, gs->num_local + 1, SORT_INT_PTR);CHKERRQ(ierr); 455827bd09bSSatish Balay 456827bd09bSSatish Balay gs->num_local_gop = t1; 457827bd09bSSatish Balay gs->num_local_total = gs->num_local; 458827bd09bSSatish Balay gs->num_local -= t1; 459827bd09bSSatish Balay gs->gop_local_reduce = gs->local_reduce; 460827bd09bSSatish Balay gs->num_gop_local_reduce = gs->num_local_reduce; 461827bd09bSSatish Balay 462827bd09bSSatish Balay for (i=0; i<t1; i++) 463827bd09bSSatish Balay { 464e32f2f54SBarry Smith if (gs->num_gop_local_reduce[i]>=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"they aren't negative?"); 465827bd09bSSatish Balay gs->num_gop_local_reduce[i] *= -1; 466827bd09bSSatish Balay gs->local_reduce++; 467827bd09bSSatish Balay gs->num_local_reduce++; 468827bd09bSSatish Balay } 469827bd09bSSatish Balay gs->local_reduce++; 470827bd09bSSatish Balay gs->num_local_reduce++; 471827bd09bSSatish Balay } 472827bd09bSSatish Balay } 473827bd09bSSatish Balay 474827bd09bSSatish Balay elms = gs->pw_elm_list; 475827bd09bSSatish Balay nel = gs->len_pw_list; 476827bd09bSSatish Balay for (i=0; i<nel; i++) 477827bd09bSSatish Balay {elms[i] = map[elms[i]];} 478827bd09bSSatish Balay 479827bd09bSSatish Balay elms = gs->tree_map_in; 480827bd09bSSatish Balay nel = gs->tree_map_sz; 481827bd09bSSatish Balay for (i=0; i<nel; i++) 482827bd09bSSatish Balay {elms[i] = map[elms[i]];} 483827bd09bSSatish Balay 484827bd09bSSatish Balay /* clean up */ 485a501084fSBarry Smith free((void*) gs->local_elms); 486a501084fSBarry Smith free((void*) gs->companion); 487a501084fSBarry Smith free((void*) gs->elms); 488a501084fSBarry Smith free((void*) gs->ngh_buf); 489827bd09bSSatish Balay gs->local_elms = gs->companion = gs->elms = gs->ngh_buf = NULL; 4903fdc5746SBarry Smith PetscFunctionReturn(0); 491827bd09bSSatish Balay } 492827bd09bSSatish Balay 493f1ed62a8SBarry Smith /******************************************************************************/ 49452f87cdaSBarry Smith static PetscErrorCode place_in_tree( PetscInt elm) 495827bd09bSSatish Balay { 49652f87cdaSBarry Smith PetscInt *tp, n; 497827bd09bSSatish Balay 4983fdc5746SBarry Smith PetscFunctionBegin; 499827bd09bSSatish Balay if (ntree==tree_buf_sz) 500827bd09bSSatish Balay { 501827bd09bSSatish Balay if (tree_buf_sz) 502827bd09bSSatish Balay { 503827bd09bSSatish Balay tp = tree_buf; 504827bd09bSSatish Balay n = tree_buf_sz; 505827bd09bSSatish Balay tree_buf_sz<<=1; 50652f87cdaSBarry Smith tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt)); 507ca8e9878SJed Brown PCTFS_ivec_copy(tree_buf,tp,n); 508a501084fSBarry Smith free(tp); 509827bd09bSSatish Balay } 510827bd09bSSatish Balay else 511827bd09bSSatish Balay { 512827bd09bSSatish Balay tree_buf_sz = TREE_BUF_SZ; 51352f87cdaSBarry Smith tree_buf = (PetscInt*)malloc(tree_buf_sz*sizeof(PetscInt)); 514827bd09bSSatish Balay } 515827bd09bSSatish Balay } 516827bd09bSSatish Balay 517827bd09bSSatish Balay tree_buf[ntree++] = elm; 5183fdc5746SBarry Smith PetscFunctionReturn(0); 519827bd09bSSatish Balay } 520827bd09bSSatish Balay 521f1ed62a8SBarry Smith /******************************************************************************/ 522ca8e9878SJed Brown static PetscErrorCode get_ngh_buf(PCTFS_gs_id *gs) 523827bd09bSSatish Balay { 52452f87cdaSBarry Smith PetscInt i, j, npw=0, ntree_map=0; 52552f87cdaSBarry Smith PetscInt p_mask_size, ngh_buf_size, buf_size; 52652f87cdaSBarry Smith PetscInt *p_mask, *sh_proc_mask, *pw_sh_proc_mask; 52752f87cdaSBarry Smith PetscInt *ngh_buf, *buf1, *buf2; 52852f87cdaSBarry Smith PetscInt offset, per_load, num_loads, or_ct, start, end; 52952f87cdaSBarry Smith PetscInt *ptr1, *ptr2, i_start, negl, nel, *elms; 53052f87cdaSBarry Smith PetscInt oper=GL_B_OR; 53152f87cdaSBarry Smith PetscInt *ptr3, *t_mask, level, ct1, ct2; 532f1ed62a8SBarry Smith PetscErrorCode ierr; 533827bd09bSSatish Balay 5343fdc5746SBarry Smith PetscFunctionBegin; 535827bd09bSSatish Balay /* to make life easier */ 536827bd09bSSatish Balay nel = gs->nel; 537827bd09bSSatish Balay elms = gs->elms; 538827bd09bSSatish Balay level = gs->level; 539827bd09bSSatish Balay 540b1c944f5SJed Brown /* det #bytes needed for processor bit masks and init w/mask cor. to PCTFS_my_id */ 541ca8e9878SJed Brown p_mask = (PetscInt*) malloc(p_mask_size=PCTFS_len_bit_mask(PCTFS_num_nodes)); 542ca8e9878SJed Brown ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr); 543827bd09bSSatish Balay 544827bd09bSSatish Balay /* allocate space for masks and info bufs */ 54552f87cdaSBarry Smith gs->nghs = sh_proc_mask = (PetscInt*) malloc(p_mask_size); 54652f87cdaSBarry Smith gs->pw_nghs = pw_sh_proc_mask = (PetscInt*) malloc(p_mask_size); 547827bd09bSSatish Balay gs->ngh_buf_sz = ngh_buf_size = p_mask_size*nel; 54852f87cdaSBarry Smith t_mask = (PetscInt*) malloc(p_mask_size); 54952f87cdaSBarry Smith gs->ngh_buf = ngh_buf = (PetscInt*) malloc(ngh_buf_size); 550827bd09bSSatish Balay 551827bd09bSSatish Balay /* comm buffer size ... memory usage bounded by ~2*msg_buf */ 552827bd09bSSatish Balay /* had thought I could exploit rendezvous threshold */ 553827bd09bSSatish Balay 554827bd09bSSatish Balay /* default is one pass */ 555827bd09bSSatish Balay per_load = negl = gs->negl; 556827bd09bSSatish Balay gs->num_loads = num_loads = 1; 557827bd09bSSatish Balay i=p_mask_size*negl; 558827bd09bSSatish Balay 559827bd09bSSatish Balay /* possible overflow on buffer size */ 560827bd09bSSatish Balay /* overflow hack */ 561827bd09bSSatish Balay if (i<0) {i=INT_MAX;} 562827bd09bSSatish Balay 56339945688SSatish Balay buf_size = PetscMin(msg_buf,i); 564827bd09bSSatish Balay 565827bd09bSSatish Balay /* can we do it? */ 566e32f2f54SBarry 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); 567827bd09bSSatish Balay 568b1c944f5SJed Brown /* get PCTFS_giop buf space ... make *only* one malloc */ 56952f87cdaSBarry Smith buf1 = (PetscInt*) malloc(buf_size<<1); 570827bd09bSSatish Balay 571827bd09bSSatish Balay /* more than one gior exchange needed? */ 572827bd09bSSatish Balay if (buf_size!=i) 573827bd09bSSatish Balay { 574827bd09bSSatish Balay per_load = buf_size/p_mask_size; 575827bd09bSSatish Balay buf_size = per_load*p_mask_size; 576827bd09bSSatish Balay gs->num_loads = num_loads = negl/per_load + (negl%per_load>0); 577827bd09bSSatish Balay } 578827bd09bSSatish Balay 579827bd09bSSatish Balay 580827bd09bSSatish Balay /* convert buf sizes from #bytes to #ints - 32 bit only! */ 581a501084fSBarry Smith p_mask_size/=sizeof(PetscInt); ngh_buf_size/=sizeof(PetscInt); buf_size/=sizeof(PetscInt); 582827bd09bSSatish Balay 583b1c944f5SJed Brown /* find PCTFS_giop work space */ 584827bd09bSSatish Balay buf2 = buf1+buf_size; 585827bd09bSSatish Balay 586827bd09bSSatish Balay /* hold #ints needed for processor masks */ 587827bd09bSSatish Balay gs->mask_sz=p_mask_size; 588827bd09bSSatish Balay 589827bd09bSSatish Balay /* init buffers */ 590ca8e9878SJed Brown ierr = PCTFS_ivec_zero(sh_proc_mask,p_mask_size);CHKERRQ(ierr); 591ca8e9878SJed Brown ierr = PCTFS_ivec_zero(pw_sh_proc_mask,p_mask_size);CHKERRQ(ierr); 592ca8e9878SJed Brown ierr = PCTFS_ivec_zero(ngh_buf,ngh_buf_size);CHKERRQ(ierr); 593827bd09bSSatish Balay 594827bd09bSSatish Balay /* HACK reset tree info */ 595827bd09bSSatish Balay tree_buf=NULL; 596827bd09bSSatish Balay tree_buf_sz=ntree=0; 597827bd09bSSatish Balay 598827bd09bSSatish Balay /* ok do it */ 599827bd09bSSatish Balay for (ptr1=ngh_buf,ptr2=elms,end=gs->gl_min,or_ct=i=0; or_ct<num_loads; or_ct++) 600827bd09bSSatish Balay { 601827bd09bSSatish Balay /* identity for bitwise or is 000...000 */ 602ca8e9878SJed Brown PCTFS_ivec_zero(buf1,buf_size); 603827bd09bSSatish Balay 604827bd09bSSatish Balay /* load msg buffer */ 605827bd09bSSatish Balay for (start=end,end+=per_load,i_start=i; (offset=*ptr2)<end; i++, ptr2++) 606827bd09bSSatish Balay { 607827bd09bSSatish Balay offset = (offset-start)*p_mask_size; 608ca8e9878SJed Brown PCTFS_ivec_copy(buf1+offset,p_mask,p_mask_size); 609827bd09bSSatish Balay } 610827bd09bSSatish Balay 611827bd09bSSatish Balay /* GLOBAL: pass buffer */ 612b1c944f5SJed Brown ierr = PCTFS_giop(buf1,buf2,buf_size,&oper);CHKERRQ(ierr); 613827bd09bSSatish Balay 614827bd09bSSatish Balay 615827bd09bSSatish Balay /* unload buffer into ngh_buf */ 616827bd09bSSatish Balay ptr2=(elms+i_start); 617827bd09bSSatish Balay for(ptr3=buf1,j=start; j<end; ptr3+=p_mask_size,j++) 618827bd09bSSatish Balay { 619827bd09bSSatish Balay /* I own it ... may have to pairwise it */ 620827bd09bSSatish Balay if (j==*ptr2) 621827bd09bSSatish Balay { 622827bd09bSSatish Balay /* do i share it w/anyone? */ 623ca8e9878SJed Brown ct1 = PCTFS_ct_bits((char *)ptr3,p_mask_size*sizeof(PetscInt)); 624827bd09bSSatish Balay /* guess not */ 625827bd09bSSatish Balay if (ct1<2) 626827bd09bSSatish Balay {ptr2++; ptr1+=p_mask_size; continue;} 627827bd09bSSatish Balay 628827bd09bSSatish Balay /* i do ... so keep info and turn off my bit */ 629ca8e9878SJed Brown PCTFS_ivec_copy(ptr1,ptr3,p_mask_size); 630ca8e9878SJed Brown ierr = PCTFS_ivec_xor(ptr1,p_mask,p_mask_size);CHKERRQ(ierr); 631ca8e9878SJed Brown ierr = PCTFS_ivec_or(sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr); 632827bd09bSSatish Balay 633827bd09bSSatish Balay /* is it to be done pairwise? */ 634827bd09bSSatish Balay if (--ct1<=level) 635827bd09bSSatish Balay { 636827bd09bSSatish Balay npw++; 637827bd09bSSatish Balay 638827bd09bSSatish Balay /* turn on high bit to indicate pw need to process */ 639827bd09bSSatish Balay *ptr2++ |= TOP_BIT; 640ca8e9878SJed Brown ierr = PCTFS_ivec_or(pw_sh_proc_mask,ptr1,p_mask_size);CHKERRQ(ierr); 641827bd09bSSatish Balay ptr1+=p_mask_size; 642827bd09bSSatish Balay continue; 643827bd09bSSatish Balay } 644827bd09bSSatish Balay 645827bd09bSSatish Balay /* get set for next and note that I have a tree contribution */ 646827bd09bSSatish Balay /* could save exact elm index for tree here -> save a search */ 647827bd09bSSatish Balay ptr2++; ptr1+=p_mask_size; ntree_map++; 648827bd09bSSatish Balay } 649827bd09bSSatish Balay /* i don't but still might be involved in tree */ 650827bd09bSSatish Balay else 651827bd09bSSatish Balay { 652827bd09bSSatish Balay 653827bd09bSSatish Balay /* shared by how many? */ 654ca8e9878SJed Brown ct1 = PCTFS_ct_bits((char *)ptr3,p_mask_size*sizeof(PetscInt)); 655827bd09bSSatish Balay 656827bd09bSSatish Balay /* none! */ 657f1ed62a8SBarry Smith if (ct1<2) continue; 658827bd09bSSatish Balay 659827bd09bSSatish Balay /* is it going to be done pairwise? but not by me of course!*/ 660f1ed62a8SBarry Smith if (--ct1<=level) continue; 661827bd09bSSatish Balay } 662827bd09bSSatish Balay /* LATER we're going to have to process it NOW */ 663827bd09bSSatish Balay /* nope ... tree it */ 664f1ed62a8SBarry Smith ierr = place_in_tree(j);CHKERRQ(ierr); 665827bd09bSSatish Balay } 666827bd09bSSatish Balay } 667827bd09bSSatish Balay 668a501084fSBarry Smith free((void*)t_mask); 669a501084fSBarry Smith free((void*)buf1); 670827bd09bSSatish Balay 671827bd09bSSatish Balay gs->len_pw_list=npw; 672ca8e9878SJed Brown gs->num_nghs = PCTFS_ct_bits((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt)); 673827bd09bSSatish Balay 674827bd09bSSatish Balay /* expand from bit mask list to int list and save ngh list */ 67552f87cdaSBarry Smith gs->nghs = (PetscInt*) malloc(gs->num_nghs * sizeof(PetscInt)); 676ca8e9878SJed Brown PCTFS_bm_to_proc((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt),gs->nghs); 677827bd09bSSatish Balay 678ca8e9878SJed Brown gs->num_pw_nghs = PCTFS_ct_bits((char *)pw_sh_proc_mask,p_mask_size*sizeof(PetscInt)); 679827bd09bSSatish Balay 680827bd09bSSatish Balay oper = GL_MAX; 681827bd09bSSatish Balay ct1 = gs->num_nghs; 682b1c944f5SJed Brown ierr = PCTFS_giop(&ct1,&ct2,1,&oper);CHKERRQ(ierr); 683827bd09bSSatish Balay gs->max_nghs = ct1; 684827bd09bSSatish Balay 685827bd09bSSatish Balay gs->tree_map_sz = ntree_map; 686827bd09bSSatish Balay gs->max_left_over=ntree; 687827bd09bSSatish Balay 688a501084fSBarry Smith free((void*)p_mask); 689a501084fSBarry Smith free((void*)sh_proc_mask); 6903fdc5746SBarry Smith PetscFunctionReturn(0); 691827bd09bSSatish Balay } 692827bd09bSSatish Balay 693f1ed62a8SBarry Smith /******************************************************************************/ 694ca8e9878SJed Brown static PetscErrorCode set_pairwise(PCTFS_gs_id *gs) 695827bd09bSSatish Balay { 69652f87cdaSBarry Smith PetscInt i, j; 69752f87cdaSBarry Smith PetscInt p_mask_size; 69852f87cdaSBarry Smith PetscInt *p_mask, *sh_proc_mask, *tmp_proc_mask; 69952f87cdaSBarry Smith PetscInt *ngh_buf, *buf2; 70052f87cdaSBarry Smith PetscInt offset; 70152f87cdaSBarry Smith PetscInt *msg_list, *msg_size, **msg_nodes, nprs; 70252f87cdaSBarry Smith PetscInt *pairwise_elm_list, len_pair_list=0; 70352f87cdaSBarry Smith PetscInt *iptr, t1, i_start, nel, *elms; 70452f87cdaSBarry Smith PetscInt ct; 705f1ed62a8SBarry Smith PetscErrorCode ierr; 706827bd09bSSatish Balay 7073fdc5746SBarry Smith PetscFunctionBegin; 708827bd09bSSatish Balay /* to make life easier */ 709827bd09bSSatish Balay nel = gs->nel; 710827bd09bSSatish Balay elms = gs->elms; 711827bd09bSSatish Balay ngh_buf = gs->ngh_buf; 712827bd09bSSatish Balay sh_proc_mask = gs->pw_nghs; 713827bd09bSSatish Balay 714827bd09bSSatish Balay /* need a few temp masks */ 715ca8e9878SJed Brown p_mask_size = PCTFS_len_bit_mask(PCTFS_num_nodes); 71652f87cdaSBarry Smith p_mask = (PetscInt*) malloc(p_mask_size); 71752f87cdaSBarry Smith tmp_proc_mask = (PetscInt*) malloc(p_mask_size); 718827bd09bSSatish Balay 719b1c944f5SJed Brown /* set mask to my PCTFS_my_id's bit mask */ 720ca8e9878SJed Brown ierr = PCTFS_set_bit_mask(p_mask,p_mask_size,PCTFS_my_id);CHKERRQ(ierr); 721827bd09bSSatish Balay 722a501084fSBarry Smith p_mask_size /= sizeof(PetscInt); 723827bd09bSSatish Balay 724827bd09bSSatish Balay len_pair_list=gs->len_pw_list; 72552f87cdaSBarry Smith gs->pw_elm_list=pairwise_elm_list=(PetscInt*)malloc((len_pair_list+1)*sizeof(PetscInt)); 726827bd09bSSatish Balay 727827bd09bSSatish Balay /* how many processors (nghs) do we have to exchange with? */ 728ca8e9878SJed Brown nprs=gs->num_pairs=PCTFS_ct_bits((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt)); 729827bd09bSSatish Balay 730827bd09bSSatish Balay 731ca8e9878SJed Brown /* allocate space for PCTFS_gs_gop() info */ 73252f87cdaSBarry Smith gs->pair_list = msg_list = (PetscInt *) malloc(sizeof(PetscInt)*nprs); 73352f87cdaSBarry Smith gs->msg_sizes = msg_size = (PetscInt *) malloc(sizeof(PetscInt)*nprs); 73452f87cdaSBarry Smith gs->node_list = msg_nodes = (PetscInt **) malloc(sizeof(PetscInt*)*(nprs+1)); 735827bd09bSSatish Balay 736827bd09bSSatish Balay /* init msg_size list */ 737ca8e9878SJed Brown ierr = PCTFS_ivec_zero(msg_size,nprs);CHKERRQ(ierr); 738827bd09bSSatish Balay 739827bd09bSSatish Balay /* expand from bit mask list to int list */ 740ca8e9878SJed Brown ierr = PCTFS_bm_to_proc((char *)sh_proc_mask,p_mask_size*sizeof(PetscInt),msg_list);CHKERRQ(ierr); 741827bd09bSSatish Balay 742827bd09bSSatish Balay /* keep list of elements being handled pairwise */ 743827bd09bSSatish Balay for (i=j=0;i<nel;i++) 744827bd09bSSatish Balay { 745827bd09bSSatish Balay if (elms[i] & TOP_BIT) 746827bd09bSSatish Balay {elms[i] ^= TOP_BIT; pairwise_elm_list[j++] = i;} 747827bd09bSSatish Balay } 748827bd09bSSatish Balay pairwise_elm_list[j] = -1; 749827bd09bSSatish Balay 750a501084fSBarry Smith gs->msg_ids_out = (MPI_Request *) malloc(sizeof(MPI_Request)*(nprs+1)); 751827bd09bSSatish Balay gs->msg_ids_out[nprs] = MPI_REQUEST_NULL; 752a501084fSBarry Smith gs->msg_ids_in = (MPI_Request *) malloc(sizeof(MPI_Request)*(nprs+1)); 753827bd09bSSatish Balay gs->msg_ids_in[nprs] = MPI_REQUEST_NULL; 754a501084fSBarry Smith gs->pw_vals = (PetscScalar *) malloc(sizeof(PetscScalar)*len_pair_list*vec_sz); 755827bd09bSSatish Balay 756827bd09bSSatish Balay /* find who goes to each processor */ 757827bd09bSSatish Balay for (i_start=i=0;i<nprs;i++) 758827bd09bSSatish Balay { 759827bd09bSSatish Balay /* processor i's mask */ 760ca8e9878SJed Brown ierr = PCTFS_set_bit_mask(p_mask,p_mask_size*sizeof(PetscInt),msg_list[i]);CHKERRQ(ierr); 761827bd09bSSatish Balay 762827bd09bSSatish Balay /* det # going to processor i */ 763827bd09bSSatish Balay for (ct=j=0;j<len_pair_list;j++) 764827bd09bSSatish Balay { 765827bd09bSSatish Balay buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size); 766ca8e9878SJed Brown ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr); 767ca8e9878SJed Brown if (PCTFS_ct_bits((char *)tmp_proc_mask,p_mask_size*sizeof(PetscInt))) 768827bd09bSSatish Balay {ct++;} 769827bd09bSSatish Balay } 770827bd09bSSatish Balay msg_size[i] = ct; 77139945688SSatish Balay i_start = PetscMax(i_start,ct); 772827bd09bSSatish Balay 773827bd09bSSatish Balay /*space to hold nodes in message to first neighbor */ 77452f87cdaSBarry Smith msg_nodes[i] = iptr = (PetscInt*) malloc(sizeof(PetscInt)*(ct+1)); 775827bd09bSSatish Balay 776827bd09bSSatish Balay for (j=0;j<len_pair_list;j++) 777827bd09bSSatish Balay { 778827bd09bSSatish Balay buf2 = ngh_buf+(pairwise_elm_list[j]*p_mask_size); 779ca8e9878SJed Brown ierr = PCTFS_ivec_and3(tmp_proc_mask,p_mask,buf2,p_mask_size);CHKERRQ(ierr); 780ca8e9878SJed Brown if (PCTFS_ct_bits((char *)tmp_proc_mask,p_mask_size*sizeof(PetscInt))) 781827bd09bSSatish Balay {*iptr++ = j;} 782827bd09bSSatish Balay } 783827bd09bSSatish Balay *iptr = -1; 784827bd09bSSatish Balay } 785827bd09bSSatish Balay msg_nodes[nprs] = NULL; 786827bd09bSSatish Balay 787827bd09bSSatish Balay j=gs->loc_node_pairs=i_start; 788827bd09bSSatish Balay t1 = GL_MAX; 789b1c944f5SJed Brown ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 790827bd09bSSatish Balay gs->max_node_pairs = i_start; 791827bd09bSSatish Balay 792827bd09bSSatish Balay i_start=j; 793827bd09bSSatish Balay t1 = GL_MIN; 794b1c944f5SJed Brown ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 795827bd09bSSatish Balay gs->min_node_pairs = i_start; 796827bd09bSSatish Balay 797827bd09bSSatish Balay i_start=j; 798827bd09bSSatish Balay t1 = GL_ADD; 799b1c944f5SJed Brown ierr = PCTFS_giop(&i_start,&offset,1,&t1);CHKERRQ(ierr); 800b1c944f5SJed Brown gs->avg_node_pairs = i_start/PCTFS_num_nodes + 1; 801827bd09bSSatish Balay 802827bd09bSSatish Balay i_start=nprs; 803827bd09bSSatish Balay t1 = GL_MAX; 804b1c944f5SJed Brown PCTFS_giop(&i_start,&offset,1,&t1); 805827bd09bSSatish Balay gs->max_pairs = i_start; 806827bd09bSSatish Balay 807827bd09bSSatish Balay 808827bd09bSSatish Balay /* remap pairwise in tail of gsi_via_bit_mask() */ 809ca8e9878SJed Brown gs->msg_total = PCTFS_ivec_sum(gs->msg_sizes,nprs); 810a501084fSBarry Smith gs->out = (PetscScalar *) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz); 811a501084fSBarry Smith gs->in = (PetscScalar *) malloc(sizeof(PetscScalar)*gs->msg_total*vec_sz); 812827bd09bSSatish Balay 813827bd09bSSatish Balay /* reset malloc pool */ 814a501084fSBarry Smith free((void*)p_mask); 815a501084fSBarry Smith free((void*)tmp_proc_mask); 8163fdc5746SBarry Smith PetscFunctionReturn(0); 817827bd09bSSatish Balay } 818827bd09bSSatish Balay 819f1ed62a8SBarry Smith /* to do pruned tree just save ngh buf copy for each one and decode here! 820827bd09bSSatish Balay ******************************************************************************/ 821ca8e9878SJed Brown static PetscErrorCode set_tree(PCTFS_gs_id *gs) 822827bd09bSSatish Balay { 82352f87cdaSBarry Smith PetscInt i, j, n, nel; 82452f87cdaSBarry Smith PetscInt *iptr_in, *iptr_out, *tree_elms, *elms; 825827bd09bSSatish Balay 8263fdc5746SBarry Smith PetscFunctionBegin; 827827bd09bSSatish Balay /* local work ptrs */ 828827bd09bSSatish Balay elms = gs->elms; 829827bd09bSSatish Balay nel = gs->nel; 830827bd09bSSatish Balay 831827bd09bSSatish Balay /* how many via tree */ 832827bd09bSSatish Balay gs->tree_nel = n = ntree; 833827bd09bSSatish Balay gs->tree_elms = tree_elms = iptr_in = tree_buf; 834a501084fSBarry Smith gs->tree_buf = (PetscScalar *) malloc(sizeof(PetscScalar)*n*vec_sz); 835a501084fSBarry Smith gs->tree_work = (PetscScalar *) malloc(sizeof(PetscScalar)*n*vec_sz); 836827bd09bSSatish Balay j=gs->tree_map_sz; 83752f87cdaSBarry Smith gs->tree_map_in = iptr_in = (PetscInt*) malloc(sizeof(PetscInt)*(j+1)); 83852f87cdaSBarry Smith gs->tree_map_out = iptr_out = (PetscInt*) malloc(sizeof(PetscInt)*(j+1)); 839827bd09bSSatish Balay 840827bd09bSSatish Balay /* search the longer of the two lists */ 841827bd09bSSatish Balay /* note ... could save this info in get_ngh_buf and save searches */ 842827bd09bSSatish Balay if (n<=nel) 843827bd09bSSatish Balay { 844827bd09bSSatish Balay /* bijective fct w/remap - search elm list */ 845827bd09bSSatish Balay for (i=0; i<n; i++) 846827bd09bSSatish Balay { 847ca8e9878SJed Brown if ((j=PCTFS_ivec_binary_search(*tree_elms++,elms,nel))>=0) 848827bd09bSSatish Balay {*iptr_in++ = j; *iptr_out++ = i;} 849827bd09bSSatish Balay } 850827bd09bSSatish Balay } 851827bd09bSSatish Balay else 852827bd09bSSatish Balay { 853827bd09bSSatish Balay for (i=0; i<nel; i++) 854827bd09bSSatish Balay { 855ca8e9878SJed Brown if ((j=PCTFS_ivec_binary_search(*elms++,tree_elms,n))>=0) 856827bd09bSSatish Balay {*iptr_in++ = i; *iptr_out++ = j;} 857827bd09bSSatish Balay } 858827bd09bSSatish Balay } 859827bd09bSSatish Balay 860827bd09bSSatish Balay /* sentinel */ 861827bd09bSSatish Balay *iptr_in = *iptr_out = -1; 8623fdc5746SBarry Smith PetscFunctionReturn(0); 863827bd09bSSatish Balay } 864827bd09bSSatish Balay 865f1ed62a8SBarry Smith /******************************************************************************/ 866ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_out( PCTFS_gs_id *gs, PetscScalar *vals) 867827bd09bSSatish Balay { 86852f87cdaSBarry Smith PetscInt *num, *map, **reduce; 869a501084fSBarry Smith PetscScalar tmp; 870827bd09bSSatish Balay 8713fdc5746SBarry Smith PetscFunctionBegin; 872827bd09bSSatish Balay num = gs->num_gop_local_reduce; 873827bd09bSSatish Balay reduce = gs->gop_local_reduce; 874827bd09bSSatish Balay while ((map = *reduce++)) 875827bd09bSSatish Balay { 876827bd09bSSatish Balay /* wall */ 877827bd09bSSatish Balay if (*num == 2) 878827bd09bSSatish Balay { 879827bd09bSSatish Balay num ++; 880827bd09bSSatish Balay vals[map[1]] = vals[map[0]]; 881827bd09bSSatish Balay } 882827bd09bSSatish Balay /* corner shared by three elements */ 883827bd09bSSatish Balay else if (*num == 3) 884827bd09bSSatish Balay { 885827bd09bSSatish Balay num ++; 886827bd09bSSatish Balay vals[map[2]] = vals[map[1]] = vals[map[0]]; 887827bd09bSSatish Balay } 888827bd09bSSatish Balay /* corner shared by four elements */ 889827bd09bSSatish Balay else if (*num == 4) 890827bd09bSSatish Balay { 891827bd09bSSatish Balay num ++; 892827bd09bSSatish Balay vals[map[3]] = vals[map[2]] = vals[map[1]] = vals[map[0]]; 893827bd09bSSatish Balay } 894827bd09bSSatish Balay /* general case ... odd geoms ... 3D*/ 895827bd09bSSatish Balay else 896827bd09bSSatish Balay { 897827bd09bSSatish Balay num++; 898827bd09bSSatish Balay tmp = *(vals + *map++); 899827bd09bSSatish Balay while (*map >= 0) 900827bd09bSSatish Balay {*(vals + *map++) = tmp;} 901827bd09bSSatish Balay } 902827bd09bSSatish Balay } 9033fdc5746SBarry Smith PetscFunctionReturn(0); 904827bd09bSSatish Balay } 905827bd09bSSatish Balay 9067b1ae94cSBarry Smith /******************************************************************************/ 907ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_plus( PCTFS_gs_id *gs, PetscScalar *vals) 908827bd09bSSatish Balay { 90952f87cdaSBarry Smith PetscInt *num, *map, **reduce; 910a501084fSBarry Smith PetscScalar tmp; 911827bd09bSSatish Balay 9123fdc5746SBarry Smith PetscFunctionBegin; 913827bd09bSSatish Balay num = gs->num_local_reduce; 914827bd09bSSatish Balay reduce = gs->local_reduce; 915827bd09bSSatish Balay while ((map = *reduce)) 916827bd09bSSatish Balay { 917827bd09bSSatish Balay /* wall */ 918827bd09bSSatish Balay if (*num == 2) 919827bd09bSSatish Balay { 920827bd09bSSatish Balay num ++; reduce++; 921827bd09bSSatish Balay vals[map[1]] = vals[map[0]] += vals[map[1]]; 922827bd09bSSatish Balay } 923827bd09bSSatish Balay /* corner shared by three elements */ 924827bd09bSSatish Balay else if (*num == 3) 925827bd09bSSatish Balay { 926827bd09bSSatish Balay num ++; reduce++; 927827bd09bSSatish Balay vals[map[2]]=vals[map[1]]=vals[map[0]]+=(vals[map[1]]+vals[map[2]]); 928827bd09bSSatish Balay } 929827bd09bSSatish Balay /* corner shared by four elements */ 930827bd09bSSatish Balay else if (*num == 4) 931827bd09bSSatish Balay { 932827bd09bSSatish Balay num ++; reduce++; 933827bd09bSSatish Balay vals[map[1]]=vals[map[2]]=vals[map[3]]=vals[map[0]] += 934827bd09bSSatish Balay (vals[map[1]] + vals[map[2]] + vals[map[3]]); 935827bd09bSSatish Balay } 936827bd09bSSatish Balay /* general case ... odd geoms ... 3D*/ 937827bd09bSSatish Balay else 938827bd09bSSatish Balay { 939827bd09bSSatish Balay num ++; 940827bd09bSSatish Balay tmp = 0.0; 941827bd09bSSatish Balay while (*map >= 0) 942827bd09bSSatish Balay {tmp += *(vals + *map++);} 943827bd09bSSatish Balay 944827bd09bSSatish Balay map = *reduce++; 945827bd09bSSatish Balay while (*map >= 0) 946827bd09bSSatish Balay {*(vals + *map++) = tmp;} 947827bd09bSSatish Balay } 948827bd09bSSatish Balay } 9493fdc5746SBarry Smith PetscFunctionReturn(0); 950827bd09bSSatish Balay } 951827bd09bSSatish Balay 9527b1ae94cSBarry Smith /******************************************************************************/ 953ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_local_in_plus( PCTFS_gs_id *gs, PetscScalar *vals) 954827bd09bSSatish Balay { 95552f87cdaSBarry Smith PetscInt *num, *map, **reduce; 956a501084fSBarry Smith PetscScalar *base; 957827bd09bSSatish Balay 9583fdc5746SBarry Smith PetscFunctionBegin; 959827bd09bSSatish Balay num = gs->num_gop_local_reduce; 960827bd09bSSatish Balay reduce = gs->gop_local_reduce; 961827bd09bSSatish Balay while ((map = *reduce++)) 962827bd09bSSatish Balay { 963827bd09bSSatish Balay /* wall */ 964827bd09bSSatish Balay if (*num == 2) 965827bd09bSSatish Balay { 966827bd09bSSatish Balay num ++; 967827bd09bSSatish Balay vals[map[0]] += vals[map[1]]; 968827bd09bSSatish Balay } 969827bd09bSSatish Balay /* corner shared by three elements */ 970827bd09bSSatish Balay else if (*num == 3) 971827bd09bSSatish Balay { 972827bd09bSSatish Balay num ++; 973827bd09bSSatish Balay vals[map[0]] += (vals[map[1]] + vals[map[2]]); 974827bd09bSSatish Balay } 975827bd09bSSatish Balay /* corner shared by four elements */ 976827bd09bSSatish Balay else if (*num == 4) 977827bd09bSSatish Balay { 978827bd09bSSatish Balay num ++; 979827bd09bSSatish Balay vals[map[0]] += (vals[map[1]] + vals[map[2]] + vals[map[3]]); 980827bd09bSSatish Balay } 981827bd09bSSatish Balay /* general case ... odd geoms ... 3D*/ 982827bd09bSSatish Balay else 983827bd09bSSatish Balay { 984827bd09bSSatish Balay num++; 985827bd09bSSatish Balay base = vals + *map++; 986827bd09bSSatish Balay while (*map >= 0) 987827bd09bSSatish Balay {*base += *(vals + *map++);} 988827bd09bSSatish Balay } 989827bd09bSSatish Balay } 9903fdc5746SBarry Smith PetscFunctionReturn(0); 991827bd09bSSatish Balay } 992827bd09bSSatish Balay 9937b1ae94cSBarry Smith /******************************************************************************/ 994ca8e9878SJed Brown PetscErrorCode PCTFS_gs_free( PCTFS_gs_id *gs) 995827bd09bSSatish Balay { 99652f87cdaSBarry Smith PetscInt i; 997827bd09bSSatish Balay 9983fdc5746SBarry Smith PetscFunctionBegin; 999a501084fSBarry Smith if (gs->nghs) {free((void*) gs->nghs);} 1000a501084fSBarry Smith if (gs->pw_nghs) {free((void*) gs->pw_nghs);} 1001827bd09bSSatish Balay 1002827bd09bSSatish Balay /* tree */ 1003827bd09bSSatish Balay if (gs->max_left_over) 1004827bd09bSSatish Balay { 1005a501084fSBarry Smith if (gs->tree_elms) {free((void*) gs->tree_elms);} 1006a501084fSBarry Smith if (gs->tree_buf) {free((void*) gs->tree_buf);} 1007a501084fSBarry Smith if (gs->tree_work) {free((void*) gs->tree_work);} 1008a501084fSBarry Smith if (gs->tree_map_in) {free((void*) gs->tree_map_in);} 1009a501084fSBarry Smith if (gs->tree_map_out) {free((void*) gs->tree_map_out);} 1010827bd09bSSatish Balay } 1011827bd09bSSatish Balay 1012827bd09bSSatish Balay /* pairwise info */ 1013827bd09bSSatish Balay if (gs->num_pairs) 1014827bd09bSSatish Balay { 1015827bd09bSSatish Balay /* should be NULL already */ 1016a501084fSBarry Smith if (gs->ngh_buf) {free((void*) gs->ngh_buf);} 1017a501084fSBarry Smith if (gs->elms) {free((void*) gs->elms);} 1018a501084fSBarry Smith if (gs->local_elms) {free((void*) gs->local_elms);} 1019a501084fSBarry Smith if (gs->companion) {free((void*) gs->companion);} 1020827bd09bSSatish Balay 1021827bd09bSSatish Balay /* only set if pairwise */ 1022a501084fSBarry Smith if (gs->vals) {free((void*) gs->vals);} 1023a501084fSBarry Smith if (gs->in) {free((void*) gs->in);} 1024a501084fSBarry Smith if (gs->out) {free((void*) gs->out);} 1025a501084fSBarry Smith if (gs->msg_ids_in) {free((void*) gs->msg_ids_in);} 1026a501084fSBarry Smith if (gs->msg_ids_out) {free((void*) gs->msg_ids_out);} 1027a501084fSBarry Smith if (gs->pw_vals) {free((void*) gs->pw_vals);} 1028a501084fSBarry Smith if (gs->pw_elm_list) {free((void*) gs->pw_elm_list);} 1029827bd09bSSatish Balay if (gs->node_list) 1030827bd09bSSatish Balay { 1031827bd09bSSatish Balay for (i=0;i<gs->num_pairs;i++) 1032a501084fSBarry Smith {if (gs->node_list[i]) {free((void*) gs->node_list[i]);}} 1033a501084fSBarry Smith free((void*) gs->node_list); 1034827bd09bSSatish Balay } 1035a501084fSBarry Smith if (gs->msg_sizes) {free((void*) gs->msg_sizes);} 1036a501084fSBarry Smith if (gs->pair_list) {free((void*) gs->pair_list);} 1037827bd09bSSatish Balay } 1038827bd09bSSatish Balay 1039827bd09bSSatish Balay /* local info */ 1040827bd09bSSatish Balay if (gs->num_local_total>=0) 1041827bd09bSSatish Balay { 1042827bd09bSSatish Balay for (i=0;i<gs->num_local_total+1;i++) 1043827bd09bSSatish Balay /* for (i=0;i<gs->num_local_total;i++) */ 1044827bd09bSSatish Balay { 1045827bd09bSSatish Balay if (gs->num_gop_local_reduce[i]) 1046a501084fSBarry Smith {free((void*) gs->gop_local_reduce[i]);} 1047827bd09bSSatish Balay } 1048827bd09bSSatish Balay } 1049827bd09bSSatish Balay 1050827bd09bSSatish Balay /* if intersection tree/pairwise and local isn't empty */ 1051a501084fSBarry Smith if (gs->gop_local_reduce) {free((void*) gs->gop_local_reduce);} 1052a501084fSBarry Smith if (gs->num_gop_local_reduce) {free((void*) gs->num_gop_local_reduce);} 1053827bd09bSSatish Balay 1054a501084fSBarry Smith free((void*) gs); 10553fdc5746SBarry Smith PetscFunctionReturn(0); 1056827bd09bSSatish Balay } 1057827bd09bSSatish Balay 10587b1ae94cSBarry Smith /******************************************************************************/ 1059ca8e9878SJed Brown PetscErrorCode PCTFS_gs_gop_vec( PCTFS_gs_id *gs, PetscScalar *vals, const char *op, PetscInt step) 1060827bd09bSSatish Balay { 1061d1528f56SBarry Smith PetscErrorCode ierr; 1062d1528f56SBarry Smith 10633fdc5746SBarry Smith PetscFunctionBegin; 1064827bd09bSSatish Balay switch (*op) { 1065827bd09bSSatish Balay case '+': 1066ca8e9878SJed Brown PCTFS_gs_gop_vec_plus(gs,vals,step); 1067827bd09bSSatish Balay break; 1068827bd09bSSatish Balay default: 1069ca8e9878SJed Brown ierr = PetscInfo1(0,"PCTFS_gs_gop_vec() :: %c is not a valid op",op[0]);CHKERRQ(ierr); 1070ca8e9878SJed Brown ierr = PetscInfo(0,"PCTFS_gs_gop_vec() :: default :: plus");CHKERRQ(ierr); 1071ca8e9878SJed Brown PCTFS_gs_gop_vec_plus(gs,vals,step); 1072827bd09bSSatish Balay break; 1073827bd09bSSatish Balay } 10743fdc5746SBarry Smith PetscFunctionReturn(0); 1075827bd09bSSatish Balay } 1076827bd09bSSatish Balay 10777b1ae94cSBarry Smith /******************************************************************************/ 1078ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_plus( PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1079827bd09bSSatish Balay { 10803fdc5746SBarry Smith PetscFunctionBegin; 1081ca8e9878SJed Brown if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_gs_gop_vec() passed NULL gs handle!!!"); 1082827bd09bSSatish Balay 1083827bd09bSSatish Balay /* local only operations!!! */ 1084827bd09bSSatish Balay if (gs->num_local) 1085ca8e9878SJed Brown {PCTFS_gs_gop_vec_local_plus(gs,vals,step);} 1086827bd09bSSatish Balay 1087827bd09bSSatish Balay /* if intersection tree/pairwise and local isn't empty */ 1088827bd09bSSatish Balay if (gs->num_local_gop) 1089827bd09bSSatish Balay { 1090ca8e9878SJed Brown PCTFS_gs_gop_vec_local_in_plus(gs,vals,step); 1091827bd09bSSatish Balay 1092827bd09bSSatish Balay /* pairwise */ 1093827bd09bSSatish Balay if (gs->num_pairs) 1094ca8e9878SJed Brown {PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step);} 1095827bd09bSSatish Balay 1096827bd09bSSatish Balay /* tree */ 1097827bd09bSSatish Balay else if (gs->max_left_over) 1098ca8e9878SJed Brown {PCTFS_gs_gop_vec_tree_plus(gs,vals,step);} 1099827bd09bSSatish Balay 1100ca8e9878SJed Brown PCTFS_gs_gop_vec_local_out(gs,vals,step); 1101827bd09bSSatish Balay } 1102827bd09bSSatish Balay /* if intersection tree/pairwise and local is empty */ 1103827bd09bSSatish Balay else 1104827bd09bSSatish Balay { 1105827bd09bSSatish Balay /* pairwise */ 1106827bd09bSSatish Balay if (gs->num_pairs) 1107ca8e9878SJed Brown {PCTFS_gs_gop_vec_pairwise_plus(gs,vals,step);} 1108827bd09bSSatish Balay 1109827bd09bSSatish Balay /* tree */ 1110827bd09bSSatish Balay else if (gs->max_left_over) 1111ca8e9878SJed Brown {PCTFS_gs_gop_vec_tree_plus(gs,vals,step);} 1112827bd09bSSatish Balay } 11133fdc5746SBarry Smith PetscFunctionReturn(0); 1114827bd09bSSatish Balay } 1115827bd09bSSatish Balay 11167b1ae94cSBarry Smith /******************************************************************************/ 1117ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_plus( PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1118827bd09bSSatish Balay { 111952f87cdaSBarry Smith PetscInt *num, *map, **reduce; 1120a501084fSBarry Smith PetscScalar *base; 1121827bd09bSSatish Balay 11223fdc5746SBarry Smith PetscFunctionBegin; 1123827bd09bSSatish Balay num = gs->num_local_reduce; 1124827bd09bSSatish Balay reduce = gs->local_reduce; 1125827bd09bSSatish Balay while ((map = *reduce)) 1126827bd09bSSatish Balay { 1127827bd09bSSatish Balay base = vals + map[0] * step; 1128827bd09bSSatish Balay 1129827bd09bSSatish Balay /* wall */ 1130827bd09bSSatish Balay if (*num == 2) 1131827bd09bSSatish Balay { 1132827bd09bSSatish Balay num++; reduce++; 1133ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[1]*step,step); 1134ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1135827bd09bSSatish Balay } 1136827bd09bSSatish Balay /* corner shared by three elements */ 1137827bd09bSSatish Balay else if (*num == 3) 1138827bd09bSSatish Balay { 1139827bd09bSSatish Balay num++; reduce++; 1140ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[1]*step,step); 1141ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[2]*step,step); 1142ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1143ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1144827bd09bSSatish Balay } 1145827bd09bSSatish Balay /* corner shared by four elements */ 1146827bd09bSSatish Balay else if (*num == 4) 1147827bd09bSSatish Balay { 1148827bd09bSSatish Balay num++; reduce++; 1149ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[1]*step,step); 1150ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[2]*step,step); 1151ca8e9878SJed Brown PCTFS_rvec_add (base,vals+map[3]*step,step); 1152ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[3]*step,base,step); 1153ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1154ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1155827bd09bSSatish Balay } 1156827bd09bSSatish Balay /* general case ... odd geoms ... 3D */ 1157827bd09bSSatish Balay else 1158827bd09bSSatish Balay { 1159827bd09bSSatish Balay num++; 1160827bd09bSSatish Balay while (*++map >= 0) 1161ca8e9878SJed Brown {PCTFS_rvec_add (base,vals+*map*step,step);} 1162827bd09bSSatish Balay 1163827bd09bSSatish Balay map = *reduce; 1164827bd09bSSatish Balay while (*++map >= 0) 1165ca8e9878SJed Brown {PCTFS_rvec_copy(vals+*map*step,base,step);} 1166827bd09bSSatish Balay 1167827bd09bSSatish Balay reduce++; 1168827bd09bSSatish Balay } 1169827bd09bSSatish Balay } 11703fdc5746SBarry Smith PetscFunctionReturn(0); 1171827bd09bSSatish Balay } 1172827bd09bSSatish Balay 11737b1ae94cSBarry Smith /******************************************************************************/ 1174ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_in_plus( PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1175827bd09bSSatish Balay { 117652f87cdaSBarry Smith PetscInt *num, *map, **reduce; 1177a501084fSBarry Smith PetscScalar *base; 11783fdc5746SBarry Smith PetscFunctionBegin; 1179827bd09bSSatish Balay num = gs->num_gop_local_reduce; 1180827bd09bSSatish Balay reduce = gs->gop_local_reduce; 1181827bd09bSSatish Balay while ((map = *reduce++)) 1182827bd09bSSatish Balay { 1183827bd09bSSatish Balay base = vals + map[0] * step; 1184827bd09bSSatish Balay 1185827bd09bSSatish Balay /* wall */ 1186827bd09bSSatish Balay if (*num == 2) 1187827bd09bSSatish Balay { 1188827bd09bSSatish Balay num ++; 1189ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[1]*step,step); 1190827bd09bSSatish Balay } 1191827bd09bSSatish Balay /* corner shared by three elements */ 1192827bd09bSSatish Balay else if (*num == 3) 1193827bd09bSSatish Balay { 1194827bd09bSSatish Balay num ++; 1195ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[1]*step,step); 1196ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[2]*step,step); 1197827bd09bSSatish Balay } 1198827bd09bSSatish Balay /* corner shared by four elements */ 1199827bd09bSSatish Balay else if (*num == 4) 1200827bd09bSSatish Balay { 1201827bd09bSSatish Balay num ++; 1202ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[1]*step,step); 1203ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[2]*step,step); 1204ca8e9878SJed Brown PCTFS_rvec_add(base,vals+map[3]*step,step); 1205827bd09bSSatish Balay } 1206827bd09bSSatish Balay /* general case ... odd geoms ... 3D*/ 1207827bd09bSSatish Balay else 1208827bd09bSSatish Balay { 1209827bd09bSSatish Balay num++; 1210827bd09bSSatish Balay while (*++map >= 0) 1211ca8e9878SJed Brown {PCTFS_rvec_add(base,vals+*map*step,step);} 1212827bd09bSSatish Balay } 1213827bd09bSSatish Balay } 12143fdc5746SBarry Smith PetscFunctionReturn(0); 1215827bd09bSSatish Balay } 1216827bd09bSSatish Balay 12177b1ae94cSBarry Smith /******************************************************************************/ 1218ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_local_out( PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1219827bd09bSSatish Balay { 122052f87cdaSBarry Smith PetscInt *num, *map, **reduce; 1221a501084fSBarry Smith PetscScalar *base; 1222827bd09bSSatish Balay 12233fdc5746SBarry Smith PetscFunctionBegin; 1224827bd09bSSatish Balay num = gs->num_gop_local_reduce; 1225827bd09bSSatish Balay reduce = gs->gop_local_reduce; 1226827bd09bSSatish Balay while ((map = *reduce++)) 1227827bd09bSSatish Balay { 1228827bd09bSSatish Balay base = vals + map[0] * step; 1229827bd09bSSatish Balay 1230827bd09bSSatish Balay /* wall */ 1231827bd09bSSatish Balay if (*num == 2) 1232827bd09bSSatish Balay { 1233827bd09bSSatish Balay num ++; 1234ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1235827bd09bSSatish Balay } 1236827bd09bSSatish Balay /* corner shared by three elements */ 1237827bd09bSSatish Balay else if (*num == 3) 1238827bd09bSSatish Balay { 1239827bd09bSSatish Balay num ++; 1240ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1241ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1242827bd09bSSatish Balay } 1243827bd09bSSatish Balay /* corner shared by four elements */ 1244827bd09bSSatish Balay else if (*num == 4) 1245827bd09bSSatish Balay { 1246827bd09bSSatish Balay num ++; 1247ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[1]*step,base,step); 1248ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[2]*step,base,step); 1249ca8e9878SJed Brown PCTFS_rvec_copy(vals+map[3]*step,base,step); 1250827bd09bSSatish Balay } 1251827bd09bSSatish Balay /* general case ... odd geoms ... 3D*/ 1252827bd09bSSatish Balay else 1253827bd09bSSatish Balay { 1254827bd09bSSatish Balay num++; 1255827bd09bSSatish Balay while (*++map >= 0) 1256ca8e9878SJed Brown {PCTFS_rvec_copy(vals+*map*step,base,step);} 1257827bd09bSSatish Balay } 1258827bd09bSSatish Balay } 12593fdc5746SBarry Smith PetscFunctionReturn(0); 1260827bd09bSSatish Balay } 1261827bd09bSSatish Balay 12627b1ae94cSBarry Smith /******************************************************************************/ 1263ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_pairwise_plus( PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt step) 1264827bd09bSSatish Balay { 1265a501084fSBarry Smith PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2; 126652f87cdaSBarry Smith PetscInt *iptr, *msg_list, *msg_size, **msg_nodes; 126752f87cdaSBarry Smith PetscInt *pw, *list, *size, **nodes; 1268827bd09bSSatish Balay MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out; 1269827bd09bSSatish Balay MPI_Status status; 12700805154bSBarry Smith PetscBLASInt i1 = 1,dstep; 12713fdc5746SBarry Smith PetscErrorCode ierr; 1272827bd09bSSatish Balay 12733fdc5746SBarry Smith PetscFunctionBegin; 1274a501084fSBarry Smith /* strip and load s */ 1275827bd09bSSatish Balay msg_list =list = gs->pair_list; 1276827bd09bSSatish Balay msg_size =size = gs->msg_sizes; 1277827bd09bSSatish Balay msg_nodes=nodes = gs->node_list; 1278827bd09bSSatish Balay iptr=pw = gs->pw_elm_list; 1279827bd09bSSatish Balay dptr1=dptr3 = gs->pw_vals; 1280827bd09bSSatish Balay msg_ids_in = ids_in = gs->msg_ids_in; 1281827bd09bSSatish Balay msg_ids_out = ids_out = gs->msg_ids_out; 1282827bd09bSSatish Balay dptr2 = gs->out; 1283827bd09bSSatish Balay in1=in2 = gs->in; 1284827bd09bSSatish Balay 1285827bd09bSSatish Balay /* post the receives */ 1286827bd09bSSatish Balay /* msg_nodes=nodes; */ 1287827bd09bSSatish Balay do 1288827bd09bSSatish Balay { 1289827bd09bSSatish Balay /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the 1290827bd09bSSatish Balay second one *list and do list++ afterwards */ 1291ca8e9878SJed Brown ierr = MPI_Irecv(in1, *size *step, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr); 12929182e22cSBarry Smith list++;msg_ids_in++; 1293827bd09bSSatish Balay in1 += *size++ *step; 1294827bd09bSSatish Balay } 1295827bd09bSSatish Balay while (*++msg_nodes); 1296827bd09bSSatish Balay msg_nodes=nodes; 1297827bd09bSSatish Balay 1298827bd09bSSatish Balay /* load gs values into in out gs buffers */ 1299827bd09bSSatish Balay while (*iptr >= 0) 1300827bd09bSSatish Balay { 1301ca8e9878SJed Brown PCTFS_rvec_copy(dptr3,in_vals + *iptr*step,step); 1302827bd09bSSatish Balay dptr3+=step; 1303827bd09bSSatish Balay iptr++; 1304827bd09bSSatish Balay } 1305827bd09bSSatish Balay 1306827bd09bSSatish Balay /* load out buffers and post the sends */ 1307827bd09bSSatish Balay while ((iptr = *msg_nodes++)) 1308827bd09bSSatish Balay { 1309827bd09bSSatish Balay dptr3 = dptr2; 1310827bd09bSSatish Balay while (*iptr >= 0) 1311827bd09bSSatish Balay { 1312ca8e9878SJed Brown PCTFS_rvec_copy(dptr2,dptr1 + *iptr*step,step); 1313827bd09bSSatish Balay dptr2+=step; 1314827bd09bSSatish Balay iptr++; 1315827bd09bSSatish Balay } 1316ca8e9878SJed 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); 13179182e22cSBarry Smith msg_size++; msg_list++;msg_ids_out++; 1318827bd09bSSatish Balay } 1319827bd09bSSatish Balay 1320827bd09bSSatish Balay /* tree */ 1321827bd09bSSatish Balay if (gs->max_left_over) 1322ca8e9878SJed Brown {PCTFS_gs_gop_vec_tree_plus(gs,in_vals,step);} 1323827bd09bSSatish Balay 1324827bd09bSSatish Balay /* process the received data */ 1325827bd09bSSatish Balay msg_nodes=nodes; 1326a501084fSBarry Smith while ((iptr = *nodes++)){ 1327a501084fSBarry Smith PetscScalar d1 = 1.0; 1328827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1329827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 13309182e22cSBarry Smith ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr); 13319182e22cSBarry Smith ids_in++; 1332a501084fSBarry Smith while (*iptr >= 0) { 13330805154bSBarry Smith dstep = PetscBLASIntCast(step); 13344a0de3f6SBarry Smith BLASaxpy_(&dstep,&d1,in2,&i1,dptr1 + *iptr*step,&i1); 1335827bd09bSSatish Balay in2+=step; 1336827bd09bSSatish Balay iptr++; 1337827bd09bSSatish Balay } 1338827bd09bSSatish Balay } 1339827bd09bSSatish Balay 1340827bd09bSSatish Balay /* replace vals */ 1341827bd09bSSatish Balay while (*pw >= 0) 1342827bd09bSSatish Balay { 1343ca8e9878SJed Brown PCTFS_rvec_copy(in_vals + *pw*step,dptr1,step); 1344827bd09bSSatish Balay dptr1+=step; 1345827bd09bSSatish Balay pw++; 1346827bd09bSSatish Balay } 1347827bd09bSSatish Balay 1348827bd09bSSatish Balay /* clear isend message handles */ 1349827bd09bSSatish Balay /* This changed for clarity though it could be the same */ 1350827bd09bSSatish Balay while (*msg_nodes++) 1351827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1352827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 13539182e22cSBarry Smith {ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr);ids_out++;} 13549182e22cSBarry Smith 1355827bd09bSSatish Balay 13563fdc5746SBarry Smith PetscFunctionReturn(0); 1357827bd09bSSatish Balay } 1358827bd09bSSatish Balay 13597b1ae94cSBarry Smith /******************************************************************************/ 1360ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_vec_tree_plus( PCTFS_gs_id *gs, PetscScalar *vals, PetscInt step) 1361827bd09bSSatish Balay { 136252f87cdaSBarry Smith PetscInt size, *in, *out; 1363a501084fSBarry Smith PetscScalar *buf, *work; 136452f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 1365a501084fSBarry Smith PetscBLASInt i1 = 1; 1366827bd09bSSatish Balay 13673fdc5746SBarry Smith PetscFunctionBegin; 1368827bd09bSSatish Balay /* copy over to local variables */ 1369827bd09bSSatish Balay in = gs->tree_map_in; 1370827bd09bSSatish Balay out = gs->tree_map_out; 1371827bd09bSSatish Balay buf = gs->tree_buf; 1372827bd09bSSatish Balay work = gs->tree_work; 1373827bd09bSSatish Balay size = gs->tree_nel*step; 1374827bd09bSSatish Balay 1375827bd09bSSatish Balay /* zero out collection buffer */ 1376ca8e9878SJed Brown PCTFS_rvec_zero(buf,size); 1377827bd09bSSatish Balay 1378827bd09bSSatish Balay 1379827bd09bSSatish Balay /* copy over my contributions */ 1380827bd09bSSatish Balay while (*in >= 0) 1381827bd09bSSatish Balay { 13820805154bSBarry Smith PetscBLASInt dstep = PetscBLASIntCast(step); 13836e4f4d19SBarry Smith BLAScopy_(&dstep,vals + *in++*step,&i1,buf + *out++*step,&i1); 1384827bd09bSSatish Balay } 1385827bd09bSSatish Balay 1386827bd09bSSatish Balay /* perform fan in/out on full buffer */ 1387b1c944f5SJed Brown /* must change PCTFS_grop to handle the blas */ 1388b1c944f5SJed Brown PCTFS_grop(buf,work,size,op); 1389827bd09bSSatish Balay 1390827bd09bSSatish Balay /* reset */ 1391827bd09bSSatish Balay in = gs->tree_map_in; 1392827bd09bSSatish Balay out = gs->tree_map_out; 1393827bd09bSSatish Balay 1394827bd09bSSatish Balay /* get the portion of the results I need */ 1395827bd09bSSatish Balay while (*in >= 0) 1396827bd09bSSatish Balay { 13970805154bSBarry Smith PetscBLASInt dstep = PetscBLASIntCast(step); 13986e4f4d19SBarry Smith BLAScopy_(&dstep,buf + *out++*step,&i1,vals + *in++*step,&i1); 1399827bd09bSSatish Balay } 14003fdc5746SBarry Smith PetscFunctionReturn(0); 1401827bd09bSSatish Balay } 1402827bd09bSSatish Balay 14037b1ae94cSBarry Smith /******************************************************************************/ 1404ca8e9878SJed Brown PetscErrorCode PCTFS_gs_gop_hc( PCTFS_gs_id *gs, PetscScalar *vals, const char *op, PetscInt dim) 1405827bd09bSSatish Balay { 1406d1528f56SBarry Smith PetscErrorCode ierr; 1407d1528f56SBarry Smith 14083fdc5746SBarry Smith PetscFunctionBegin; 1409827bd09bSSatish Balay switch (*op) { 1410827bd09bSSatish Balay case '+': 1411ca8e9878SJed Brown PCTFS_gs_gop_plus_hc(gs,vals,dim); 1412827bd09bSSatish Balay break; 1413827bd09bSSatish Balay default: 1414ca8e9878SJed Brown ierr = PetscInfo1(0,"PCTFS_gs_gop_hc() :: %c is not a valid op",op[0]);CHKERRQ(ierr); 1415ca8e9878SJed Brown ierr = PetscInfo(0,"PCTFS_gs_gop_hc() :: default :: plus\n");CHKERRQ(ierr); 1416ca8e9878SJed Brown PCTFS_gs_gop_plus_hc(gs,vals,dim); 1417827bd09bSSatish Balay break; 1418827bd09bSSatish Balay } 14193fdc5746SBarry Smith PetscFunctionReturn(0); 1420827bd09bSSatish Balay } 1421827bd09bSSatish Balay 14227b1ae94cSBarry Smith /******************************************************************************/ 1423ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_plus_hc( PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim) 1424827bd09bSSatish Balay { 14253fdc5746SBarry Smith PetscFunctionBegin; 1426827bd09bSSatish Balay /* if there's nothing to do return */ 1427827bd09bSSatish Balay if (dim<=0) 14283fdc5746SBarry Smith { PetscFunctionReturn(0);} 1429827bd09bSSatish Balay 1430827bd09bSSatish Balay /* can't do more dimensions then exist */ 1431b1c944f5SJed Brown dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 1432827bd09bSSatish Balay 1433827bd09bSSatish Balay /* local only operations!!! */ 1434827bd09bSSatish Balay if (gs->num_local) 1435ca8e9878SJed Brown {PCTFS_gs_gop_local_plus(gs,vals);} 1436827bd09bSSatish Balay 1437827bd09bSSatish Balay /* if intersection tree/pairwise and local isn't empty */ 1438827bd09bSSatish Balay if (gs->num_local_gop) 1439827bd09bSSatish Balay { 1440ca8e9878SJed Brown PCTFS_gs_gop_local_in_plus(gs,vals); 1441827bd09bSSatish Balay 1442827bd09bSSatish Balay /* pairwise will do tree inside ... */ 1443827bd09bSSatish Balay if (gs->num_pairs) 1444ca8e9878SJed Brown {PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim);} 1445827bd09bSSatish Balay 1446827bd09bSSatish Balay /* tree only */ 1447827bd09bSSatish Balay else if (gs->max_left_over) 1448ca8e9878SJed Brown {PCTFS_gs_gop_tree_plus_hc(gs,vals,dim);} 1449827bd09bSSatish Balay 1450ca8e9878SJed Brown PCTFS_gs_gop_local_out(gs,vals); 1451827bd09bSSatish Balay } 1452827bd09bSSatish Balay /* if intersection tree/pairwise and local is empty */ 1453827bd09bSSatish Balay else 1454827bd09bSSatish Balay { 1455827bd09bSSatish Balay /* pairwise will do tree inside */ 1456827bd09bSSatish Balay if (gs->num_pairs) 1457ca8e9878SJed Brown {PCTFS_gs_gop_pairwise_plus_hc(gs,vals,dim);} 1458827bd09bSSatish Balay 1459827bd09bSSatish Balay /* tree */ 1460827bd09bSSatish Balay else if (gs->max_left_over) 1461ca8e9878SJed Brown {PCTFS_gs_gop_tree_plus_hc(gs,vals,dim);} 1462827bd09bSSatish Balay } 14633fdc5746SBarry Smith PetscFunctionReturn(0); 1464827bd09bSSatish Balay } 1465827bd09bSSatish Balay 14667b1ae94cSBarry Smith /******************************************************************************/ 1467ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_pairwise_plus_hc( PCTFS_gs_id *gs, PetscScalar *in_vals, PetscInt dim) 1468827bd09bSSatish Balay { 1469a501084fSBarry Smith PetscScalar *dptr1, *dptr2, *dptr3, *in1, *in2; 147052f87cdaSBarry Smith PetscInt *iptr, *msg_list, *msg_size, **msg_nodes; 147152f87cdaSBarry Smith PetscInt *pw, *list, *size, **nodes; 1472827bd09bSSatish Balay MPI_Request *msg_ids_in, *msg_ids_out, *ids_in, *ids_out; 1473827bd09bSSatish Balay MPI_Status status; 147452f87cdaSBarry Smith PetscInt i, mask=1; 14753fdc5746SBarry Smith PetscErrorCode ierr; 1476827bd09bSSatish Balay 14773fdc5746SBarry Smith PetscFunctionBegin; 1478827bd09bSSatish Balay for (i=1; i<dim; i++) 1479827bd09bSSatish Balay {mask<<=1; mask++;} 1480827bd09bSSatish Balay 1481827bd09bSSatish Balay 1482a501084fSBarry Smith /* strip and load s */ 1483827bd09bSSatish Balay msg_list =list = gs->pair_list; 1484827bd09bSSatish Balay msg_size =size = gs->msg_sizes; 1485827bd09bSSatish Balay msg_nodes=nodes = gs->node_list; 1486827bd09bSSatish Balay iptr=pw = gs->pw_elm_list; 1487827bd09bSSatish Balay dptr1=dptr3 = gs->pw_vals; 1488827bd09bSSatish Balay msg_ids_in = ids_in = gs->msg_ids_in; 1489827bd09bSSatish Balay msg_ids_out = ids_out = gs->msg_ids_out; 1490827bd09bSSatish Balay dptr2 = gs->out; 1491827bd09bSSatish Balay in1=in2 = gs->in; 1492827bd09bSSatish Balay 1493827bd09bSSatish Balay /* post the receives */ 1494827bd09bSSatish Balay /* msg_nodes=nodes; */ 1495827bd09bSSatish Balay do 1496827bd09bSSatish Balay { 1497827bd09bSSatish Balay /* Should MPI_ANY_SOURCE be replaced by *list ? In that case do the 1498827bd09bSSatish Balay second one *list and do list++ afterwards */ 1499b1c944f5SJed Brown if ((PCTFS_my_id|mask)==(*list|mask)) 1500827bd09bSSatish Balay { 1501ca8e9878SJed Brown ierr = MPI_Irecv(in1, *size, MPIU_SCALAR, MPI_ANY_SOURCE, MSGTAG1 + *list, gs->PCTFS_gs_comm, msg_ids_in);CHKERRQ(ierr); 15029182e22cSBarry Smith list++; msg_ids_in++;in1 += *size++; 1503827bd09bSSatish Balay } 1504827bd09bSSatish Balay else 1505827bd09bSSatish Balay {list++; size++;} 1506827bd09bSSatish Balay } 1507827bd09bSSatish Balay while (*++msg_nodes); 1508827bd09bSSatish Balay 1509827bd09bSSatish Balay /* load gs values into in out gs buffers */ 1510827bd09bSSatish Balay while (*iptr >= 0) 1511827bd09bSSatish Balay {*dptr3++ = *(in_vals + *iptr++);} 1512827bd09bSSatish Balay 1513827bd09bSSatish Balay /* load out buffers and post the sends */ 1514827bd09bSSatish Balay msg_nodes=nodes; 1515827bd09bSSatish Balay list = msg_list; 1516827bd09bSSatish Balay while ((iptr = *msg_nodes++)) 1517827bd09bSSatish Balay { 1518b1c944f5SJed Brown if ((PCTFS_my_id|mask)==(*list|mask)) 1519827bd09bSSatish Balay { 1520827bd09bSSatish Balay dptr3 = dptr2; 1521827bd09bSSatish Balay while (*iptr >= 0) 1522827bd09bSSatish Balay {*dptr2++ = *(dptr1 + *iptr++);} 1523827bd09bSSatish Balay /* CHECK PERSISTENT COMMS MODE FOR ALL THIS STUFF */ 1524827bd09bSSatish Balay /* is msg_ids_out++ correct? */ 1525ca8e9878SJed Brown ierr = MPI_Isend(dptr3, *msg_size, MPIU_SCALAR, *list, MSGTAG1+PCTFS_my_id, gs->PCTFS_gs_comm, msg_ids_out);CHKERRQ(ierr); 15269182e22cSBarry Smith msg_size++;list++;msg_ids_out++; 1527827bd09bSSatish Balay } 1528827bd09bSSatish Balay else 1529827bd09bSSatish Balay {list++; msg_size++;} 1530827bd09bSSatish Balay } 1531827bd09bSSatish Balay 1532827bd09bSSatish Balay /* do the tree while we're waiting */ 1533827bd09bSSatish Balay if (gs->max_left_over) 1534ca8e9878SJed Brown {PCTFS_gs_gop_tree_plus_hc(gs,in_vals,dim);} 1535827bd09bSSatish Balay 1536827bd09bSSatish Balay /* process the received data */ 1537827bd09bSSatish Balay msg_nodes=nodes; 1538827bd09bSSatish Balay list = msg_list; 1539827bd09bSSatish Balay while ((iptr = *nodes++)) 1540827bd09bSSatish Balay { 1541b1c944f5SJed Brown if ((PCTFS_my_id|mask)==(*list|mask)) 1542827bd09bSSatish Balay { 1543827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1544827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 15459182e22cSBarry Smith ierr = MPI_Wait(ids_in, &status);CHKERRQ(ierr); 15469182e22cSBarry Smith ids_in++; 1547827bd09bSSatish Balay while (*iptr >= 0) 1548827bd09bSSatish Balay {*(dptr1 + *iptr++) += *in2++;} 1549827bd09bSSatish Balay } 1550827bd09bSSatish Balay list++; 1551827bd09bSSatish Balay } 1552827bd09bSSatish Balay 1553827bd09bSSatish Balay /* replace vals */ 1554827bd09bSSatish Balay while (*pw >= 0) 1555827bd09bSSatish Balay {*(in_vals + *pw++) = *dptr1++;} 1556827bd09bSSatish Balay 1557827bd09bSSatish Balay /* clear isend message handles */ 1558827bd09bSSatish Balay /* This changed for clarity though it could be the same */ 1559827bd09bSSatish Balay while (*msg_nodes++) 1560827bd09bSSatish Balay { 1561b1c944f5SJed Brown if ((PCTFS_my_id|mask)==(*msg_list|mask)) 1562827bd09bSSatish Balay { 1563827bd09bSSatish Balay /* Should I check the return value of MPI_Wait() or status? */ 1564827bd09bSSatish Balay /* Can this loop be replaced by a call to MPI_Waitall()? */ 15659182e22cSBarry Smith ierr = MPI_Wait(ids_out, &status);CHKERRQ(ierr); 15669182e22cSBarry Smith ids_out++; 1567827bd09bSSatish Balay } 1568827bd09bSSatish Balay msg_list++; 1569827bd09bSSatish Balay } 1570827bd09bSSatish Balay 15713fdc5746SBarry Smith PetscFunctionReturn(0); 1572827bd09bSSatish Balay } 1573827bd09bSSatish Balay 15747b1ae94cSBarry Smith /******************************************************************************/ 1575ca8e9878SJed Brown static PetscErrorCode PCTFS_gs_gop_tree_plus_hc(PCTFS_gs_id *gs, PetscScalar *vals, PetscInt dim) 1576827bd09bSSatish Balay { 157752f87cdaSBarry Smith PetscInt size; 157852f87cdaSBarry Smith PetscInt *in, *out; 1579a501084fSBarry Smith PetscScalar *buf, *work; 158052f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 1581827bd09bSSatish Balay 15823fdc5746SBarry Smith PetscFunctionBegin; 1583827bd09bSSatish Balay in = gs->tree_map_in; 1584827bd09bSSatish Balay out = gs->tree_map_out; 1585827bd09bSSatish Balay buf = gs->tree_buf; 1586827bd09bSSatish Balay work = gs->tree_work; 1587827bd09bSSatish Balay size = gs->tree_nel; 1588827bd09bSSatish Balay 1589ca8e9878SJed Brown PCTFS_rvec_zero(buf,size); 1590827bd09bSSatish Balay 1591827bd09bSSatish Balay while (*in >= 0) 1592827bd09bSSatish Balay {*(buf + *out++) = *(vals + *in++);} 1593827bd09bSSatish Balay 1594827bd09bSSatish Balay in = gs->tree_map_in; 1595827bd09bSSatish Balay out = gs->tree_map_out; 1596827bd09bSSatish Balay 1597b1c944f5SJed Brown PCTFS_grop_hc(buf,work,size,op,dim); 1598827bd09bSSatish Balay 1599827bd09bSSatish Balay while (*in >= 0) 1600827bd09bSSatish Balay {*(vals + *in++) = *(buf + *out++);} 16013fdc5746SBarry Smith PetscFunctionReturn(0); 1602827bd09bSSatish Balay } 1603827bd09bSSatish Balay 1604827bd09bSSatish Balay 1605827bd09bSSatish Balay 1606