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 if (is->map->bs > 0 && bs != is->map->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot change blocksize %D (to %D) if ISType is ISBLOCK",is->map->bs,bs); 359 ierr = PetscLayoutSetBlockSize(is->map, bs);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 static PetscErrorCode ISToGeneral_Block(IS inis) 364 { 365 IS_Block *sub = (IS_Block*)inis->data; 366 PetscInt bs,n; 367 const PetscInt *idx; 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 ierr = ISGetBlockSize(inis,&bs);CHKERRQ(ierr); 372 ierr = ISGetLocalSize(inis,&n);CHKERRQ(ierr); 373 ierr = ISGetIndices(inis,&idx);CHKERRQ(ierr); 374 if (bs == 1) { 375 PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER; 376 sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/ 377 ierr = ISSetType(inis,ISGENERAL);CHKERRQ(ierr); 378 ierr = ISGeneralSetIndices(inis,n,idx,mode);CHKERRQ(ierr); 379 } else { 380 ierr = ISSetType(inis,ISGENERAL);CHKERRQ(ierr); 381 ierr = ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);CHKERRQ(ierr); 382 } 383 PetscFunctionReturn(0); 384 } 385 386 387 static struct _ISOps myops = { ISGetIndices_Block, 388 ISRestoreIndices_Block, 389 ISInvertPermutation_Block, 390 ISSort_Block, 391 ISSortRemoveDups_Block, 392 ISSorted_Block, 393 ISDuplicate_Block, 394 ISDestroy_Block, 395 ISView_Block, 396 ISLoad_Default, 397 ISCopy_Block, 398 ISToGeneral_Block, 399 ISOnComm_Block, 400 ISSetBlockSize_Block, 401 0, 402 ISLocate_Block, 403 /* we can have specialized local routines for determining properties, 404 * but unless the block size is the same on each process (which is not guaranteed at 405 * the moment), then trying to do something specialized for global properties is too 406 * complicated */ 407 ISSortedLocal_Block, 408 NULL, 409 ISUniqueLocal_Block, 410 NULL, 411 ISPermutationLocal_Block, 412 NULL, 413 ISIntervalLocal_Block, 414 NULL}; 415 416 /*@ 417 ISBlockSetIndices - The indices are relative to entries, not blocks. 418 419 Collective on IS 420 421 Input Parameters: 422 + is - the index set 423 . bs - number of elements in each block, one for each block and count of block not indices 424 . n - the length of the index set (the number of blocks) 425 . idx - the list of integers, these are by block, not by location 426 - mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported 427 428 429 Notes: 430 When the communicator is not MPI_COMM_SELF, the operations on the 431 index sets, IS, are NOT conceptually the same as MPI_Group operations. 432 The index sets are then distributed sets of indices and thus certain operations 433 on them are collective. 434 435 Example: 436 If you wish to index the values {0,1,4,5}, then use 437 a block size of 2 and idx of {0,2}. 438 439 Level: beginner 440 441 .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather() 442 @*/ 443 PetscErrorCode ISBlockSetIndices(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode) 444 { 445 PetscErrorCode ierr; 446 447 PetscFunctionBegin; 448 ierr = ISClearInfoCache(is,PETSC_FALSE);CHKERRQ(ierr); 449 ierr = PetscUseMethod(is,"ISBlockSetIndices_C",(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode),(is,bs,n,idx,mode));CHKERRQ(ierr); 450 PetscFunctionReturn(0); 451 } 452 453 static PetscErrorCode ISBlockSetIndices_Block(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode) 454 { 455 PetscErrorCode ierr; 456 PetscInt i,min,max; 457 IS_Block *sub = (IS_Block*)is->data; 458 PetscLayout map; 459 460 PetscFunctionBegin; 461 if (bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1"); 462 if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0"); 463 if (n) PetscValidIntPointer(idx,3); 464 465 ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n*bs,is->map->N,bs,&map);CHKERRQ(ierr); 466 ierr = PetscLayoutDestroy(&is->map);CHKERRQ(ierr); 467 is->map = map; 468 469 if (sub->allocated) {ierr = PetscFree(sub->idx);CHKERRQ(ierr);} 470 if (mode == PETSC_COPY_VALUES) { 471 ierr = PetscMalloc1(n,&sub->idx);CHKERRQ(ierr); 472 ierr = PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));CHKERRQ(ierr); 473 ierr = PetscArraycpy(sub->idx,idx,n);CHKERRQ(ierr); 474 sub->allocated = PETSC_TRUE; 475 } else if (mode == PETSC_OWN_POINTER) { 476 sub->idx = (PetscInt*) idx; 477 ierr = PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));CHKERRQ(ierr); 478 sub->allocated = PETSC_TRUE; 479 } else if (mode == PETSC_USE_POINTER) { 480 sub->idx = (PetscInt*) idx; 481 sub->allocated = PETSC_FALSE; 482 } 483 484 if (n) { 485 min = max = idx[0]; 486 for (i=1; i<n; i++) { 487 if (idx[i] < min) min = idx[i]; 488 if (idx[i] > max) max = idx[i]; 489 } 490 is->min = bs*min; 491 is->max = bs*max+bs-1; 492 } else { 493 is->min = PETSC_MAX_INT; 494 is->max = PETSC_MIN_INT; 495 } 496 PetscFunctionReturn(0); 497 } 498 499 /*@ 500 ISCreateBlock - Creates a data structure for an index set containing 501 a list of integers. The indices are relative to entries, not blocks. 502 503 Collective 504 505 Input Parameters: 506 + comm - the MPI communicator 507 . bs - number of elements in each block 508 . n - the length of the index set (the number of blocks) 509 . idx - the list of integers, one for each block and count of block not indices 510 - mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine 511 512 Output Parameter: 513 . is - the new index set 514 515 Notes: 516 When the communicator is not MPI_COMM_SELF, the operations on the 517 index sets, IS, are NOT conceptually the same as MPI_Group operations. 518 The index sets are then distributed sets of indices and thus certain operations 519 on them are collective. 520 521 Example: 522 If you wish to index the values {0,1,6,7}, then use 523 a block size of 2 and idx of {0,3}. 524 525 Level: beginner 526 527 .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather() 528 @*/ 529 PetscErrorCode ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode,IS *is) 530 { 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 PetscValidPointer(is,5); 535 if (bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1"); 536 if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0"); 537 if (n) PetscValidIntPointer(idx,4); 538 539 ierr = ISCreate(comm,is);CHKERRQ(ierr); 540 ierr = ISSetType(*is,ISBLOCK);CHKERRQ(ierr); 541 ierr = ISBlockSetIndices(*is,bs,n,idx,mode);CHKERRQ(ierr); 542 PetscFunctionReturn(0); 543 } 544 545 static PetscErrorCode ISBlockGetIndices_Block(IS is,const PetscInt *idx[]) 546 { 547 IS_Block *sub = (IS_Block*)is->data; 548 549 PetscFunctionBegin; 550 *idx = sub->idx; 551 PetscFunctionReturn(0); 552 } 553 554 static PetscErrorCode ISBlockRestoreIndices_Block(IS is,const PetscInt *idx[]) 555 { 556 PetscFunctionBegin; 557 PetscFunctionReturn(0); 558 } 559 560 /*@C 561 ISBlockGetIndices - Gets the indices associated with each block. 562 563 Not Collective 564 565 Input Parameter: 566 . is - the index set 567 568 Output Parameter: 569 . idx - the integer indices, one for each block and count of block not indices 570 571 Level: intermediate 572 573 .seealso: ISGetIndices(), ISBlockRestoreIndices() 574 @*/ 575 PetscErrorCode ISBlockGetIndices(IS is,const PetscInt *idx[]) 576 { 577 PetscErrorCode ierr; 578 579 PetscFunctionBegin; 580 ierr = PetscUseMethod(is,"ISBlockGetIndices_C",(IS,const PetscInt*[]),(is,idx));CHKERRQ(ierr); 581 PetscFunctionReturn(0); 582 } 583 584 /*@C 585 ISBlockRestoreIndices - Restores the indices associated with each block. 586 587 Not Collective 588 589 Input Parameter: 590 . is - the index set 591 592 Output Parameter: 593 . idx - the integer indices 594 595 Level: intermediate 596 597 .seealso: ISRestoreIndices(), ISBlockGetIndices() 598 @*/ 599 PetscErrorCode ISBlockRestoreIndices(IS is,const PetscInt *idx[]) 600 { 601 PetscErrorCode ierr; 602 603 PetscFunctionBegin; 604 ierr = PetscUseMethod(is,"ISBlockRestoreIndices_C",(IS,const PetscInt*[]),(is,idx));CHKERRQ(ierr); 605 PetscFunctionReturn(0); 606 } 607 608 /*@ 609 ISBlockGetLocalSize - Returns the local number of blocks in the index set. 610 611 Not Collective 612 613 Input Parameter: 614 . is - the index set 615 616 Output Parameter: 617 . size - the local number of blocks 618 619 Level: intermediate 620 621 622 .seealso: ISGetBlockSize(), ISBlockGetSize(), ISGetSize(), ISCreateBlock() 623 @*/ 624 PetscErrorCode ISBlockGetLocalSize(IS is,PetscInt *size) 625 { 626 PetscErrorCode ierr; 627 628 PetscFunctionBegin; 629 ierr = PetscUseMethod(is,"ISBlockGetLocalSize_C",(IS,PetscInt*),(is,size));CHKERRQ(ierr); 630 PetscFunctionReturn(0); 631 } 632 633 static PetscErrorCode ISBlockGetLocalSize_Block(IS is,PetscInt *size) 634 { 635 PetscInt bs, n; 636 PetscErrorCode ierr; 637 638 PetscFunctionBegin; 639 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 640 ierr = PetscLayoutGetLocalSize(is->map, &n);CHKERRQ(ierr); 641 *size = n/bs; 642 PetscFunctionReturn(0); 643 } 644 645 /*@ 646 ISBlockGetSize - Returns the global number of blocks in the index set. 647 648 Not Collective 649 650 Input Parameter: 651 . is - the index set 652 653 Output Parameter: 654 . size - the global number of blocks 655 656 Level: intermediate 657 658 659 .seealso: ISGetBlockSize(), ISBlockGetLocalSize(), ISGetSize(), ISCreateBlock() 660 @*/ 661 PetscErrorCode ISBlockGetSize(IS is,PetscInt *size) 662 { 663 PetscErrorCode ierr; 664 665 PetscFunctionBegin; 666 ierr = PetscUseMethod(is,"ISBlockGetSize_C",(IS,PetscInt*),(is,size));CHKERRQ(ierr); 667 PetscFunctionReturn(0); 668 } 669 670 static PetscErrorCode ISBlockGetSize_Block(IS is,PetscInt *size) 671 { 672 PetscInt bs, N; 673 PetscErrorCode ierr; 674 675 PetscFunctionBegin; 676 ierr = PetscLayoutGetBlockSize(is->map, &bs);CHKERRQ(ierr); 677 ierr = PetscLayoutGetSize(is->map, &N);CHKERRQ(ierr); 678 *size = N/bs; 679 PetscFunctionReturn(0); 680 } 681 682 PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is) 683 { 684 PetscErrorCode ierr; 685 IS_Block *sub; 686 687 PetscFunctionBegin; 688 ierr = PetscNewLog(is,&sub);CHKERRQ(ierr); 689 is->data = (void *) sub; 690 ierr = PetscMemcpy(is->ops,&myops,sizeof(myops));CHKERRQ(ierr); 691 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",ISBlockSetIndices_Block);CHKERRQ(ierr); 692 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",ISBlockGetIndices_Block);CHKERRQ(ierr); 693 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",ISBlockRestoreIndices_Block);CHKERRQ(ierr); 694 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",ISBlockGetSize_Block);CHKERRQ(ierr); 695 ierr = PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",ISBlockGetLocalSize_Block);CHKERRQ(ierr); 696 PetscFunctionReturn(0); 697 } 698