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 extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm); 11 extern 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_WORLD);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 = PetscSynchronizedPrintf(comm," [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n",rank,psubcomm->color,subsize,subrank,duprank);CHKERRQ(ierr); 122 ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 123 } 124 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet"); 125 PetscFunctionReturn(0); 126 } 127 128 /*@C 129 PetscSubcommSetNumber - Set total number of subcommunicators. 130 131 Collective on MPI_Comm 132 133 Input Parameter: 134 + psubcomm - PetscSubcomm context 135 - nsubcomm - the total number of subcommunicators in psubcomm 136 137 Level: advanced 138 139 .keywords: communicator 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 .keywords: communicator 171 172 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 173 @*/ 174 PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype) 175 { 176 PetscErrorCode ierr; 177 178 PetscFunctionBegin; 179 if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 180 if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 181 182 if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) { 183 ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 184 } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) { 185 ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 186 } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[subcommtype]); 187 PetscFunctionReturn(0); 188 } 189 190 /*@C 191 PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications 192 193 Collective on MPI_Comm 194 195 Input Parameter: 196 + psubcomm - PetscSubcomm context 197 . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 198 - subrank - rank in the subcommunicator 199 200 Level: advanced 201 202 .keywords: communicator, create 203 204 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 205 @*/ 206 PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank) 207 { 208 PetscErrorCode ierr; 209 MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 210 PetscMPIInt size,icolor,duprank,*recvbuf,sendbuf[3],mysubsize,rank,*subsize; 211 PetscMPIInt i,nsubcomm=psubcomm->n; 212 213 PetscFunctionBegin; 214 if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 215 if (nsubcomm < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",nsubcomm); 216 217 ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 218 219 /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 220 /* TODO: this can be done in an ostensibly scalale way (i.e., without allocating an array of size 'size') as is done in PetscObjectsCreateGlobalOrdering(). */ 221 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 222 ierr = PetscMalloc1(2*size,&recvbuf);CHKERRQ(ierr); 223 224 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 225 ierr = MPI_Comm_size(subcomm,&mysubsize);CHKERRQ(ierr); 226 227 sendbuf[0] = color; 228 sendbuf[1] = mysubsize; 229 ierr = MPI_Allgather(sendbuf,2,MPI_INT,recvbuf,2,MPI_INT,comm);CHKERRQ(ierr); 230 231 ierr = PetscCalloc1(nsubcomm,&subsize);CHKERRQ(ierr); 232 for (i=0; i<2*size; i+=2) { 233 subsize[recvbuf[i]] = recvbuf[i+1]; 234 } 235 ierr = PetscFree(recvbuf);CHKERRQ(ierr); 236 237 duprank = 0; 238 for (icolor=0; icolor<nsubcomm; icolor++) { 239 if (icolor != color) { /* not color of this process */ 240 duprank += subsize[icolor]; 241 } else { 242 duprank += subrank; 243 break; 244 } 245 } 246 ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 247 248 ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 249 ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 250 ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 251 ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 252 253 psubcomm->color = color; 254 psubcomm->subsize = subsize; 255 psubcomm->type = PETSC_SUBCOMM_GENERAL; 256 PetscFunctionReturn(0); 257 } 258 259 /*@C 260 PetscSubcommDestroy - Destroys a PetscSubcomm object 261 262 Collective on PetscSubcomm 263 264 Input Parameter: 265 . psubcomm - the PetscSubcomm context 266 267 Level: advanced 268 269 .seealso: PetscSubcommCreate(),PetscSubcommSetType() 270 @*/ 271 PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 272 { 273 PetscErrorCode ierr; 274 275 PetscFunctionBegin; 276 if (!*psubcomm) PetscFunctionReturn(0); 277 ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 278 ierr = PetscCommDestroy(&(*psubcomm)->child);CHKERRQ(ierr); 279 ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr); 280 if ((*psubcomm)->subcommprefix) { ierr = PetscFree((*psubcomm)->subcommprefix);CHKERRQ(ierr); } 281 ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 /*@C 286 PetscSubcommCreate - Create a PetscSubcomm context. 287 288 Collective on MPI_Comm 289 290 Input Parameter: 291 . comm - MPI communicator 292 293 Output Parameter: 294 . psubcomm - location to store the PetscSubcomm context 295 296 Level: advanced 297 298 .keywords: communicator, create 299 300 .seealso: PetscSubcommDestroy() 301 @*/ 302 PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 303 { 304 PetscErrorCode ierr; 305 PetscMPIInt rank,size; 306 307 PetscFunctionBegin; 308 ierr = PetscNew(psubcomm);CHKERRQ(ierr); 309 310 /* set defaults */ 311 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 312 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 313 314 (*psubcomm)->parent = comm; 315 (*psubcomm)->dupparent = comm; 316 (*psubcomm)->child = PETSC_COMM_SELF; 317 (*psubcomm)->n = size; 318 (*psubcomm)->color = rank; 319 (*psubcomm)->subsize = NULL; 320 (*psubcomm)->type = PETSC_SUBCOMM_INTERLACED; 321 PetscFunctionReturn(0); 322 } 323 324 PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 325 { 326 PetscErrorCode ierr; 327 PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 328 PetscMPIInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 329 MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 330 331 PetscFunctionBegin; 332 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 333 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 334 335 /* get size of each subcommunicator */ 336 ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr); 337 338 np_subcomm = size/nsubcomm; 339 nleftover = size - nsubcomm*np_subcomm; 340 for (i=0; i<nsubcomm; i++) { 341 subsize[i] = np_subcomm; 342 if (i<nleftover) subsize[i]++; 343 } 344 345 /* get color and subrank of this proc */ 346 rankstart = 0; 347 for (i=0; i<nsubcomm; i++) { 348 if (rank >= rankstart && rank < rankstart+subsize[i]) { 349 color = i; 350 subrank = rank - rankstart; 351 duprank = rank; 352 break; 353 } else rankstart += subsize[i]; 354 } 355 356 ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 357 358 /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 359 ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 360 ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 361 ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 362 ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 363 ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 364 365 psubcomm->color = color; 366 psubcomm->subsize = subsize; 367 psubcomm->type = PETSC_SUBCOMM_CONTIGUOUS; 368 PetscFunctionReturn(0); 369 } 370 371 /* 372 Note: 373 In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 374 by iteratively taking a process into a subcommunicator. 375 Example: size=4, nsubcomm=(*psubcomm)->n=3 376 comm=(*psubcomm)->parent: 377 rank: [0] [1] [2] [3] 378 color: 0 1 2 0 379 380 subcomm=(*psubcomm)->comm: 381 subrank: [0] [0] [0] [1] 382 383 dupcomm=(*psubcomm)->dupparent: 384 duprank: [0] [2] [3] [1] 385 386 Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 387 subcomm[color = 1] has subsize=1, owns process [1] 388 subcomm[color = 2] has subsize=1, owns process [2] 389 dupcomm has same number of processes as comm, and its duprank maps 390 processes in subcomm contiguously into a 1d array: 391 duprank: [0] [1] [2] [3] 392 rank: [0] [3] [1] [2] 393 subcomm[0] subcomm[1] subcomm[2] 394 */ 395 396 PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 397 { 398 PetscErrorCode ierr; 399 PetscMPIInt rank,size,*subsize,duprank,subrank; 400 PetscMPIInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 401 MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 402 403 PetscFunctionBegin; 404 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 405 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 406 407 /* get size of each subcommunicator */ 408 ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr); 409 410 np_subcomm = size/nsubcomm; 411 nleftover = size - nsubcomm*np_subcomm; 412 for (i=0; i<nsubcomm; i++) { 413 subsize[i] = np_subcomm; 414 if (i<nleftover) subsize[i]++; 415 } 416 417 /* find color for this proc */ 418 color = rank%nsubcomm; 419 subrank = rank/nsubcomm; 420 421 ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 422 423 j = 0; duprank = 0; 424 for (i=0; i<nsubcomm; i++) { 425 if (j == color) { 426 duprank += subrank; 427 break; 428 } 429 duprank += subsize[i]; j++; 430 } 431 432 /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 433 ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 434 ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 435 ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 436 ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 437 ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 438 439 psubcomm->color = color; 440 psubcomm->subsize = subsize; 441 psubcomm->type = PETSC_SUBCOMM_INTERLACED; 442 PetscFunctionReturn(0); 443 } 444 445 446