1 2 /* 3 Provides utility routines for split MPI communicator. 4 */ 5 #include <petscsys.h> /*I "petscsys.h" I*/ 6 #include <petscviewer.h> 7 8 const char *const PetscSubcommTypes[] = {"GENERAL","CONTIGUOUS","INTERLACED","PetscSubcommType","PETSC_SUBCOMM_",0}; 9 10 static PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm); 11 static PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm); 12 13 /*@C 14 PetscSubcommSetFromOptions - Allows setting options from a PetscSubcomm 15 16 Collective on PetscSubcomm 17 18 Input Parameter: 19 . psubcomm - PetscSubcomm context 20 21 Level: beginner 22 23 @*/ 24 PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm) 25 { 26 PetscErrorCode ierr; 27 PetscSubcommType type; 28 PetscBool flg; 29 30 PetscFunctionBegin; 31 if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt"); 32 33 ierr = PetscOptionsBegin(psubcomm->parent,psubcomm->subcommprefix,"Options for PetscSubcomm",NULL);CHKERRQ(ierr); 34 ierr = PetscOptionsEnum("-psubcomm_type",NULL,NULL,PetscSubcommTypes,(PetscEnum)psubcomm->type,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 35 if (flg && psubcomm->type != type) { 36 /* free old structures */ 37 ierr = PetscCommDestroy(&(psubcomm)->dupparent);CHKERRQ(ierr); 38 ierr = PetscCommDestroy(&(psubcomm)->child);CHKERRQ(ierr); 39 ierr = PetscFree((psubcomm)->subsize);CHKERRQ(ierr); 40 switch (type) { 41 case PETSC_SUBCOMM_GENERAL: 42 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Runtime option PETSC_SUBCOMM_GENERAL is not supported, use PetscSubcommSetTypeGeneral()"); 43 case PETSC_SUBCOMM_CONTIGUOUS: 44 ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 45 break; 46 case PETSC_SUBCOMM_INTERLACED: 47 ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 48 break; 49 default: 50 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[type]); 51 } 52 } 53 54 ierr = PetscOptionsName("-psubcomm_view","Triggers display of PetscSubcomm context","PetscSubcommView",&flg);CHKERRQ(ierr); 55 if (flg) { 56 ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_(psubcomm->parent));CHKERRQ(ierr); 57 } 58 ierr = PetscOptionsEnd();CHKERRQ(ierr); 59 PetscFunctionReturn(0); 60 } 61 62 /*@C 63 PetscSubcommSetOptionsPrefix - Sets the prefix used for searching for all 64 PetscSubcomm items in the options database. 65 66 Logically collective on PetscSubcomm. 67 68 Level: Intermediate 69 70 Input Parameters: 71 + psubcomm - PetscSubcomm context 72 - prefix - the prefix to prepend all PetscSubcomm item names with. 73 74 @*/ 75 PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm psubcomm,const char pre[]) 76 { 77 PetscErrorCode ierr; 78 79 PetscFunctionBegin; 80 if (!pre) { 81 ierr = PetscFree(psubcomm->subcommprefix);CHKERRQ(ierr); 82 } else { 83 if (pre[0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Options prefix should not begin with a hypen"); 84 ierr = PetscFree(psubcomm->subcommprefix);CHKERRQ(ierr); 85 ierr = PetscStrallocpy(pre,&(psubcomm->subcommprefix));CHKERRQ(ierr); 86 } 87 PetscFunctionReturn(0); 88 } 89 90 /*@C 91 PetscSubcommView - Views a PetscSubcomm of values as either ASCII text or a binary file 92 93 Collective on PetscSubcomm 94 95 Input Parameter: 96 + psubcomm - PetscSubcomm context 97 - viewer - location to view the values 98 99 Level: beginner 100 @*/ 101 PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer) 102 { 103 PetscErrorCode ierr; 104 PetscBool iascii; 105 PetscViewerFormat format; 106 107 PetscFunctionBegin; 108 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 109 if (iascii) { 110 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 111 if (format == PETSC_VIEWER_DEFAULT) { 112 MPI_Comm comm=psubcomm->parent; 113 PetscMPIInt rank,size,subsize,subrank,duprank; 114 115 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 116 ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %d MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr); 117 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 118 ierr = MPI_Comm_size(psubcomm->child,&subsize);CHKERRQ(ierr); 119 ierr = MPI_Comm_rank(psubcomm->child,&subrank);CHKERRQ(ierr); 120 ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr); 121 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 122 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n",rank,psubcomm->color,subsize,subrank,duprank);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 124 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 125 } 126 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet"); 127 PetscFunctionReturn(0); 128 } 129 130 /*@C 131 PetscSubcommSetNumber - Set total number of subcommunicators. 132 133 Collective on MPI_Comm 134 135 Input Parameter: 136 + psubcomm - PetscSubcomm context 137 - nsubcomm - the total number of subcommunicators in psubcomm 138 139 Level: advanced 140 141 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral() 142 @*/ 143 PetscErrorCode PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm) 144 { 145 PetscErrorCode ierr; 146 MPI_Comm comm=psubcomm->parent; 147 PetscMPIInt msub,size; 148 149 PetscFunctionBegin; 150 if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first"); 151 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 152 ierr = PetscMPIIntCast(nsubcomm,&msub);CHKERRQ(ierr); 153 if (msub < 1 || msub > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %d cannot be < 1 or > input comm size %d",msub,size); 154 155 psubcomm->n = msub; 156 PetscFunctionReturn(0); 157 } 158 159 /*@C 160 PetscSubcommSetType - Set type of subcommunicators. 161 162 Collective on MPI_Comm 163 164 Input Parameter: 165 + psubcomm - PetscSubcomm context 166 - subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED 167 168 Level: advanced 169 170 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 171 @*/ 172 PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype) 173 { 174 PetscErrorCode ierr; 175 176 PetscFunctionBegin; 177 if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 178 if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 179 180 if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) { 181 ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 182 } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) { 183 ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 184 } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[subcommtype]); 185 PetscFunctionReturn(0); 186 } 187 188 /*@C 189 PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications 190 191 Collective on MPI_Comm 192 193 Input Parameter: 194 + psubcomm - PetscSubcomm context 195 . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 196 - subrank - rank in the subcommunicator 197 198 Level: advanced 199 200 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 201 @*/ 202 PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank) 203 { 204 PetscErrorCode ierr; 205 MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 206 PetscMPIInt size,icolor,duprank,*recvbuf,sendbuf[3],mysubsize,rank,*subsize; 207 PetscMPIInt i,nsubcomm=psubcomm->n; 208 209 PetscFunctionBegin; 210 if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 211 if (nsubcomm < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",nsubcomm); 212 213 ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 214 215 /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 216 /* TODO: this can be done in an ostensibly scalale way (i.e., without allocating an array of size 'size') as is done in PetscObjectsCreateGlobalOrdering(). */ 217 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 218 ierr = PetscMalloc1(2*size,&recvbuf);CHKERRQ(ierr); 219 220 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 221 ierr = MPI_Comm_size(subcomm,&mysubsize);CHKERRQ(ierr); 222 223 sendbuf[0] = color; 224 sendbuf[1] = mysubsize; 225 ierr = MPI_Allgather(sendbuf,2,MPI_INT,recvbuf,2,MPI_INT,comm);CHKERRQ(ierr); 226 227 ierr = PetscCalloc1(nsubcomm,&subsize);CHKERRQ(ierr); 228 for (i=0; i<2*size; i+=2) { 229 subsize[recvbuf[i]] = recvbuf[i+1]; 230 } 231 ierr = PetscFree(recvbuf);CHKERRQ(ierr); 232 233 duprank = 0; 234 for (icolor=0; icolor<nsubcomm; icolor++) { 235 if (icolor != color) { /* not color of this process */ 236 duprank += subsize[icolor]; 237 } else { 238 duprank += subrank; 239 break; 240 } 241 } 242 ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 243 244 ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 245 ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 246 ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 247 ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 248 249 psubcomm->color = color; 250 psubcomm->subsize = subsize; 251 psubcomm->type = PETSC_SUBCOMM_GENERAL; 252 PetscFunctionReturn(0); 253 } 254 255 /*@C 256 PetscSubcommDestroy - Destroys a PetscSubcomm object 257 258 Collective on PetscSubcomm 259 260 Input Parameter: 261 . psubcomm - the PetscSubcomm context 262 263 Level: advanced 264 265 .seealso: PetscSubcommCreate(),PetscSubcommSetType() 266 @*/ 267 PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 268 { 269 PetscErrorCode ierr; 270 271 PetscFunctionBegin; 272 if (!*psubcomm) PetscFunctionReturn(0); 273 ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 274 ierr = PetscCommDestroy(&(*psubcomm)->child);CHKERRQ(ierr); 275 ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr); 276 if ((*psubcomm)->subcommprefix) { ierr = PetscFree((*psubcomm)->subcommprefix);CHKERRQ(ierr); } 277 ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 278 PetscFunctionReturn(0); 279 } 280 281 /*@C 282 PetscSubcommCreate - Create a PetscSubcomm context. 283 284 Collective on MPI_Comm 285 286 Input Parameter: 287 . comm - MPI communicator 288 289 Output Parameter: 290 . psubcomm - location to store the PetscSubcomm context 291 292 Level: advanced 293 294 .seealso: PetscSubcommDestroy() 295 @*/ 296 PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 297 { 298 PetscErrorCode ierr; 299 PetscMPIInt rank,size; 300 301 PetscFunctionBegin; 302 ierr = PetscNew(psubcomm);CHKERRQ(ierr); 303 304 /* set defaults */ 305 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 306 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 307 308 (*psubcomm)->parent = comm; 309 (*psubcomm)->dupparent = comm; 310 (*psubcomm)->child = PETSC_COMM_SELF; 311 (*psubcomm)->n = size; 312 (*psubcomm)->color = rank; 313 (*psubcomm)->subsize = NULL; 314 (*psubcomm)->type = PETSC_SUBCOMM_INTERLACED; 315 PetscFunctionReturn(0); 316 } 317 318 static PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 319 { 320 PetscErrorCode ierr; 321 PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 322 PetscMPIInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 323 MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 324 325 PetscFunctionBegin; 326 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 327 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 328 329 /* get size of each subcommunicator */ 330 ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr); 331 332 np_subcomm = size/nsubcomm; 333 nleftover = size - nsubcomm*np_subcomm; 334 for (i=0; i<nsubcomm; i++) { 335 subsize[i] = np_subcomm; 336 if (i<nleftover) subsize[i]++; 337 } 338 339 /* get color and subrank of this proc */ 340 rankstart = 0; 341 for (i=0; i<nsubcomm; i++) { 342 if (rank >= rankstart && rank < rankstart+subsize[i]) { 343 color = i; 344 subrank = rank - rankstart; 345 duprank = rank; 346 break; 347 } else rankstart += subsize[i]; 348 } 349 350 ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 351 352 /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 353 ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 354 ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 355 ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 356 ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 357 ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 358 359 psubcomm->color = color; 360 psubcomm->subsize = subsize; 361 psubcomm->type = PETSC_SUBCOMM_CONTIGUOUS; 362 PetscFunctionReturn(0); 363 } 364 365 /* 366 Note: 367 In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 368 by iteratively taking a process into a subcommunicator. 369 Example: size=4, nsubcomm=(*psubcomm)->n=3 370 comm=(*psubcomm)->parent: 371 rank: [0] [1] [2] [3] 372 color: 0 1 2 0 373 374 subcomm=(*psubcomm)->comm: 375 subrank: [0] [0] [0] [1] 376 377 dupcomm=(*psubcomm)->dupparent: 378 duprank: [0] [2] [3] [1] 379 380 Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 381 subcomm[color = 1] has subsize=1, owns process [1] 382 subcomm[color = 2] has subsize=1, owns process [2] 383 dupcomm has same number of processes as comm, and its duprank maps 384 processes in subcomm contiguously into a 1d array: 385 duprank: [0] [1] [2] [3] 386 rank: [0] [3] [1] [2] 387 subcomm[0] subcomm[1] subcomm[2] 388 */ 389 390 static PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 391 { 392 PetscErrorCode ierr; 393 PetscMPIInt rank,size,*subsize,duprank,subrank; 394 PetscMPIInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 395 MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 396 397 PetscFunctionBegin; 398 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 399 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 400 401 /* get size of each subcommunicator */ 402 ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr); 403 404 np_subcomm = size/nsubcomm; 405 nleftover = size - nsubcomm*np_subcomm; 406 for (i=0; i<nsubcomm; i++) { 407 subsize[i] = np_subcomm; 408 if (i<nleftover) subsize[i]++; 409 } 410 411 /* find color for this proc */ 412 color = rank%nsubcomm; 413 subrank = rank/nsubcomm; 414 415 ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 416 417 j = 0; duprank = 0; 418 for (i=0; i<nsubcomm; i++) { 419 if (j == color) { 420 duprank += subrank; 421 break; 422 } 423 duprank += subsize[i]; j++; 424 } 425 426 /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 427 ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 428 ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 429 ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 430 ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 431 ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 432 433 psubcomm->color = color; 434 psubcomm->subsize = subsize; 435 psubcomm->type = PETSC_SUBCOMM_INTERLACED; 436 PetscFunctionReturn(0); 437 } 438 439 440