1 2 /* 3 Provides the functions for index sets (IS) defined by a list of integers. 4 These are for blocks of data, each block is indicated with a single integer. 5 */ 6 #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 7 #include <petscviewer.h> 8 9 typedef struct { 10 PetscBool sorted; /* are the blocks sorted? */ 11 PetscBool allocated; /* did we allocate the index array ourselves? */ 12 PetscInt *idx; 13 } IS_Block; 14 15 static PetscErrorCode ISDestroy_Block(IS is) 16 { 17 IS_Block *sub = (IS_Block*)is->data; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 if (sub->allocated) {ierr = PetscFree(sub->idx);CHKERRQ(ierr);} 22 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",NULL);CHKERRQ(ierr); 23 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",NULL);CHKERRQ(ierr); 24 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",NULL);CHKERRQ(ierr); 25 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",NULL);CHKERRQ(ierr); 26 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",NULL);CHKERRQ(ierr); 27 ierr = PetscFree(is->data);CHKERRQ(ierr); 28 PetscFunctionReturn(0); 29 } 30 31 static PetscErrorCode ISLocate_Block(IS is,PetscInt key,PetscInt *location) 32 { 33 IS_Block *sub = (IS_Block*)is->data; 34 PetscInt numIdx, i, bs, bkey, mkey; 35 PetscBool sorted; 36 PetscErrorCode ierr; 37 38 PetscFunctionBegin; 39 ierr = PetscLayoutGetBlockSize(is->map,&bs);CHKERRQ(ierr); 40 ierr = PetscLayoutGetSize(is->map,&numIdx);CHKERRQ(ierr); 41 numIdx /= bs; 42 bkey = key / bs; 43 mkey = key % bs; 44 if (mkey < 0) { 45 bkey--; 46 mkey += bs; 47 } 48 ierr = ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted);CHKERRQ(ierr); 49 if (sorted) { 50 ierr = PetscFindInt(bkey,numIdx,sub->idx,location);CHKERRQ(ierr); 51 } else { 52 const PetscInt *idx = sub->idx; 53 54 *location = -1; 55 for (i = 0; i < numIdx; i++) { 56 if (idx[i] == bkey) { 57 *location = i; 58 break; 59 } 60 } 61 } 62 if (*location >= 0) { 63 *location = *location * bs + mkey; 64 } 65 PetscFunctionReturn(0); 66 } 67 68 static PetscErrorCode ISGetIndices_Block(IS in,const PetscInt *idx[]) 69 { 70 IS_Block *sub = (IS_Block*)in->data; 71 PetscErrorCode ierr; 72 PetscInt i,j,k,bs,n,*ii,*jj; 73 74 PetscFunctionBegin; 75 ierr = PetscLayoutGetBlockSize(in->map, &bs);CHKERRQ(ierr); 76 ierr = PetscLayoutGetLocalSize(in->map, &n);CHKERRQ(ierr); 77 n /= bs; 78 if (bs == 1) *idx = sub->idx; 79 else { 80 ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr); 81 *idx = jj; 82 k = 0; 83 ii = sub->idx; 84 for (i=0; i<n; i++) 85 for (j=0; j<bs; j++) 86 jj[k++] = bs*ii[i] + j; 87 } 88 PetscFunctionReturn(0); 89 } 90 91 static PetscErrorCode ISRestoreIndices_Block(IS is,const PetscInt *idx[]) 92 { 93 IS_Block *sub = (IS_Block*)is->data; 94 PetscInt bs; 95 PetscErrorCode ierr; 96 97 PetscFunctionBegin; 98 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 99 if (bs != 1) { 100 ierr = PetscFree(*(void**)idx);CHKERRQ(ierr); 101 } else { 102 if (*idx != sub->idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must restore with value from ISGetIndices()"); 103 } 104 PetscFunctionReturn(0); 105 } 106 107 static PetscErrorCode ISInvertPermutation_Block(IS is,PetscInt nlocal,IS *isout) 108 { 109 IS_Block *sub = (IS_Block*)is->data; 110 PetscInt i,*ii,bs,n,*idx = sub->idx; 111 PetscMPIInt size; 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);CHKERRQ(ierr); 116 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 117 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 118 n /= bs; 119 if (size == 1) { 120 ierr = PetscMalloc1(n,&ii);CHKERRQ(ierr); 121 for (i=0; i<n; i++) ii[idx[i]] = i; 122 ierr = ISCreateBlock(PETSC_COMM_SELF,bs,n,ii,PETSC_OWN_POINTER,isout);CHKERRQ(ierr); 123 ierr = ISSetPermutation(*isout);CHKERRQ(ierr); 124 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No inversion written yet for block IS"); 125 PetscFunctionReturn(0); 126 } 127 128 static PetscErrorCode ISView_Block(IS is, PetscViewer viewer) 129 { 130 IS_Block *sub = (IS_Block*)is->data; 131 PetscErrorCode ierr; 132 PetscInt i,bs,n,*idx = sub->idx; 133 PetscBool iascii; 134 135 PetscFunctionBegin; 136 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 137 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 138 n /= bs; 139 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 140 if (iascii) { 141 PetscViewerFormat fmt; 142 143 ierr = PetscViewerGetFormat(viewer,&fmt);CHKERRQ(ierr); 144 if (fmt == PETSC_VIEWER_ASCII_MATLAB) { 145 IS ist; 146 const char *name; 147 const PetscInt *idx; 148 PetscInt n; 149 150 ierr = PetscObjectGetName((PetscObject)is,&name);CHKERRQ(ierr); 151 ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 152 ierr = ISGetIndices(is,&idx);CHKERRQ(ierr); 153 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idx,PETSC_USE_POINTER,&ist);CHKERRQ(ierr); 154 ierr = PetscObjectSetName((PetscObject)ist,name);CHKERRQ(ierr); 155 ierr = ISView(ist,viewer);CHKERRQ(ierr); 156 ierr = ISDestroy(&ist);CHKERRQ(ierr); 157 ierr = ISRestoreIndices(is,&idx);CHKERRQ(ierr); 158 } else { 159 PetscBool isperm; 160 161 ierr = ISGetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_FALSE,&isperm);CHKERRQ(ierr); 162 if (isperm) {ierr = PetscViewerASCIIPrintf(viewer,"Block Index set is permutation\n");CHKERRQ(ierr);} 163 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 164 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Block size %D\n",bs);CHKERRQ(ierr); 165 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of block indices in set %D\n",n);CHKERRQ(ierr); 166 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"The first indices of each block are\n");CHKERRQ(ierr); 167 for (i=0; i<n; i++) { 168 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Block %D Index %D\n",i,idx[i]);CHKERRQ(ierr); 169 } 170 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 171 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 172 } 173 } 174 PetscFunctionReturn(0); 175 } 176 177 static PetscErrorCode ISSort_Block(IS is) 178 { 179 IS_Block *sub = (IS_Block*)is->data; 180 PetscInt bs, n; 181 PetscErrorCode ierr; 182 183 PetscFunctionBegin; 184 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 185 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 186 ierr = PetscSortInt(n/bs,sub->idx);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 190 static PetscErrorCode ISSortRemoveDups_Block(IS is) 191 { 192 IS_Block *sub = (IS_Block*)is->data; 193 PetscInt bs, n, nb; 194 PetscBool sorted; 195 PetscErrorCode ierr; 196 197 PetscFunctionBegin; 198 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 199 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 200 nb = n/bs; 201 ierr = ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted);CHKERRQ(ierr); 202 if (sorted) { 203 ierr = PetscSortedRemoveDupsInt(&nb,sub->idx);CHKERRQ(ierr); 204 } else { 205 ierr = PetscSortRemoveDupsInt(&nb,sub->idx);CHKERRQ(ierr); 206 } 207 ierr = PetscLayoutDestroy(&is->map);CHKERRQ(ierr); 208 ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), nb*bs, PETSC_DECIDE, bs, &is->map);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 static PetscErrorCode ISSorted_Block(IS is,PetscBool *flg) 213 { 214 PetscErrorCode ierr; 215 216 PetscFunctionBegin; 217 ierr = ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,flg);CHKERRQ(ierr); 218 PetscFunctionReturn(0); 219 } 220 221 static PetscErrorCode ISSortedLocal_Block(IS is,PetscBool *flg) 222 { 223 IS_Block *sub = (IS_Block*)is->data; 224 PetscInt n, bs, i, *idx; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 229 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 230 n /= bs; 231 idx = sub->idx; 232 for (i = 1; i < n; i++) if (idx[i] < idx[i - 1]) break; 233 if (i < n) *flg = PETSC_FALSE; 234 else *flg = PETSC_TRUE; 235 PetscFunctionReturn(0); 236 } 237 238 static PetscErrorCode ISUniqueLocal_Block(IS is,PetscBool *flg) 239 { 240 IS_Block *sub = (IS_Block*)is->data; 241 PetscInt n, bs, i, *idx, *idxcopy = NULL; 242 PetscBool sortedLocal; 243 PetscErrorCode ierr; 244 245 PetscFunctionBegin; 246 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 247 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 248 n /= bs; 249 idx = sub->idx; 250 ierr = ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal);CHKERRQ(ierr); 251 if (!sortedLocal) { 252 ierr = PetscMalloc1(n, &idxcopy);CHKERRQ(ierr); 253 ierr = PetscArraycpy(idxcopy, idx, n);CHKERRQ(ierr); 254 ierr = PetscSortInt(n, idxcopy);CHKERRQ(ierr); 255 idx = idxcopy; 256 } 257 for (i = 1; i < n; i++) if (idx[i] == idx[i - 1]) break; 258 if (i < n) *flg = PETSC_FALSE; 259 else *flg = PETSC_TRUE; 260 ierr = PetscFree(idxcopy);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 static PetscErrorCode ISPermutationLocal_Block(IS is,PetscBool *flg) 265 { 266 IS_Block *sub = (IS_Block*)is->data; 267 PetscInt n, bs, i, *idx, *idxcopy = NULL; 268 PetscBool sortedLocal; 269 PetscErrorCode ierr; 270 271 PetscFunctionBegin; 272 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 273 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 274 n /= bs; 275 idx = sub->idx; 276 ierr = ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal);CHKERRQ(ierr); 277 if (!sortedLocal) { 278 ierr = PetscMalloc1(n, &idxcopy);CHKERRQ(ierr); 279 ierr = PetscArraycpy(idxcopy, idx, n);CHKERRQ(ierr); 280 ierr = PetscSortInt(n, idxcopy);CHKERRQ(ierr); 281 idx = idxcopy; 282 } 283 for (i = 0; i < n; i++) if (idx[i] != i) break; 284 if (i < n) *flg = PETSC_FALSE; 285 else *flg = PETSC_TRUE; 286 ierr = PetscFree(idxcopy);CHKERRQ(ierr); 287 PetscFunctionReturn(0); 288 } 289 290 static PetscErrorCode ISIntervalLocal_Block(IS is,PetscBool *flg) 291 { 292 IS_Block *sub = (IS_Block*)is->data; 293 PetscInt n, bs, i, *idx; 294 PetscErrorCode ierr; 295 296 PetscFunctionBegin; 297 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 298 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 299 n /= bs; 300 idx = sub->idx; 301 for (i = 1; i < n; i++) if (idx[i] != idx[i - 1] + 1) break; 302 if (i < n) *flg = PETSC_FALSE; 303 else *flg = PETSC_TRUE; 304 PetscFunctionReturn(0); 305 } 306 307 static PetscErrorCode ISDuplicate_Block(IS is,IS *newIS) 308 { 309 PetscErrorCode ierr; 310 IS_Block *sub = (IS_Block*)is->data; 311 PetscInt bs, n; 312 313 PetscFunctionBegin; 314 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 315 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 316 n /= bs; 317 ierr = ISCreateBlock(PetscObjectComm((PetscObject)is),bs,n,sub->idx,PETSC_COPY_VALUES,newIS);CHKERRQ(ierr); 318 PetscFunctionReturn(0); 319 } 320 321 static PetscErrorCode ISCopy_Block(IS is,IS isy) 322 { 323 IS_Block *is_block = (IS_Block*)is->data,*isy_block = (IS_Block*)isy->data; 324 PetscInt bs, n, N, bsy, ny, Ny; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 329 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 330 ierr = PetscLayoutGetSize(is->map, &N);CHKERRQ(ierr); 331 ierr = PetscLayoutGetBlockSize(isy->map, &bsy);CHKERRQ(ierr); 332 ierr = PetscLayoutGetLocalSize(isy->map, &ny);CHKERRQ(ierr); 333 ierr = PetscLayoutGetSize(isy->map, &Ny);CHKERRQ(ierr); 334 if (n != ny || N != Ny || bs != bsy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Index sets incompatible"); 335 ierr = PetscArraycpy(isy_block->idx,is_block->idx,(n/bs));CHKERRQ(ierr); 336 PetscFunctionReturn(0); 337 } 338 339 static PetscErrorCode ISOnComm_Block(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis) 340 { 341 PetscErrorCode ierr; 342 IS_Block *sub = (IS_Block*)is->data; 343 PetscInt bs, n; 344 345 PetscFunctionBegin; 346 if (mode == PETSC_OWN_POINTER) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Cannot use PETSC_OWN_POINTER"); 347 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 348 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 349 ierr = ISCreateBlock(comm,bs,n/bs,sub->idx,mode,newis);CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352 353 static PetscErrorCode ISSetBlockSize_Block(IS is,PetscInt bs) 354 { 355 PetscErrorCode ierr; 356 357 PetscFunctionBegin; 358 ierr = PetscLayoutSetBlockSize(is->map, bs);CHKERRQ(ierr); 359 PetscFunctionReturn(0); 360 } 361 362 static PetscErrorCode ISToGeneral_Block(IS inis) 363 { 364 IS_Block *sub = (IS_Block*)inis->data; 365 PetscInt bs,n; 366 const PetscInt *idx; 367 PetscErrorCode ierr; 368 369 PetscFunctionBegin; 370 ierr = ISGetBlockSize(inis,&bs);CHKERRQ(ierr); 371 ierr = ISGetLocalSize(inis,&n);CHKERRQ(ierr); 372 ierr = ISGetIndices(inis,&idx);CHKERRQ(ierr); 373 if (bs == 1) { 374 PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER; 375 sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/ 376 ierr = ISSetType(inis,ISGENERAL);CHKERRQ(ierr); 377 ierr = ISGeneralSetIndices(inis,n,idx,mode);CHKERRQ(ierr); 378 } else { 379 ierr = ISSetType(inis,ISGENERAL);CHKERRQ(ierr); 380 ierr = ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);CHKERRQ(ierr); 381 } 382 PetscFunctionReturn(0); 383 } 384 385 386 static struct _ISOps myops = { ISGetIndices_Block, 387 ISRestoreIndices_Block, 388 ISInvertPermutation_Block, 389 ISSort_Block, 390 ISSortRemoveDups_Block, 391 ISSorted_Block, 392 ISDuplicate_Block, 393 ISDestroy_Block, 394 ISView_Block, 395 ISLoad_Default, 396 ISCopy_Block, 397 ISToGeneral_Block, 398 ISOnComm_Block, 399 ISSetBlockSize_Block, 400 0, 401 ISLocate_Block, 402 /* we can have specialized local routines for determining properties, 403 * but unless the block size is the same on each process (which is not guaranteed at 404 * the moment), then trying to do something specialized for global properties is too 405 * complicated */ 406 ISSortedLocal_Block, 407 NULL, 408 ISUniqueLocal_Block, 409 NULL, 410 ISPermutationLocal_Block, 411 NULL, 412 ISIntervalLocal_Block, 413 NULL}; 414 415 /*@ 416 ISBlockSetIndices - The indices are relative to entries, not blocks. 417 418 Collective on IS 419 420 Input Parameters: 421 + is - the index set 422 . bs - number of elements in each block, one for each block and count of block not indices 423 . n - the length of the index set (the number of blocks) 424 . idx - the list of integers, these are by block, not by location 425 - mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported 426 427 428 Notes: 429 When the communicator is not MPI_COMM_SELF, the operations on the 430 index sets, IS, are NOT conceptually the same as MPI_Group operations. 431 The index sets are then distributed sets of indices and thus certain operations 432 on them are collective. 433 434 Example: 435 If you wish to index the values {0,1,4,5}, then use 436 a block size of 2 and idx of {0,2}. 437 438 Level: beginner 439 440 .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather() 441 @*/ 442 PetscErrorCode ISBlockSetIndices(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode) 443 { 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 ierr = ISClearInfoCache(is,PETSC_FALSE);CHKERRQ(ierr); 448 ierr = PetscUseMethod(is,"ISBlockSetIndices_C",(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode),(is,bs,n,idx,mode));CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 452 static PetscErrorCode ISBlockSetIndices_Block(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode) 453 { 454 PetscErrorCode ierr; 455 PetscInt i,min,max; 456 IS_Block *sub = (IS_Block*)is->data; 457 PetscLayout map; 458 459 PetscFunctionBegin; 460 if (bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1"); 461 if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0"); 462 if (n) PetscValidIntPointer(idx,3); 463 464 ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n*bs,is->map->N,bs,&map);CHKERRQ(ierr); 465 ierr = PetscLayoutDestroy(&is->map);CHKERRQ(ierr); 466 is->map = map; 467 468 if (sub->allocated) {ierr = PetscFree(sub->idx);CHKERRQ(ierr);} 469 if (mode == PETSC_COPY_VALUES) { 470 ierr = PetscMalloc1(n,&sub->idx);CHKERRQ(ierr); 471 ierr = PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));CHKERRQ(ierr); 472 ierr = PetscArraycpy(sub->idx,idx,n);CHKERRQ(ierr); 473 sub->allocated = PETSC_TRUE; 474 } else if (mode == PETSC_OWN_POINTER) { 475 sub->idx = (PetscInt*) idx; 476 ierr = PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));CHKERRQ(ierr); 477 sub->allocated = PETSC_TRUE; 478 } else if (mode == PETSC_USE_POINTER) { 479 sub->idx = (PetscInt*) idx; 480 sub->allocated = PETSC_FALSE; 481 } 482 483 if (n) { 484 min = max = idx[0]; 485 for (i=1; i<n; i++) { 486 if (idx[i] < min) min = idx[i]; 487 if (idx[i] > max) max = idx[i]; 488 } 489 is->min = bs*min; 490 is->max = bs*max+bs-1; 491 } else { 492 is->min = PETSC_MAX_INT; 493 is->max = PETSC_MIN_INT; 494 } 495 PetscFunctionReturn(0); 496 } 497 498 /*@ 499 ISCreateBlock - Creates a data structure for an index set containing 500 a list of integers. The indices are relative to entries, not blocks. 501 502 Collective 503 504 Input Parameters: 505 + comm - the MPI communicator 506 . bs - number of elements in each block 507 . n - the length of the index set (the number of blocks) 508 . idx - the list of integers, one for each block and count of block not indices 509 - mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine 510 511 Output Parameter: 512 . is - the new index set 513 514 Notes: 515 When the communicator is not MPI_COMM_SELF, the operations on the 516 index sets, IS, are NOT conceptually the same as MPI_Group operations. 517 The index sets are then distributed sets of indices and thus certain operations 518 on them are collective. 519 520 Example: 521 If you wish to index the values {0,1,6,7}, then use 522 a block size of 2 and idx of {0,3}. 523 524 Level: beginner 525 526 .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather() 527 @*/ 528 PetscErrorCode ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode,IS *is) 529 { 530 PetscErrorCode ierr; 531 532 PetscFunctionBegin; 533 PetscValidPointer(is,5); 534 if (bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1"); 535 if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0"); 536 if (n) PetscValidIntPointer(idx,4); 537 538 ierr = ISCreate(comm,is);CHKERRQ(ierr); 539 ierr = ISSetType(*is,ISBLOCK);CHKERRQ(ierr); 540 ierr = ISBlockSetIndices(*is,bs,n,idx,mode);CHKERRQ(ierr); 541 PetscFunctionReturn(0); 542 } 543 544 static PetscErrorCode ISBlockGetIndices_Block(IS is,const PetscInt *idx[]) 545 { 546 IS_Block *sub = (IS_Block*)is->data; 547 548 PetscFunctionBegin; 549 *idx = sub->idx; 550 PetscFunctionReturn(0); 551 } 552 553 static PetscErrorCode ISBlockRestoreIndices_Block(IS is,const PetscInt *idx[]) 554 { 555 PetscFunctionBegin; 556 PetscFunctionReturn(0); 557 } 558 559 /*@C 560 ISBlockGetIndices - Gets the indices associated with each block. 561 562 Not Collective 563 564 Input Parameter: 565 . is - the index set 566 567 Output Parameter: 568 . idx - the integer indices, one for each block and count of block not indices 569 570 Level: intermediate 571 572 .seealso: ISGetIndices(), ISBlockRestoreIndices() 573 @*/ 574 PetscErrorCode ISBlockGetIndices(IS is,const PetscInt *idx[]) 575 { 576 PetscErrorCode ierr; 577 578 PetscFunctionBegin; 579 ierr = PetscUseMethod(is,"ISBlockGetIndices_C",(IS,const PetscInt*[]),(is,idx));CHKERRQ(ierr); 580 PetscFunctionReturn(0); 581 } 582 583 /*@C 584 ISBlockRestoreIndices - Restores the indices associated with each block. 585 586 Not Collective 587 588 Input Parameter: 589 . is - the index set 590 591 Output Parameter: 592 . idx - the integer indices 593 594 Level: intermediate 595 596 .seealso: ISRestoreIndices(), ISBlockGetIndices() 597 @*/ 598 PetscErrorCode ISBlockRestoreIndices(IS is,const PetscInt *idx[]) 599 { 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 ierr = PetscUseMethod(is,"ISBlockRestoreIndices_C",(IS,const PetscInt*[]),(is,idx));CHKERRQ(ierr); 604 PetscFunctionReturn(0); 605 } 606 607 /*@ 608 ISBlockGetLocalSize - Returns the local number of blocks in the index set. 609 610 Not Collective 611 612 Input Parameter: 613 . is - the index set 614 615 Output Parameter: 616 . size - the local number of blocks 617 618 Level: intermediate 619 620 621 .seealso: ISGetBlockSize(), ISBlockGetSize(), ISGetSize(), ISCreateBlock() 622 @*/ 623 PetscErrorCode ISBlockGetLocalSize(IS is,PetscInt *size) 624 { 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 ierr = PetscUseMethod(is,"ISBlockGetLocalSize_C",(IS,PetscInt*),(is,size));CHKERRQ(ierr); 629 PetscFunctionReturn(0); 630 } 631 632 static PetscErrorCode ISBlockGetLocalSize_Block(IS is,PetscInt *size) 633 { 634 PetscInt bs, n; 635 PetscErrorCode ierr; 636 637 PetscFunctionBegin; 638 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 639 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 640 *size = n/bs; 641 PetscFunctionReturn(0); 642 } 643 644 /*@ 645 ISBlockGetSize - Returns the global number of blocks in the index set. 646 647 Not Collective 648 649 Input Parameter: 650 . is - the index set 651 652 Output Parameter: 653 . size - the global number of blocks 654 655 Level: intermediate 656 657 658 .seealso: ISGetBlockSize(), ISBlockGetLocalSize(), ISGetSize(), ISCreateBlock() 659 @*/ 660 PetscErrorCode ISBlockGetSize(IS is,PetscInt *size) 661 { 662 PetscErrorCode ierr; 663 664 PetscFunctionBegin; 665 ierr = PetscUseMethod(is,"ISBlockGetSize_C",(IS,PetscInt*),(is,size));CHKERRQ(ierr); 666 PetscFunctionReturn(0); 667 } 668 669 static PetscErrorCode ISBlockGetSize_Block(IS is,PetscInt *size) 670 { 671 PetscInt bs, N; 672 PetscErrorCode ierr; 673 674 PetscFunctionBegin; 675 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 676 ierr = PetscLayoutGetSize(is->map, &N);CHKERRQ(ierr); 677 *size = N/bs; 678 PetscFunctionReturn(0); 679 } 680 681 PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is) 682 { 683 PetscErrorCode ierr; 684 IS_Block *sub; 685 686 PetscFunctionBegin; 687 ierr = PetscNewLog(is,&sub);CHKERRQ(ierr); 688 is->data = (void *) sub; 689 ierr = PetscMemcpy(is->ops,&myops,sizeof(myops));CHKERRQ(ierr); 690 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",ISBlockSetIndices_Block);CHKERRQ(ierr); 691 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",ISBlockGetIndices_Block);CHKERRQ(ierr); 692 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",ISBlockRestoreIndices_Block);CHKERRQ(ierr); 693 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",ISBlockGetSize_Block);CHKERRQ(ierr); 694 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",ISBlockGetLocalSize_Block);CHKERRQ(ierr); 695 PetscFunctionReturn(0); 696 } 697