1dba47a55SKris Buschelman #define PETSCKSP_DLL 2827bd09bSSatish Balay 3827bd09bSSatish Balay /***********************************comm.c************************************* 4827bd09bSSatish Balay 5827bd09bSSatish Balay Author: Henry M. Tufo III 6827bd09bSSatish Balay 7827bd09bSSatish Balay e-mail: hmt@cs.brown.edu 8827bd09bSSatish Balay 9827bd09bSSatish Balay snail-mail: 10827bd09bSSatish Balay Division of Applied Mathematics 11827bd09bSSatish Balay Brown University 12827bd09bSSatish Balay Providence, RI 02912 13827bd09bSSatish Balay 14827bd09bSSatish Balay Last Modification: 15827bd09bSSatish Balay 11.21.97 16827bd09bSSatish Balay ***********************************comm.c*************************************/ 177c4f633dSBarry Smith #include "../src/ksp/pc/impls/tfs/tfs.h" 18827bd09bSSatish Balay 19827bd09bSSatish Balay 20827bd09bSSatish Balay /* global program control variables - explicitly exported */ 2152f87cdaSBarry Smith PetscMPIInt my_id = 0; 2252f87cdaSBarry Smith PetscMPIInt num_nodes = 1; 236e4f4d19SBarry Smith PetscMPIInt floor_num_nodes = 0; 246e4f4d19SBarry Smith PetscMPIInt i_log2_num_nodes = 0; 25827bd09bSSatish Balay 26827bd09bSSatish Balay /* global program control variables */ 2752f87cdaSBarry Smith static PetscInt p_init = 0; 2852f87cdaSBarry Smith static PetscInt modfl_num_nodes; 2952f87cdaSBarry Smith static PetscInt edge_not_pow_2; 30827bd09bSSatish Balay 3152f87cdaSBarry Smith static PetscInt edge_node[sizeof(PetscInt)*32]; 32827bd09bSSatish Balay 337b1ae94cSBarry Smith /***********************************comm.c*************************************/ 340924e98cSBarry Smith PetscErrorCode comm_init (void) 35827bd09bSSatish Balay { 36827bd09bSSatish Balay 373fdc5746SBarry Smith if (p_init++) PetscFunctionReturn(0); 38827bd09bSSatish Balay 39827bd09bSSatish Balay MPI_Comm_size(MPI_COMM_WORLD,&num_nodes); 40827bd09bSSatish Balay MPI_Comm_rank(MPI_COMM_WORLD,&my_id); 41827bd09bSSatish Balay 42e7e72b3dSBarry Smith if (num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!"); 43827bd09bSSatish Balay 443fdc5746SBarry Smith ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32); 45827bd09bSSatish Balay 46827bd09bSSatish Balay floor_num_nodes = 1; 47827bd09bSSatish Balay i_log2_num_nodes = modfl_num_nodes = 0; 48827bd09bSSatish Balay while (floor_num_nodes <= num_nodes) 49827bd09bSSatish Balay { 50827bd09bSSatish Balay edge_node[i_log2_num_nodes] = my_id ^ floor_num_nodes; 51827bd09bSSatish Balay floor_num_nodes <<= 1; 52827bd09bSSatish Balay i_log2_num_nodes++; 53827bd09bSSatish Balay } 54827bd09bSSatish Balay 55827bd09bSSatish Balay i_log2_num_nodes--; 56827bd09bSSatish Balay floor_num_nodes >>= 1; 57827bd09bSSatish Balay modfl_num_nodes = (num_nodes - floor_num_nodes); 58827bd09bSSatish Balay 59827bd09bSSatish Balay if ((my_id > 0) && (my_id <= modfl_num_nodes)) 60827bd09bSSatish Balay {edge_not_pow_2=((my_id|floor_num_nodes)-1);} 61827bd09bSSatish Balay else if (my_id >= floor_num_nodes) 62827bd09bSSatish Balay {edge_not_pow_2=((my_id^floor_num_nodes)+1); 63827bd09bSSatish Balay } 64827bd09bSSatish Balay else 65827bd09bSSatish Balay {edge_not_pow_2 = 0;} 663fdc5746SBarry Smith PetscFunctionReturn(0); 67827bd09bSSatish Balay } 68827bd09bSSatish Balay 697b1ae94cSBarry Smith /***********************************comm.c*************************************/ 700924e98cSBarry Smith PetscErrorCode giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs) 71827bd09bSSatish Balay { 723fdc5746SBarry Smith PetscInt mask, edge; 733fdc5746SBarry Smith PetscInt type, dest; 74827bd09bSSatish Balay vfp fp; 75827bd09bSSatish Balay MPI_Status status; 763fdc5746SBarry Smith PetscInt ierr; 77827bd09bSSatish Balay 783fdc5746SBarry Smith PetscFunctionBegin; 79827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 80*c1235816SBarry Smith if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 81827bd09bSSatish Balay 82827bd09bSSatish Balay /* non-uniform should have at least two entries */ 83e7e72b3dSBarry Smith if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: non_uniform and n=0,1?"); 84827bd09bSSatish Balay 85827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 86827bd09bSSatish Balay if (!p_init) 87827bd09bSSatish Balay {comm_init();} 88827bd09bSSatish Balay 89827bd09bSSatish Balay /* if there's nothing to do return */ 90827bd09bSSatish Balay if ((num_nodes<2)||(!n)) 91827bd09bSSatish Balay { 923fdc5746SBarry Smith PetscFunctionReturn(0); 93827bd09bSSatish Balay } 94827bd09bSSatish Balay 9571a0148aSBarry Smith 96827bd09bSSatish Balay /* a negative number if items to send ==> fatal */ 97*c1235816SBarry Smith if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop() :: n=%D<0?",n); 98827bd09bSSatish Balay 99827bd09bSSatish Balay /* advance to list of n operations for custom */ 100827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 101827bd09bSSatish Balay {oprs++;} 102827bd09bSSatish Balay 103827bd09bSSatish Balay /* major league hack */ 104d890fc11SSatish Balay if (!(fp = (vfp) ivec_fct_addr(type))) { 105f1ed62a8SBarry Smith ierr = PetscInfo(0,"giop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 106827bd09bSSatish Balay fp = (vfp) oprs; 107827bd09bSSatish Balay } 108827bd09bSSatish Balay 109827bd09bSSatish Balay /* all msgs will be of the same length */ 110827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 111827bd09bSSatish Balay if (edge_not_pow_2) 112827bd09bSSatish Balay { 113827bd09bSSatish Balay if (my_id >= floor_num_nodes) 1143fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 115827bd09bSSatish Balay else 116827bd09bSSatish Balay { 1173fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr); 118827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 119827bd09bSSatish Balay } 120827bd09bSSatish Balay } 121827bd09bSSatish Balay 122827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 123827bd09bSSatish Balay if (my_id<floor_num_nodes) 124827bd09bSSatish Balay { 125827bd09bSSatish Balay for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1) 126827bd09bSSatish Balay { 127827bd09bSSatish Balay dest = my_id^mask; 128827bd09bSSatish Balay if (my_id > dest) 1293fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 130827bd09bSSatish Balay else 131827bd09bSSatish Balay { 1323fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 133827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 134827bd09bSSatish Balay } 135827bd09bSSatish Balay } 136827bd09bSSatish Balay 137827bd09bSSatish Balay mask=floor_num_nodes>>1; 138827bd09bSSatish Balay for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1) 139827bd09bSSatish Balay { 140827bd09bSSatish Balay if (my_id%mask) 141827bd09bSSatish Balay {continue;} 142827bd09bSSatish Balay 143827bd09bSSatish Balay dest = my_id^mask; 144827bd09bSSatish Balay if (my_id < dest) 1453fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 146827bd09bSSatish Balay else 147827bd09bSSatish Balay { 1483fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 149827bd09bSSatish Balay } 150827bd09bSSatish Balay } 151827bd09bSSatish Balay } 152827bd09bSSatish Balay 153827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 154827bd09bSSatish Balay if (edge_not_pow_2) 155827bd09bSSatish Balay { 156827bd09bSSatish Balay if (my_id >= floor_num_nodes) 157827bd09bSSatish Balay { 1583fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 159827bd09bSSatish Balay } 160827bd09bSSatish Balay else 1613fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 162827bd09bSSatish Balay } 1633fdc5746SBarry Smith PetscFunctionReturn(0); 164827bd09bSSatish Balay } 165827bd09bSSatish Balay 1667b1ae94cSBarry Smith /***********************************comm.c*************************************/ 1676e4f4d19SBarry Smith PetscErrorCode grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs) 168827bd09bSSatish Balay { 1693fdc5746SBarry Smith PetscInt mask, edge; 1703fdc5746SBarry Smith PetscInt type, dest; 171827bd09bSSatish Balay vfp fp; 172827bd09bSSatish Balay MPI_Status status; 1733fdc5746SBarry Smith PetscErrorCode ierr; 174827bd09bSSatish Balay 1753fdc5746SBarry Smith PetscFunctionBegin; 176827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 177*c1235816SBarry Smith if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 178827bd09bSSatish Balay 179827bd09bSSatish Balay /* non-uniform should have at least two entries */ 180e7e72b3dSBarry Smith if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop() :: non_uniform and n=0,1?"); 181827bd09bSSatish Balay 182827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 183827bd09bSSatish Balay if (!p_init) 184827bd09bSSatish Balay {comm_init();} 185827bd09bSSatish Balay 186827bd09bSSatish Balay /* if there's nothing to do return */ 187827bd09bSSatish Balay if ((num_nodes<2)||(!n)) 1883fdc5746SBarry Smith { PetscFunctionReturn(0);} 189827bd09bSSatish Balay 190827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 191*c1235816SBarry Smith if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n); 192827bd09bSSatish Balay 193827bd09bSSatish Balay /* advance to list of n operations for custom */ 194827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 195827bd09bSSatish Balay {oprs++;} 196827bd09bSSatish Balay 197d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))) { 198f1ed62a8SBarry Smith ierr = PetscInfo(0,"grop() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 199827bd09bSSatish Balay fp = (vfp) oprs; 200827bd09bSSatish Balay } 201827bd09bSSatish Balay 202827bd09bSSatish Balay /* all msgs will be of the same length */ 203827bd09bSSatish Balay /* if not a hypercube must colapse partial dim */ 204827bd09bSSatish Balay if (edge_not_pow_2) 205827bd09bSSatish Balay { 206827bd09bSSatish Balay if (my_id >= floor_num_nodes) 2073fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 208827bd09bSSatish Balay else 209827bd09bSSatish Balay { 2103fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 211827bd09bSSatish Balay (*fp)(vals,work,n,oprs); 212827bd09bSSatish Balay } 213827bd09bSSatish Balay } 214827bd09bSSatish Balay 215827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 216827bd09bSSatish Balay if (my_id<floor_num_nodes) 217827bd09bSSatish Balay { 218827bd09bSSatish Balay for (mask=1,edge=0; edge<i_log2_num_nodes; edge++,mask<<=1) 219827bd09bSSatish Balay { 220827bd09bSSatish Balay dest = my_id^mask; 221827bd09bSSatish Balay if (my_id > dest) 2223fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 223827bd09bSSatish Balay else 224827bd09bSSatish Balay { 2253fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 226827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 227827bd09bSSatish Balay } 228827bd09bSSatish Balay } 229827bd09bSSatish Balay 230827bd09bSSatish Balay mask=floor_num_nodes>>1; 231827bd09bSSatish Balay for (edge=0; edge<i_log2_num_nodes; edge++,mask>>=1) 232827bd09bSSatish Balay { 233827bd09bSSatish Balay if (my_id%mask) 234827bd09bSSatish Balay {continue;} 235827bd09bSSatish Balay 236827bd09bSSatish Balay dest = my_id^mask; 237827bd09bSSatish Balay if (my_id < dest) 2383fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 239827bd09bSSatish Balay else 240827bd09bSSatish Balay { 2413fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 242827bd09bSSatish Balay } 243827bd09bSSatish Balay } 244827bd09bSSatish Balay } 245827bd09bSSatish Balay 246827bd09bSSatish Balay /* if not a hypercube must expand to partial dim */ 247827bd09bSSatish Balay if (edge_not_pow_2) 248827bd09bSSatish Balay { 249827bd09bSSatish Balay if (my_id >= floor_num_nodes) 250827bd09bSSatish Balay { 2513fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr); 252827bd09bSSatish Balay } 253827bd09bSSatish Balay else 2543fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 255827bd09bSSatish Balay } 2563fdc5746SBarry Smith PetscFunctionReturn(0); 257827bd09bSSatish Balay } 258827bd09bSSatish Balay 2597b1ae94cSBarry Smith /***********************************comm.c*************************************/ 26052f87cdaSBarry Smith PetscErrorCode grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim) 261827bd09bSSatish Balay { 2623fdc5746SBarry Smith PetscInt mask, edge; 2633fdc5746SBarry Smith PetscInt type, dest; 264827bd09bSSatish Balay vfp fp; 265827bd09bSSatish Balay MPI_Status status; 2663fdc5746SBarry Smith PetscErrorCode ierr; 267827bd09bSSatish Balay 2683fdc5746SBarry Smith PetscFunctionBegin; 269827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 270*c1235816SBarry Smith if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 271827bd09bSSatish Balay 272827bd09bSSatish Balay /* non-uniform should have at least two entries */ 273e7e72b3dSBarry Smith if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: non_uniform and n=0,1?"); 274827bd09bSSatish Balay 275827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 276827bd09bSSatish Balay if (!p_init) 277827bd09bSSatish Balay {comm_init();} 278827bd09bSSatish Balay 279827bd09bSSatish Balay /* if there's nothing to do return */ 280827bd09bSSatish Balay if ((num_nodes<2)||(!n)||(dim<=0)) 2810924e98cSBarry Smith {PetscFunctionReturn(0);} 282827bd09bSSatish Balay 283827bd09bSSatish Balay /* the error msg says it all!!! */ 284e7e72b3dSBarry Smith if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: num_nodes not a power of 2!?!"); 285827bd09bSSatish Balay 286827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 287*c1235816SBarry Smith if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"grop_hc() :: n=%D<0?",n); 288827bd09bSSatish Balay 289827bd09bSSatish Balay /* can't do more dimensions then exist */ 29039945688SSatish Balay dim = PetscMin(dim,i_log2_num_nodes); 291827bd09bSSatish Balay 292827bd09bSSatish Balay /* advance to list of n operations for custom */ 293827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 294827bd09bSSatish Balay {oprs++;} 295827bd09bSSatish Balay 296d890fc11SSatish Balay if (!(fp = (vfp) rvec_fct_addr(type))) { 297f1ed62a8SBarry Smith ierr = PetscInfo(0,"grop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 298827bd09bSSatish Balay fp = (vfp) oprs; 299827bd09bSSatish Balay } 300827bd09bSSatish Balay 301827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 302827bd09bSSatish Balay { 303827bd09bSSatish Balay dest = my_id^mask; 304827bd09bSSatish Balay if (my_id > dest) 3053fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 306827bd09bSSatish Balay else 307827bd09bSSatish Balay { 3083fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 309827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 310827bd09bSSatish Balay } 311827bd09bSSatish Balay } 312827bd09bSSatish Balay 313827bd09bSSatish Balay if (edge==dim) 314827bd09bSSatish Balay {mask>>=1;} 315827bd09bSSatish Balay else 316827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 317827bd09bSSatish Balay 318827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 319827bd09bSSatish Balay { 320827bd09bSSatish Balay if (my_id%mask) 321827bd09bSSatish Balay {continue;} 322827bd09bSSatish Balay 323827bd09bSSatish Balay dest = my_id^mask; 324827bd09bSSatish Balay if (my_id < dest) 3253fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 326827bd09bSSatish Balay else 327827bd09bSSatish Balay { 3283fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 329827bd09bSSatish Balay } 330827bd09bSSatish Balay } 3313fdc5746SBarry Smith PetscFunctionReturn(0); 332827bd09bSSatish Balay } 333827bd09bSSatish Balay 3347b1ae94cSBarry Smith /******************************************************************************/ 3350924e98cSBarry Smith PetscErrorCode ssgl_radd( PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs) 336827bd09bSSatish Balay { 3373fdc5746SBarry Smith PetscInt edge, type, dest, mask; 3383fdc5746SBarry Smith PetscInt stage_n; 339827bd09bSSatish Balay MPI_Status status; 3403fdc5746SBarry Smith PetscErrorCode ierr; 341827bd09bSSatish Balay 3423fdc5746SBarry Smith PetscFunctionBegin; 343827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 344827bd09bSSatish Balay if (!p_init) 345827bd09bSSatish Balay {comm_init();} 346827bd09bSSatish Balay 347827bd09bSSatish Balay 348827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 349827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 350827bd09bSSatish Balay for (mask=0, edge=0; edge<level; edge++, mask++) 351827bd09bSSatish Balay { 352827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 353827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 354827bd09bSSatish Balay { 355827bd09bSSatish Balay dest = edge_node[edge]; 356827bd09bSSatish Balay type = MSGTAG3 + my_id + (num_nodes*edge); 357827bd09bSSatish Balay if (my_id>dest) 3583fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);} 359827bd09bSSatish Balay else 360827bd09bSSatish Balay { 361827bd09bSSatish Balay type = type - my_id + dest; 3623fdc5746SBarry Smith ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 363827bd09bSSatish Balay rvec_add(vals+segs[edge], work, stage_n); 364827bd09bSSatish Balay } 365827bd09bSSatish Balay } 366827bd09bSSatish Balay mask <<= 1; 367827bd09bSSatish Balay } 368827bd09bSSatish Balay mask>>=1; 369827bd09bSSatish Balay for (edge=0; edge<level; edge++) 370827bd09bSSatish Balay { 371827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 372827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 373827bd09bSSatish Balay { 374827bd09bSSatish Balay dest = edge_node[level-edge-1]; 375827bd09bSSatish Balay type = MSGTAG6 + my_id + (num_nodes*edge); 376827bd09bSSatish Balay if (my_id<dest) 3773fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);} 378827bd09bSSatish Balay else 379827bd09bSSatish Balay { 380827bd09bSSatish Balay type = type - my_id + dest; 3813fdc5746SBarry Smith ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 382827bd09bSSatish Balay } 383827bd09bSSatish Balay } 384827bd09bSSatish Balay mask >>= 1; 385827bd09bSSatish Balay } 3863fdc5746SBarry Smith PetscFunctionReturn(0); 387827bd09bSSatish Balay } 388827bd09bSSatish Balay 3897b1ae94cSBarry Smith /******************************************************************************/ 39052f87cdaSBarry Smith PetscErrorCode new_ssgl_radd( PetscScalar *vals, PetscScalar *work, PetscInt level, PetscInt *segs) 391827bd09bSSatish Balay { 39252f87cdaSBarry Smith PetscInt edge, type, dest, mask; 39352f87cdaSBarry Smith PetscInt stage_n; 394827bd09bSSatish Balay MPI_Status status; 3953fdc5746SBarry Smith PetscErrorCode ierr; 396827bd09bSSatish Balay 3973fdc5746SBarry Smith PetscFunctionBegin; 398827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 399827bd09bSSatish Balay if (!p_init) 400827bd09bSSatish Balay {comm_init();} 401827bd09bSSatish Balay 402827bd09bSSatish Balay /* all msgs are *NOT* the same length */ 403827bd09bSSatish Balay /* implement the mesh fan in/out exchange algorithm */ 404827bd09bSSatish Balay for (mask=0, edge=0; edge<level; edge++, mask++) 405827bd09bSSatish Balay { 406827bd09bSSatish Balay stage_n = (segs[level] - segs[edge]); 407827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 408827bd09bSSatish Balay { 409827bd09bSSatish Balay dest = edge_node[edge]; 410827bd09bSSatish Balay type = MSGTAG3 + my_id + (num_nodes*edge); 411827bd09bSSatish Balay if (my_id>dest) 4123fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);} 413827bd09bSSatish Balay else 414827bd09bSSatish Balay { 415827bd09bSSatish Balay type = type - my_id + dest; 4163fdc5746SBarry Smith ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type, MPI_COMM_WORLD,&status);CHKERRQ(ierr); 417827bd09bSSatish Balay rvec_add(vals+segs[edge], work, stage_n); 418827bd09bSSatish Balay } 419827bd09bSSatish Balay } 420827bd09bSSatish Balay mask <<= 1; 421827bd09bSSatish Balay } 422827bd09bSSatish Balay mask>>=1; 423827bd09bSSatish Balay for (edge=0; edge<level; edge++) 424827bd09bSSatish Balay { 425827bd09bSSatish Balay stage_n = (segs[level] - segs[level-1-edge]); 426827bd09bSSatish Balay if (stage_n && !(my_id & mask)) 427827bd09bSSatish Balay { 428827bd09bSSatish Balay dest = edge_node[level-edge-1]; 429827bd09bSSatish Balay type = MSGTAG6 + my_id + (num_nodes*edge); 430827bd09bSSatish Balay if (my_id<dest) 4313fdc5746SBarry Smith {ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);} 432827bd09bSSatish Balay else 433827bd09bSSatish Balay { 434827bd09bSSatish Balay type = type - my_id + dest; 4353fdc5746SBarry Smith ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 436827bd09bSSatish Balay } 437827bd09bSSatish Balay } 438827bd09bSSatish Balay mask >>= 1; 439827bd09bSSatish Balay } 4403fdc5746SBarry Smith PetscFunctionReturn(0); 441827bd09bSSatish Balay } 442827bd09bSSatish Balay 4437b1ae94cSBarry Smith /***********************************comm.c*************************************/ 44452f87cdaSBarry Smith PetscErrorCode giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim) 445827bd09bSSatish Balay { 44652f87cdaSBarry Smith PetscInt mask, edge; 44752f87cdaSBarry Smith PetscInt type, dest; 448827bd09bSSatish Balay vfp fp; 449827bd09bSSatish Balay MPI_Status status; 4503fdc5746SBarry Smith PetscErrorCode ierr; 451827bd09bSSatish Balay 4523fdc5746SBarry Smith PetscFunctionBegin; 453827bd09bSSatish Balay /* ok ... should have some data, work, and operator(s) */ 454*c1235816SBarry Smith if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs); 455827bd09bSSatish Balay 456827bd09bSSatish Balay /* non-uniform should have at least two entries */ 457e7e72b3dSBarry Smith if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: non_uniform and n=0,1?"); 458827bd09bSSatish Balay 459827bd09bSSatish Balay /* check to make sure comm package has been initialized */ 460827bd09bSSatish Balay if (!p_init) 461827bd09bSSatish Balay {comm_init();} 462827bd09bSSatish Balay 463827bd09bSSatish Balay /* if there's nothing to do return */ 464827bd09bSSatish Balay if ((num_nodes<2)||(!n)||(dim<=0)) 4653fdc5746SBarry Smith { PetscFunctionReturn(0);} 466827bd09bSSatish Balay 467827bd09bSSatish Balay /* the error msg says it all!!! */ 468e7e72b3dSBarry Smith if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: num_nodes not a power of 2!?!"); 469827bd09bSSatish Balay 470827bd09bSSatish Balay /* a negative number of items to send ==> fatal */ 471*c1235816SBarry Smith if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"giop_hc() :: n=%D<0?",n); 472827bd09bSSatish Balay 473827bd09bSSatish Balay /* can't do more dimensions then exist */ 47439945688SSatish Balay dim = PetscMin(dim,i_log2_num_nodes); 475827bd09bSSatish Balay 476827bd09bSSatish Balay /* advance to list of n operations for custom */ 477827bd09bSSatish Balay if ((type=oprs[0])==NON_UNIFORM) 478827bd09bSSatish Balay {oprs++;} 479827bd09bSSatish Balay 480d890fc11SSatish Balay if (!(fp = (vfp) ivec_fct_addr(type))){ 481f1ed62a8SBarry Smith ierr = PetscInfo(0,"giop_hc() :: hope you passed in a rbfp!\n");CHKERRQ(ierr); 482827bd09bSSatish Balay fp = (vfp) oprs; 483827bd09bSSatish Balay } 484827bd09bSSatish Balay 485827bd09bSSatish Balay for (mask=1,edge=0; edge<dim; edge++,mask<<=1) 486827bd09bSSatish Balay { 487827bd09bSSatish Balay dest = my_id^mask; 488827bd09bSSatish Balay if (my_id > dest) 4893fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 490827bd09bSSatish Balay else 491827bd09bSSatish Balay { 4923fdc5746SBarry Smith ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr); 493827bd09bSSatish Balay (*fp)(vals, work, n, oprs); 494827bd09bSSatish Balay } 495827bd09bSSatish Balay } 496827bd09bSSatish Balay 497827bd09bSSatish Balay if (edge==dim) 498827bd09bSSatish Balay {mask>>=1;} 499827bd09bSSatish Balay else 500827bd09bSSatish Balay {while (++edge<dim) {mask<<=1;}} 501827bd09bSSatish Balay 502827bd09bSSatish Balay for (edge=0; edge<dim; edge++,mask>>=1) 503827bd09bSSatish Balay { 504827bd09bSSatish Balay if (my_id%mask) 505827bd09bSSatish Balay {continue;} 506827bd09bSSatish Balay 507827bd09bSSatish Balay dest = my_id^mask; 508827bd09bSSatish Balay if (my_id < dest) 5093fdc5746SBarry Smith {ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+my_id,MPI_COMM_WORLD);CHKERRQ(ierr);} 510827bd09bSSatish Balay else 511827bd09bSSatish Balay { 5123fdc5746SBarry Smith ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr); 513827bd09bSSatish Balay } 514827bd09bSSatish Balay } 5153fdc5746SBarry Smith PetscFunctionReturn(0); 516827bd09bSSatish Balay } 517