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