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