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 /* global program control variables - explicitly exported */ 19b1c944f5SJed Brown PetscMPIInt PCTFS_my_id = 0; 20b1c944f5SJed Brown PetscMPIInt PCTFS_num_nodes = 1; 21b1c944f5SJed Brown PetscMPIInt PCTFS_floor_num_nodes = 0; 22b1c944f5SJed Brown PetscMPIInt PCTFS_i_log2_num_nodes = 0; 23827bd09bSSatish Balay 24827bd09bSSatish Balay /* global program control variables */ 2552f87cdaSBarry Smith static PetscInt p_init = 0; 2652f87cdaSBarry Smith static PetscInt modfl_num_nodes; 2752f87cdaSBarry Smith static PetscInt edge_not_pow_2; 28827bd09bSSatish Balay 2952f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt)*32]; 30827bd09bSSatish Balay 317b1ae94cSBarry Smith /***********************************comm.c*************************************/ 32b1c944f5SJed Brown PetscErrorCode PCTFS_comm_init(void) 33827bd09bSSatish Balay { 34*362febeeSStefano Zampini PetscFunctionBegin; 353fdc5746SBarry Smith if (p_init++) PetscFunctionReturn(0); 36827bd09bSSatish Balay 37b1c944f5SJed Brown MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes); 38b1c944f5SJed Brown MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id); 39827bd09bSSatish Balay 40b1c944f5SJed Brown if (PCTFS_num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!"); 41827bd09bSSatish Balay 42ca8e9878SJed Brown PCTFS_ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32); 43827bd09bSSatish Balay 44b1c944f5SJed Brown PCTFS_floor_num_nodes = 1; 45b1c944f5SJed Brown PCTFS_i_log2_num_nodes = modfl_num_nodes = 0; 46db4deed7SKarl Rupp while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) { 47b1c944f5SJed Brown edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes; 48b1c944f5SJed Brown PCTFS_floor_num_nodes <<= 1; 49b1c944f5SJed Brown PCTFS_i_log2_num_nodes++; 50827bd09bSSatish Balay } 51827bd09bSSatish Balay 52b1c944f5SJed Brown PCTFS_i_log2_num_nodes--; 53b1c944f5SJed Brown PCTFS_floor_num_nodes >>= 1; 54b1c944f5SJed Brown modfl_num_nodes = (PCTFS_num_nodes - PCTFS_floor_num_nodes); 55827bd09bSSatish Balay 562fa5cd67SKarl Rupp if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1); 572fa5cd67SKarl Rupp else if (PCTFS_my_id >= PCTFS_floor_num_nodes) edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1); 582fa5cd67SKarl Rupp else edge_not_pow_2 = 0; 593fdc5746SBarry Smith PetscFunctionReturn(0); 60827bd09bSSatish Balay } 61827bd09bSSatish Balay 627b1ae94cSBarry Smith /***********************************comm.c*************************************/ 63b1c944f5SJed Brown PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs) 64827bd09bSSatish Balay { 653fdc5746SBarry Smith PetscInt mask, edge; 663fdc5746SBarry Smith PetscInt type, dest; 67827bd09bSSatish Balay vfp fp; 68827bd09bSSatish Balay MPI_Status status; 693fdc5746SBarry Smith PetscInt ierr; 70827bd09bSSatish Balay 713fdc5746SBarry Smith PetscFunctionBegin; 72827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 73b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 74827bd09bSSatish Balay 75827bd09bSSatish Balay /* non-uniform should have at least two entries */ 76b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?"); 77827bd09bSSatish Balay 78827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 792fa5cd67SKarl Rupp if (!p_init) PCTFS_comm_init(); 80827bd09bSSatish Balay 81827bd09bSSatish Balay /* if there's nothing to do return */ 822fa5cd67SKarl Rupp if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0); 8371a0148aSBarry Smith 84827bd09bSSatish Balay /* a negative number if items to send ==> fatal */ 85b1c944f5SJed Brown if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n); 86827bd09bSSatish Balay 87827bd09bSSatish Balay /* advance to list of n operations for custom */ 882fa5cd67SKarl Rupp if ((type=oprs[0])==NON_UNIFORM) oprs++; 89827bd09bSSatish Balay 90827bd09bSSatish Balay /* major league hack */ 912fa5cd67SKarl Rupp if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: Could not retrieve function pointer!\n"); 92827bd09bSSatish Balay 93827bd09bSSatish Balay /* all msgs will be of the same length */ 94827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 95db4deed7SKarl Rupp if (edge_not_pow_2) { 962fa5cd67SKarl Rupp if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 97ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 982fa5cd67SKarl Rupp } else { 99ffc4695bSBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 100827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 101827bd09bSSatish Balay } 102827bd09bSSatish Balay } 103827bd09bSSatish Balay 104827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 105db4deed7SKarl Rupp if (PCTFS_my_id<PCTFS_floor_num_nodes) { 106db4deed7SKarl Rupp for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) { 107b1c944f5SJed Brown dest = PCTFS_my_id^mask; 1082fa5cd67SKarl Rupp if (PCTFS_my_id > dest) { 109ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 1102fa5cd67SKarl Rupp } else { 111ffc4695bSBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 112827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 113827bd09bSSatish Balay } 114827bd09bSSatish Balay } 115827bd09bSSatish Balay 116b1c944f5SJed Brown mask=PCTFS_floor_num_nodes>>1; 117db4deed7SKarl Rupp for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) { 1182fa5cd67SKarl Rupp if (PCTFS_my_id%mask) continue; 119827bd09bSSatish Balay 120b1c944f5SJed Brown dest = PCTFS_my_id^mask; 1212fa5cd67SKarl Rupp if (PCTFS_my_id < dest) { 122ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 1232fa5cd67SKarl Rupp } else { 124ffc4695bSBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 125827bd09bSSatish Balay } 126827bd09bSSatish Balay } 127827bd09bSSatish Balay } 128827bd09bSSatish Balay 129827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 130db4deed7SKarl Rupp if (edge_not_pow_2) { 1312fa5cd67SKarl Rupp if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 132ffc4695bSBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 1332fa5cd67SKarl Rupp } else { 134ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 135827bd09bSSatish Balay } 136827bd09bSSatish Balay } 1373fdc5746SBarry Smith PetscFunctionReturn(0); 138827bd09bSSatish Balay } 139827bd09bSSatish Balay 1407b1ae94cSBarry Smith /***********************************comm.c*************************************/ 141b1c944f5SJed Brown PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs) 142827bd09bSSatish Balay { 1433fdc5746SBarry Smith PetscInt mask, edge; 1443fdc5746SBarry Smith PetscInt type, dest; 145827bd09bSSatish Balay vfp fp; 146827bd09bSSatish Balay MPI_Status status; 1473fdc5746SBarry Smith PetscErrorCode ierr; 148827bd09bSSatish Balay 1493fdc5746SBarry Smith PetscFunctionBegin; 150827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 151b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 152827bd09bSSatish Balay 153827bd09bSSatish Balay /* non-uniform should have at least two entries */ 154b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?"); 155827bd09bSSatish Balay 156827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 1572fa5cd67SKarl Rupp if (!p_init) PCTFS_comm_init(); 158827bd09bSSatish Balay 159827bd09bSSatish Balay /* if there's nothing to do return */ 1602fa5cd67SKarl Rupp if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0); 161827bd09bSSatish Balay 162827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 163c1235816SBarry Smith if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n); 164827bd09bSSatish Balay 165827bd09bSSatish Balay /* advance to list of n operations for custom */ 1662fa5cd67SKarl Rupp if ((type=oprs[0])==NON_UNIFORM) oprs++; 167827bd09bSSatish Balay 1682fa5cd67SKarl Rupp if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: Could not retrieve function pointer!\n"); 169827bd09bSSatish Balay 170827bd09bSSatish Balay /* all msgs will be of the same length */ 171827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 1722fa5cd67SKarl Rupp if (edge_not_pow_2) { 1732fa5cd67SKarl Rupp if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 174ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 1752fa5cd67SKarl Rupp } else { 176ffc4695bSBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 177827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 178827bd09bSSatish Balay } 179827bd09bSSatish Balay } 180827bd09bSSatish Balay 181827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 182db4deed7SKarl Rupp if (PCTFS_my_id<PCTFS_floor_num_nodes) { 183db4deed7SKarl Rupp for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) { 184b1c944f5SJed Brown dest = PCTFS_my_id^mask; 1852fa5cd67SKarl Rupp if (PCTFS_my_id > dest) { 186ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 1872fa5cd67SKarl Rupp } else { 188ffc4695bSBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 189827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 190827bd09bSSatish Balay } 191827bd09bSSatish Balay } 192827bd09bSSatish Balay 193b1c944f5SJed Brown mask=PCTFS_floor_num_nodes>>1; 194db4deed7SKarl Rupp for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) { 1952fa5cd67SKarl Rupp if (PCTFS_my_id%mask) continue; 196827bd09bSSatish Balay 197b1c944f5SJed Brown dest = PCTFS_my_id^mask; 1982fa5cd67SKarl Rupp if (PCTFS_my_id < dest) { 199ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 2002fa5cd67SKarl Rupp } else { 201ffc4695bSBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 202827bd09bSSatish Balay } 203827bd09bSSatish Balay } 204827bd09bSSatish Balay } 205827bd09bSSatish Balay 206827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 207db4deed7SKarl Rupp if (edge_not_pow_2) { 208db4deed7SKarl Rupp if (PCTFS_my_id >= PCTFS_floor_num_nodes) { 209ffc4695bSBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 210db4deed7SKarl Rupp } else { 211ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 212827bd09bSSatish Balay } 213827bd09bSSatish Balay } 2143fdc5746SBarry Smith PetscFunctionReturn(0); 215827bd09bSSatish Balay } 216827bd09bSSatish Balay 2177b1ae94cSBarry Smith /***********************************comm.c*************************************/ 218b1c944f5SJed Brown PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim) 219827bd09bSSatish Balay { 2203fdc5746SBarry Smith PetscInt mask, edge; 2213fdc5746SBarry Smith PetscInt type, dest; 222827bd09bSSatish Balay vfp fp; 223827bd09bSSatish Balay MPI_Status status; 2243fdc5746SBarry Smith PetscErrorCode ierr; 225827bd09bSSatish Balay 2263fdc5746SBarry Smith PetscFunctionBegin; 227827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 228b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 229827bd09bSSatish Balay 230827bd09bSSatish Balay /* non-uniform should have at least two entries */ 231b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?"); 232827bd09bSSatish Balay 233827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 2342fa5cd67SKarl Rupp if (!p_init) PCTFS_comm_init(); 235827bd09bSSatish Balay 236827bd09bSSatish Balay /* if there's nothing to do return */ 2372fa5cd67SKarl Rupp if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0); 238827bd09bSSatish Balay 239827bd09bSSatish Balay /* the error msg says it all!!! */ 240b1c944f5SJed Brown if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!"); 241827bd09bSSatish Balay 242827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 243b1c944f5SJed Brown if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n); 244827bd09bSSatish Balay 245827bd09bSSatish Balay /* can't do more dimensions then exist */ 246b1c944f5SJed Brown dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 247827bd09bSSatish Balay 248827bd09bSSatish Balay /* advance to list of n operations for custom */ 2492fa5cd67SKarl Rupp if ((type=oprs[0])==NON_UNIFORM) oprs++; 250827bd09bSSatish Balay 2512fa5cd67SKarl Rupp if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: Could not retrieve function pointer!\n"); 252827bd09bSSatish Balay 253db4deed7SKarl Rupp for (mask=1,edge=0; edge<dim; edge++,mask<<=1) { 254b1c944f5SJed Brown dest = PCTFS_my_id^mask; 2552fa5cd67SKarl Rupp if (PCTFS_my_id > dest) { 256ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 2572fa5cd67SKarl Rupp } else { 258ffc4695bSBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 259827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 260827bd09bSSatish Balay } 261827bd09bSSatish Balay } 262827bd09bSSatish Balay 2632fa5cd67SKarl Rupp if (edge==dim) mask>>=1; 264db4deed7SKarl Rupp else { 2652fa5cd67SKarl Rupp while (++edge<dim) mask<<=1; 266db4deed7SKarl Rupp } 267827bd09bSSatish Balay 268db4deed7SKarl Rupp for (edge=0; edge<dim; edge++,mask>>=1) { 2692fa5cd67SKarl Rupp if (PCTFS_my_id%mask) continue; 270827bd09bSSatish Balay 271b1c944f5SJed Brown dest = PCTFS_my_id^mask; 2722fa5cd67SKarl Rupp if (PCTFS_my_id < dest) { 273ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 2742fa5cd67SKarl Rupp } else { 275ffc4695bSBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 276827bd09bSSatish Balay } 277827bd09bSSatish Balay } 2783fdc5746SBarry Smith PetscFunctionReturn(0); 279827bd09bSSatish Balay } 280827bd09bSSatish Balay 2817b1ae94cSBarry Smith /******************************************************************************/ 282b1c944f5SJed Brown PetscErrorCode PCTFS_ssgl_radd(PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs) 283827bd09bSSatish Balay { 2843fdc5746SBarry Smith PetscInt edge, type, dest, mask; 2853fdc5746SBarry Smith PetscInt stage_n; 286827bd09bSSatish Balay MPI_Status status; 2873fdc5746SBarry Smith PetscErrorCode ierr; 288827bd09bSSatish Balay 2893fdc5746SBarry Smith PetscFunctionBegin; 290827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 2912fa5cd67SKarl Rupp if (!p_init) PCTFS_comm_init(); 292827bd09bSSatish Balay 293827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 294827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 295db4deed7SKarl Rupp for (mask=0, edge=0; edge<level; edge++, mask++) { 296827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 297db4deed7SKarl Rupp if (stage_n && !(PCTFS_my_id & mask)) { 298827bd09bSSatish Balay dest = edge_node[edge]; 299b1c944f5SJed Brown type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge); 3002fa5cd67SKarl Rupp if (PCTFS_my_id>dest) { 301ffc4695bSBarry Smith ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRMPI(ierr); 3022fa5cd67SKarl Rupp } else { 303b1c944f5SJed Brown type = type - PCTFS_my_id + dest; 304ffc4695bSBarry Smith ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 305ca8e9878SJed Brown PCTFS_rvec_add(vals+segs[edge], work, stage_n); 306827bd09bSSatish Balay } 307827bd09bSSatish Balay } 308827bd09bSSatish Balay mask <<= 1; 309827bd09bSSatish Balay } 310827bd09bSSatish Balay mask>>=1; 311db4deed7SKarl Rupp for (edge=0; edge<level; edge++) { 312827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 313db4deed7SKarl Rupp if (stage_n && !(PCTFS_my_id & mask)) { 314827bd09bSSatish Balay dest = edge_node[level-edge-1]; 315b1c944f5SJed Brown type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge); 3162fa5cd67SKarl Rupp if (PCTFS_my_id<dest) { 317ffc4695bSBarry Smith ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRMPI(ierr); 3182fa5cd67SKarl Rupp } else { 319b1c944f5SJed Brown type = type - PCTFS_my_id + dest; 320ffc4695bSBarry Smith ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 321827bd09bSSatish Balay } 322827bd09bSSatish Balay } 323827bd09bSSatish Balay mask >>= 1; 324827bd09bSSatish Balay } 3253fdc5746SBarry Smith PetscFunctionReturn(0); 326827bd09bSSatish Balay } 327827bd09bSSatish Balay 3287b1ae94cSBarry Smith /***********************************comm.c*************************************/ 329b1c944f5SJed Brown PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim) 330827bd09bSSatish Balay { 33152f87cdaSBarry Smith PetscInt mask, edge; 33252f87cdaSBarry Smith PetscInt type, dest; 333827bd09bSSatish Balay vfp fp; 334827bd09bSSatish Balay MPI_Status status; 3353fdc5746SBarry Smith PetscErrorCode ierr; 336827bd09bSSatish Balay 3373fdc5746SBarry Smith PetscFunctionBegin; 338827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 339b1c944f5SJed Brown if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 340827bd09bSSatish Balay 341827bd09bSSatish Balay /* non-uniform should have at least two entries */ 342b1c944f5SJed Brown if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?"); 343827bd09bSSatish Balay 344827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 3452fa5cd67SKarl Rupp if (!p_init) PCTFS_comm_init(); 346827bd09bSSatish Balay 347827bd09bSSatish Balay /* if there's nothing to do return */ 3482fa5cd67SKarl Rupp if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0); 349827bd09bSSatish Balay 350827bd09bSSatish Balay /* the error msg says it all!!! */ 351b1c944f5SJed Brown if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!"); 352827bd09bSSatish Balay 353827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 354b1c944f5SJed Brown if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n); 355827bd09bSSatish Balay 356827bd09bSSatish Balay /* can't do more dimensions then exist */ 357b1c944f5SJed Brown dim = PetscMin(dim,PCTFS_i_log2_num_nodes); 358827bd09bSSatish Balay 359827bd09bSSatish Balay /* advance to list of n operations for custom */ 3602fa5cd67SKarl Rupp if ((type=oprs[0])==NON_UNIFORM) oprs++; 361827bd09bSSatish Balay 3622fa5cd67SKarl Rupp if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: Could not retrieve function pointer!\n"); 363827bd09bSSatish Balay 364db4deed7SKarl Rupp for (mask=1,edge=0; edge<dim; edge++,mask<<=1) { 365b1c944f5SJed Brown dest = PCTFS_my_id^mask; 3662fa5cd67SKarl Rupp if (PCTFS_my_id > dest) { 367ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 3682fa5cd67SKarl Rupp } else { 369ffc4695bSBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRMPI(ierr); 370827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 371827bd09bSSatish Balay } 372827bd09bSSatish Balay } 373827bd09bSSatish Balay 3742fa5cd67SKarl Rupp if (edge==dim) mask>>=1; 3752fa5cd67SKarl Rupp else { 3762fa5cd67SKarl Rupp while (++edge<dim) mask<<=1; 3772fa5cd67SKarl Rupp } 378827bd09bSSatish Balay 379db4deed7SKarl Rupp for (edge=0; edge<dim; edge++,mask>>=1) { 3802fa5cd67SKarl Rupp if (PCTFS_my_id%mask) continue; 381827bd09bSSatish Balay 382b1c944f5SJed Brown dest = PCTFS_my_id^mask; 3832fa5cd67SKarl Rupp if (PCTFS_my_id < dest) { 384ffc4695bSBarry Smith ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRMPI(ierr); 3852fa5cd67SKarl Rupp } else { 386ffc4695bSBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRMPI(ierr); 387827bd09bSSatish Balay } 388827bd09bSSatish Balay } 3893fdc5746SBarry Smith PetscFunctionReturn(0); 390827bd09bSSatish Balay } 391