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