1827bd09bSSatish Balay 2827bd09bSSatish Balay /***********************************comm.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 11.21.97 15827bd09bSSatish Balay ***********************************comm.c*************************************/ 16c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h> 17827bd09bSSatish Balay 18827bd09bSSatish Balay 19827bd09bSSatish Balay /* global program control variables - explicitly exported */ 20*b1c944f5SJed Brown PetscMPIInt PCTFS_my_id = 0; 21*b1c944f5SJed Brown PetscMPIInt PCTFS_num_nodes = 1; 22*b1c944f5SJed Brown PetscMPIInt PCTFS_floor_num_nodes = 0; 23*b1c944f5SJed Brown PetscMPIInt PCTFS_i_log2_num_nodes = 0; 24827bd09bSSatish Balay 25827bd09bSSatish Balay /* global program control variables */ 2652f87cdaSBarry Smith static PetscInt p_init = 0; 2752f87cdaSBarry Smith static PetscInt modfl_num_nodes; 2852f87cdaSBarry Smith static PetscInt edge_not_pow_2; 29827bd09bSSatish Balay 3052f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt)*32]; 31827bd09bSSatish Balay 327b1ae94cSBarry Smith /***********************************comm.c*************************************/ 33*b1c944f5SJed Brown PetscErrorCode PCTFS_comm_init (void) 34827bd09bSSatish Balay { 35827bd09bSSatish Balay 363fdc5746SBarry Smith if (p_init++) PetscFunctionReturn(0); 37827bd09bSSatish Balay 38*b1c944f5SJed Brown MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes); 39*b1c944f5SJed Brown MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id); 40827bd09bSSatish Balay 41*b1c944f5SJed Brown if (PCTFS_num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!"); 42827bd09bSSatish Balay 433fdc5746SBarry Smith ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32); 44827bd09bSSatish Balay 45*b1c944f5SJed Brown PCTFS_floor_num_nodes = 1; 46*b1c944f5SJed Brown PCTFS_i_log2_num_nodes = modfl_num_nodes = 0; 47*b1c944f5SJed Brown while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) 48827bd09bSSatish Balay { 49*b1c944f5SJed Brown edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes; 50*b1c944f5SJed Brown PCTFS_floor_num_nodes <<= 1; 51*b1c944f5SJed Brown PCTFS_i_log2_num_nodes++; 52827bd09bSSatish Balay } 53827bd09bSSatish Balay 54*b1c944f5SJed Brown PCTFS_i_log2_num_nodes--; 55*b1c944f5SJed Brown PCTFS_floor_num_nodes >>= 1; 56*b1c944f5SJed Brown modfl_num_nodes = (PCTFS_num_nodes - PCTFS_floor_num_nodes); 57827bd09bSSatish Balay 58*b1c944f5SJed Brown if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) 59*b1c944f5SJed Brown {edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1);} 60*b1c944f5SJed Brown else if (PCTFS_my_id >= PCTFS_floor_num_nodes) 61*b1c944f5SJed Brown {edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1); 62827bd09bSSatish Balay } 63827bd09bSSatish Balay else 64827bd09bSSatish Balay {edge_not_pow_2 = 0;} 653fdc5746SBarry Smith PetscFunctionReturn(0); 66827bd09bSSatish Balay } 67827bd09bSSatish Balay 687b1ae94cSBarry Smith /***********************************comm.c*************************************/ 69*b1c944f5SJed Brown PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs) 70827bd09bSSatish Balay { 713fdc5746SBarry Smith PetscInt mask, edge; 723fdc5746SBarry Smith PetscInt type, dest; 73827bd09bSSatish Balay vfp fp; 74827bd09bSSatish Balay MPI_Status status; 753fdc5746SBarry Smith PetscInt ierr; 76827bd09bSSatish Balay 773fdc5746SBarry Smith PetscFunctionBegin; 78827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 79*b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 80827bd09bSSatish Balay 81827bd09bSSatish Balay /* non-uniform should have at least two entries */ 82*b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?"); 83827bd09bSSatish Balay 84827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 85827bd09bSSatish Balay if (!p_init) 86*b1c944f5SJed Brown {PCTFS_comm_init();} 87827bd09bSSatish Balay 88827bd09bSSatish Balay /* if there's nothing to do return */ 89*b1c944f5SJed Brown if ((PCTFS_num_nodes<2)||(!n)) 90827bd09bSSatish Balay { 913fdc5746SBarry Smith PetscFunctionReturn(0); 92827bd09bSSatish Balay } 93827bd09bSSatish Balay 9471a0148aSBarry Smith 95827bd09bSSatish Balay /* a negative number if items to send ==> fatal */ 96*b1c944f5SJed Brown if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n); 97827bd09bSSatish Balay 98827bd09bSSatish Balay /* advance to list of n operations for custom */ 99827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 100827bd09bSSatish Balay {oprs++;} 101827bd09bSSatish Balay 102827bd09bSSatish Balay /* major league hack */ 103d890fc11SSatish Balay if (!(fp = (vfp) ivec_fct_addr(type))) { 104*b1c944f5SJed Brown ierr = PetscInfo(0,"PCTFS_giop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 105827bd09bSSatish Balay fp = (vfp) oprs; 106827bd09bSSatish Balay } 107827bd09bSSatish Balay 108827bd09bSSatish Balay /* all msgs will be of the same length */ 109827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 110827bd09bSSatish Balay if (edge_not_pow_2) 111827bd09bSSatish Balay { 112*b1c944f5SJed Brown if (PCTFS_my_id >= PCTFS_floor_num_nodes) 113*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 114827bd09bSSatish Balay else 115827bd09bSSatish Balay { 1163fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr); 117827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 118827bd09bSSatish Balay } 119827bd09bSSatish Balay } 120827bd09bSSatish Balay 121827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 122*b1c944f5SJed Brown if (PCTFS_my_id<PCTFS_floor_num_nodes) 123827bd09bSSatish Balay { 124*b1c944f5SJed Brown for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) 125827bd09bSSatish Balay { 126*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 127*b1c944f5SJed Brown if (PCTFS_my_id > dest) 128*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 129827bd09bSSatish Balay else 130827bd09bSSatish Balay { 1313fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 132827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 133827bd09bSSatish Balay } 134827bd09bSSatish Balay } 135827bd09bSSatish Balay 136*b1c944f5SJed Brown mask=PCTFS_floor_num_nodes>>1; 137*b1c944f5SJed Brown for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) 138827bd09bSSatish Balay { 139*b1c944f5SJed Brown if (PCTFS_my_id%mask) 140827bd09bSSatish Balay {continue;} 141827bd09bSSatish Balay 142*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 143*b1c944f5SJed Brown if (PCTFS_my_id < dest) 144*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 145827bd09bSSatish Balay else 146827bd09bSSatish Balay { 1473fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 148827bd09bSSatish Balay } 149827bd09bSSatish Balay } 150827bd09bSSatish Balay } 151827bd09bSSatish Balay 152827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 153827bd09bSSatish Balay if (edge_not_pow_2) 154827bd09bSSatish Balay { 155*b1c944f5SJed Brown if (PCTFS_my_id >= PCTFS_floor_num_nodes) 156827bd09bSSatish Balay { 1573fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 158827bd09bSSatish Balay } 159827bd09bSSatish Balay else 160*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 161827bd09bSSatish Balay } 1623fdc5746SBarry Smith PetscFunctionReturn(0); 163827bd09bSSatish Balay } 164827bd09bSSatish Balay 1657b1ae94cSBarry Smith /***********************************comm.c*************************************/ 166*b1c944f5SJed Brown PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs) 167827bd09bSSatish Balay { 1683fdc5746SBarry Smith PetscInt mask, edge; 1693fdc5746SBarry Smith PetscInt type, dest; 170827bd09bSSatish Balay vfp fp; 171827bd09bSSatish Balay MPI_Status status; 1723fdc5746SBarry Smith PetscErrorCode ierr; 173827bd09bSSatish Balay 1743fdc5746SBarry Smith PetscFunctionBegin; 175827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 176*b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 177827bd09bSSatish Balay 178827bd09bSSatish Balay /* non-uniform should have at least two entries */ 179*b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?"); 180827bd09bSSatish Balay 181827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 182827bd09bSSatish Balay if (!p_init) 183*b1c944f5SJed Brown {PCTFS_comm_init();} 184827bd09bSSatish Balay 185827bd09bSSatish Balay /* if there's nothing to do return */ 186*b1c944f5SJed Brown if ((PCTFS_num_nodes<2)||(!n)) 1873fdc5746SBarry Smith { PetscFunctionReturn(0);} 188827bd09bSSatish Balay 189827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 190c1235816SBarry Smith if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n); 191827bd09bSSatish Balay 192827bd09bSSatish Balay /* advance to list of n operations for custom */ 193827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 194827bd09bSSatish Balay {oprs++;} 195827bd09bSSatish Balay 196d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))) { 197*b1c944f5SJed Brown ierr = PetscInfo(0,"PCTFS_grop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 198827bd09bSSatish Balay fp = (vfp) oprs; 199827bd09bSSatish Balay } 200827bd09bSSatish Balay 201827bd09bSSatish Balay /* all msgs will be of the same length */ 202827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 203827bd09bSSatish Balay if (edge_not_pow_2) 204827bd09bSSatish Balay { 205*b1c944f5SJed Brown if (PCTFS_my_id >= PCTFS_floor_num_nodes) 206*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 207827bd09bSSatish Balay else 208827bd09bSSatish Balay { 2093fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 210827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 211827bd09bSSatish Balay } 212827bd09bSSatish Balay } 213827bd09bSSatish Balay 214827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 215*b1c944f5SJed Brown if (PCTFS_my_id<PCTFS_floor_num_nodes) 216827bd09bSSatish Balay { 217*b1c944f5SJed Brown for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) 218827bd09bSSatish Balay { 219*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 220*b1c944f5SJed Brown if (PCTFS_my_id > dest) 221*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 222827bd09bSSatish Balay else 223827bd09bSSatish Balay { 2243fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 225827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 226827bd09bSSatish Balay } 227827bd09bSSatish Balay } 228827bd09bSSatish Balay 229*b1c944f5SJed Brown mask=PCTFS_floor_num_nodes>>1; 230*b1c944f5SJed Brown for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) 231827bd09bSSatish Balay { 232*b1c944f5SJed Brown if (PCTFS_my_id%mask) 233827bd09bSSatish Balay {continue;} 234827bd09bSSatish Balay 235*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 236*b1c944f5SJed Brown if (PCTFS_my_id < dest) 237*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 238827bd09bSSatish Balay else 239827bd09bSSatish Balay { 2403fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 241827bd09bSSatish Balay } 242827bd09bSSatish Balay } 243827bd09bSSatish Balay } 244827bd09bSSatish Balay 245827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 246827bd09bSSatish Balay if (edge_not_pow_2) 247827bd09bSSatish Balay { 248*b1c944f5SJed Brown if (PCTFS_my_id >= PCTFS_floor_num_nodes) 249827bd09bSSatish Balay { 2503fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr); 251827bd09bSSatish Balay } 252827bd09bSSatish Balay else 253*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 254827bd09bSSatish Balay } 2553fdc5746SBarry Smith PetscFunctionReturn(0); 256827bd09bSSatish Balay } 257827bd09bSSatish Balay 2587b1ae94cSBarry Smith /***********************************comm.c*************************************/ 259*b1c944f5SJed Brown PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim) 260827bd09bSSatish Balay { 2613fdc5746SBarry Smith PetscInt mask, edge; 2623fdc5746SBarry Smith PetscInt type, dest; 263827bd09bSSatish Balay vfp fp; 264827bd09bSSatish Balay MPI_Status status; 2653fdc5746SBarry Smith PetscErrorCode ierr; 266827bd09bSSatish Balay 2673fdc5746SBarry Smith PetscFunctionBegin; 268827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 269*b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 270827bd09bSSatish Balay 271827bd09bSSatish Balay /* non-uniform should have at least two entries */ 272*b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?"); 273827bd09bSSatish Balay 274827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 275827bd09bSSatish Balay if (!p_init) 276*b1c944f5SJed Brown {PCTFS_comm_init();} 277827bd09bSSatish Balay 278827bd09bSSatish Balay /* if there's nothing to do return */ 279*b1c944f5SJed Brown if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) 2800924e98cSBarry Smith {PetscFunctionReturn(0);} 281827bd09bSSatish Balay 282827bd09bSSatish Balay /* the error msg says it all!!! */ 283*b1c944f5SJed Brown if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!"); 284827bd09bSSatish Balay 285827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 286*b1c944f5SJed Brown if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n); 287827bd09bSSatish Balay 288827bd09bSSatish Balay /* can't do more dimensions then exist */ 289*b1c944f5SJed Brown dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 290827bd09bSSatish Balay 291827bd09bSSatish Balay /* advance to list of n operations for custom */ 292827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 293827bd09bSSatish Balay {oprs++;} 294827bd09bSSatish Balay 295d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))) { 296*b1c944f5SJed Brown ierr = PetscInfo(0,"PCTFS_grop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 297827bd09bSSatish Balay fp = (vfp) oprs; 298827bd09bSSatish Balay } 299827bd09bSSatish Balay 300827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 301827bd09bSSatish Balay { 302*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 303*b1c944f5SJed Brown if (PCTFS_my_id > dest) 304*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 305827bd09bSSatish Balay else 306827bd09bSSatish Balay { 3073fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 308827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 309827bd09bSSatish Balay } 310827bd09bSSatish Balay } 311827bd09bSSatish Balay 312827bd09bSSatish Balay if (edge==dim) 313827bd09bSSatish Balay {mask>>=1;} 314827bd09bSSatish Balay else 315827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 316827bd09bSSatish Balay 317827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 318827bd09bSSatish Balay { 319*b1c944f5SJed Brown if (PCTFS_my_id%mask) 320827bd09bSSatish Balay {continue;} 321827bd09bSSatish Balay 322*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 323*b1c944f5SJed Brown if (PCTFS_my_id < dest) 324*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 325827bd09bSSatish Balay else 326827bd09bSSatish Balay { 3273fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 328827bd09bSSatish Balay } 329827bd09bSSatish Balay } 3303fdc5746SBarry Smith PetscFunctionReturn(0); 331827bd09bSSatish Balay } 332827bd09bSSatish Balay 3337b1ae94cSBarry Smith /******************************************************************************/ 334*b1c944f5SJed Brown PetscErrorCode PCTFS_ssgl_radd( PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs) 335827bd09bSSatish Balay { 3363fdc5746SBarry Smith PetscInt edge, type, dest, mask; 3373fdc5746SBarry Smith PetscInt stage_n; 338827bd09bSSatish Balay MPI_Status status; 3393fdc5746SBarry Smith PetscErrorCode ierr; 340827bd09bSSatish Balay 3413fdc5746SBarry Smith PetscFunctionBegin; 342827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 343827bd09bSSatish Balay if (!p_init) 344*b1c944f5SJed Brown {PCTFS_comm_init();} 345827bd09bSSatish Balay 346827bd09bSSatish Balay 347827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 348827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 349827bd09bSSatish Balay for (mask=0, edge=0; edge<level; edge++, mask++) 350827bd09bSSatish Balay { 351827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 352*b1c944f5SJed Brown if (stage_n && !(PCTFS_my_id & mask)) 353827bd09bSSatish Balay { 354827bd09bSSatish Balay dest = edge_node[edge]; 355*b1c944f5SJed Brown type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge); 356*b1c944f5SJed Brown if (PCTFS_my_id>dest) 3573fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);} 358827bd09bSSatish Balay else 359827bd09bSSatish Balay { 360*b1c944f5SJed Brown type = type - PCTFS_my_id + dest; 3613fdc5746SBarry Smith ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 362827bd09bSSatish Balay rvec_add(vals+segs[edge], work, stage_n); 363827bd09bSSatish Balay } 364827bd09bSSatish Balay } 365827bd09bSSatish Balay mask <<= 1; 366827bd09bSSatish Balay } 367827bd09bSSatish Balay mask>>=1; 368827bd09bSSatish Balay for (edge=0; edge<level; edge++) 369827bd09bSSatish Balay { 370827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 371*b1c944f5SJed Brown if (stage_n && !(PCTFS_my_id & mask)) 372827bd09bSSatish Balay { 373827bd09bSSatish Balay dest = edge_node[level-edge-1]; 374*b1c944f5SJed Brown type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge); 375*b1c944f5SJed Brown if (PCTFS_my_id<dest) 3763fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);} 377827bd09bSSatish Balay else 378827bd09bSSatish Balay { 379*b1c944f5SJed Brown type = type - PCTFS_my_id + dest; 3803fdc5746SBarry Smith ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 381827bd09bSSatish Balay } 382827bd09bSSatish Balay } 383827bd09bSSatish Balay mask >>= 1; 384827bd09bSSatish Balay } 3853fdc5746SBarry Smith PetscFunctionReturn(0); 386827bd09bSSatish Balay } 387827bd09bSSatish Balay 3887b1ae94cSBarry Smith /******************************************************************************/ 38952f87cdaSBarry Smith PetscErrorCode new_ssgl_radd( PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs) 390827bd09bSSatish Balay { 39152f87cdaSBarry Smith PetscInt edge, type, dest, mask; 39252f87cdaSBarry Smith PetscInt stage_n; 393827bd09bSSatish Balay MPI_Status status; 3943fdc5746SBarry Smith PetscErrorCode ierr; 395827bd09bSSatish Balay 3963fdc5746SBarry Smith PetscFunctionBegin; 397827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 398827bd09bSSatish Balay if (!p_init) 399*b1c944f5SJed Brown {PCTFS_comm_init();} 400827bd09bSSatish Balay 401827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 402827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 403827bd09bSSatish Balay for (mask=0, edge=0; edge<level; edge++, mask++) 404827bd09bSSatish Balay { 405827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 406*b1c944f5SJed Brown if (stage_n && !(PCTFS_my_id & mask)) 407827bd09bSSatish Balay { 408827bd09bSSatish Balay dest = edge_node[edge]; 409*b1c944f5SJed Brown type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge); 410*b1c944f5SJed Brown if (PCTFS_my_id>dest) 4113fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);} 412827bd09bSSatish Balay else 413827bd09bSSatish Balay { 414*b1c944f5SJed Brown type = type - PCTFS_my_id + dest; 4153fdc5746SBarry Smith ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type, MPI_COMM_WORLD,&status);CHKERRQ(ierr); 416827bd09bSSatish Balay rvec_add(vals+segs[edge], work, stage_n); 417827bd09bSSatish Balay } 418827bd09bSSatish Balay } 419827bd09bSSatish Balay mask <<= 1; 420827bd09bSSatish Balay } 421827bd09bSSatish Balay mask>>=1; 422827bd09bSSatish Balay for (edge=0; edge<level; edge++) 423827bd09bSSatish Balay { 424827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 425*b1c944f5SJed Brown if (stage_n && !(PCTFS_my_id & mask)) 426827bd09bSSatish Balay { 427827bd09bSSatish Balay dest = edge_node[level-edge-1]; 428*b1c944f5SJed Brown type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge); 429*b1c944f5SJed Brown if (PCTFS_my_id<dest) 4303fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);} 431827bd09bSSatish Balay else 432827bd09bSSatish Balay { 433*b1c944f5SJed Brown type = type - PCTFS_my_id + dest; 4343fdc5746SBarry Smith ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 435827bd09bSSatish Balay } 436827bd09bSSatish Balay } 437827bd09bSSatish Balay mask >>= 1; 438827bd09bSSatish Balay } 4393fdc5746SBarry Smith PetscFunctionReturn(0); 440827bd09bSSatish Balay } 441827bd09bSSatish Balay 4427b1ae94cSBarry Smith /***********************************comm.c*************************************/ 443*b1c944f5SJed Brown PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim) 444827bd09bSSatish Balay { 44552f87cdaSBarry Smith PetscInt mask, edge; 44652f87cdaSBarry Smith PetscInt type, dest; 447827bd09bSSatish Balay vfp fp; 448827bd09bSSatish Balay MPI_Status status; 4493fdc5746SBarry Smith PetscErrorCode ierr; 450827bd09bSSatish Balay 4513fdc5746SBarry Smith PetscFunctionBegin; 452827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 453*b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 454827bd09bSSatish Balay 455827bd09bSSatish Balay /* non-uniform should have at least two entries */ 456*b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?"); 457827bd09bSSatish Balay 458827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 459827bd09bSSatish Balay if (!p_init) 460*b1c944f5SJed Brown {PCTFS_comm_init();} 461827bd09bSSatish Balay 462827bd09bSSatish Balay /* if there's nothing to do return */ 463*b1c944f5SJed Brown if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) 4643fdc5746SBarry Smith { PetscFunctionReturn(0);} 465827bd09bSSatish Balay 466827bd09bSSatish Balay /* the error msg says it all!!! */ 467*b1c944f5SJed Brown if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!"); 468827bd09bSSatish Balay 469827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 470*b1c944f5SJed Brown if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n); 471827bd09bSSatish Balay 472827bd09bSSatish Balay /* can't do more dimensions then exist */ 473*b1c944f5SJed Brown dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 474827bd09bSSatish Balay 475827bd09bSSatish Balay /* advance to list of n operations for custom */ 476827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 477827bd09bSSatish Balay {oprs++;} 478827bd09bSSatish Balay 479d890fc11SSatish Balay if (!(fp = (vfp) ivec_fct_addr(type))){ 480*b1c944f5SJed Brown ierr = PetscInfo(0,"PCTFS_giop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 481827bd09bSSatish Balay fp = (vfp) oprs; 482827bd09bSSatish Balay } 483827bd09bSSatish Balay 484827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 485827bd09bSSatish Balay { 486*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 487*b1c944f5SJed Brown if (PCTFS_my_id > dest) 488*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 489827bd09bSSatish Balay else 490827bd09bSSatish Balay { 4913fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 492827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 493827bd09bSSatish Balay } 494827bd09bSSatish Balay } 495827bd09bSSatish Balay 496827bd09bSSatish Balay if (edge==dim) 497827bd09bSSatish Balay {mask>>=1;} 498827bd09bSSatish Balay else 499827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 500827bd09bSSatish Balay 501827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 502827bd09bSSatish Balay { 503*b1c944f5SJed Brown if (PCTFS_my_id%mask) 504827bd09bSSatish Balay {continue;} 505827bd09bSSatish Balay 506*b1c944f5SJed Brown dest = PCTFS_my_id^mask; 507*b1c944f5SJed Brown if (PCTFS_my_id < dest) 508*b1c944f5SJed Brown {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 509827bd09bSSatish Balay else 510827bd09bSSatish Balay { 5113fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 512827bd09bSSatish Balay } 513827bd09bSSatish Balay } 5143fdc5746SBarry Smith PetscFunctionReturn(0); 515827bd09bSSatish Balay } 516