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