Actual source code: block.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  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>
  7:  #include <petscviewer.h>

  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;

 15: static PetscErrorCode ISDestroy_Block(IS is)
 16: {
 17:   IS_Block       *sub = (IS_Block*)is->data;

 21:   if (sub->allocated) {PetscFree(sub->idx);}
 22:   PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",NULL);
 23:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",NULL);
 24:   PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",NULL);
 25:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",NULL);
 26:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",NULL);
 27:   PetscFree(is->data);
 28:   return(0);
 29: }

 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;

 39:   PetscLayoutGetBlockSize(is->map,&bs);
 40:   PetscLayoutGetSize(is->map,&numIdx);
 41:   numIdx /= bs;
 42:   bkey    = key / bs;
 43:   mkey    = key % bs;
 44:   if (mkey < 0) {
 45:     bkey--;
 46:     mkey += bs;
 47:   }
 48:   ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted);
 49:   if (sorted) {
 50:     PetscFindInt(bkey,numIdx,sub->idx,location);
 51:   } else {
 52:     const PetscInt *idx = sub->idx;

 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:   return(0);
 66: }

 68: static PetscErrorCode ISGetIndices_Block(IS in,const PetscInt *idx[])
 69: {
 70:   IS_Block       *sub = (IS_Block*)in->data;
 72:   PetscInt       i,j,k,bs,n,*ii,*jj;

 75:   PetscLayoutGetBlockSize(in->map, &bs);
 76:   PetscLayoutGetLocalSize(in->map, &n);
 77:   n   /= bs;
 78:   if (bs == 1) *idx = sub->idx;
 79:   else {
 80:     PetscMalloc1(bs*n,&jj);
 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:   return(0);
 89: }

 91: static PetscErrorCode ISRestoreIndices_Block(IS is,const PetscInt *idx[])
 92: {
 93:   IS_Block       *sub = (IS_Block*)is->data;
 94:   PetscInt       bs;

 98:   PetscLayoutGetBlockSize(is->map, &bs);
 99:   if (bs != 1) {
100:     PetscFree(*(void**)idx);
101:   } else {
102:     if (*idx != sub->idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must restore with value from ISGetIndices()");
103:   }
104:   return(0);
105: }

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;

115:   MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
116:   PetscLayoutGetBlockSize(is->map, &bs);
117:   PetscLayoutGetLocalSize(is->map, &n);
118:   n   /= bs;
119:   if (size == 1) {
120:     PetscMalloc1(n,&ii);
121:     for (i=0; i<n; i++) ii[idx[i]] = i;
122:     ISCreateBlock(PETSC_COMM_SELF,bs,n,ii,PETSC_OWN_POINTER,isout);
123:     ISSetPermutation(*isout);
124:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No inversion written yet for block IS");
125:   return(0);
126: }

128: static PetscErrorCode ISView_Block(IS is, PetscViewer viewer)
129: {
130:   IS_Block       *sub = (IS_Block*)is->data;
132:   PetscInt       i,bs,n,*idx = sub->idx;
133:   PetscBool      iascii,ibinary;

136:   PetscLayoutGetBlockSize(is->map, &bs);
137:   PetscLayoutGetLocalSize(is->map, &n);
138:   n   /= bs;
139:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
140:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary);
141:   if (iascii) {
142:     PetscViewerFormat fmt;

144:     PetscViewerGetFormat(viewer,&fmt);
145:     if (fmt == PETSC_VIEWER_ASCII_MATLAB) {
146:       IS             ist;
147:       const char     *name;
148:       const PetscInt *idx;
149:       PetscInt       n;

151:       PetscObjectGetName((PetscObject)is,&name);
152:       ISGetLocalSize(is,&n);
153:       ISGetIndices(is,&idx);
154:       ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idx,PETSC_USE_POINTER,&ist);
155:       PetscObjectSetName((PetscObject)ist,name);
156:       ISView(ist,viewer);
157:       ISDestroy(&ist);
158:       ISRestoreIndices(is,&idx);
159:     } else {
160:       PetscBool isperm;

162:       ISGetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_FALSE,&isperm);
163:       if (isperm) {PetscViewerASCIIPrintf(viewer,"Block Index set is permutation\n");}
164:       PetscViewerASCIIPushSynchronized(viewer);
165:       PetscViewerASCIISynchronizedPrintf(viewer,"Block size %D\n",bs);
166:       PetscViewerASCIISynchronizedPrintf(viewer,"Number of block indices in set %D\n",n);
167:       PetscViewerASCIISynchronizedPrintf(viewer,"The first indices of each block are\n");
168:       for (i=0; i<n; i++) {
169:         PetscViewerASCIISynchronizedPrintf(viewer,"Block %D Index %D\n",i,idx[i]);
170:       }
171:       PetscViewerFlush(viewer);
172:       PetscViewerASCIIPopSynchronized(viewer);
173:     }
174:   } else if (ibinary) {
175:     ISView_Binary(is,viewer);
176:   }
177:   return(0);
178: }

180: static PetscErrorCode ISSort_Block(IS is)
181: {
182:   IS_Block       *sub = (IS_Block*)is->data;
183:   PetscInt       bs, n;

187:   PetscLayoutGetBlockSize(is->map, &bs);
188:   PetscLayoutGetLocalSize(is->map, &n);
189:   PetscSortInt(n/bs,sub->idx);
190:   return(0);
191: }

193: static PetscErrorCode ISSortRemoveDups_Block(IS is)
194: {
195:   IS_Block       *sub = (IS_Block*)is->data;
196:   PetscInt       bs, n, nb;
197:   PetscBool      sorted;

201:   PetscLayoutGetBlockSize(is->map, &bs);
202:   PetscLayoutGetLocalSize(is->map, &n);
203:   nb   = n/bs;
204:   ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sorted);
205:   if (sorted) {
206:     PetscSortedRemoveDupsInt(&nb,sub->idx);
207:   } else {
208:     PetscSortRemoveDupsInt(&nb,sub->idx);
209:   }
210:   PetscLayoutDestroy(&is->map);
211:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), nb*bs, PETSC_DECIDE, bs, &is->map);
212:   return(0);
213: }

215: static PetscErrorCode ISSorted_Block(IS is,PetscBool  *flg)
216: {

220:   ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,flg);
221:   return(0);
222: }

224: static PetscErrorCode ISSortedLocal_Block(IS is,PetscBool *flg)
225: {
226:   IS_Block       *sub = (IS_Block*)is->data;
227:   PetscInt       n, bs, i, *idx;

231:   PetscLayoutGetLocalSize(is->map, &n);
232:   PetscLayoutGetBlockSize(is->map, &bs);
233:   n   /= bs;
234:   idx  = sub->idx;
235:   for (i = 1; i < n; i++) if (idx[i] < idx[i - 1]) break;
236:   if (i < n) *flg = PETSC_FALSE;
237:   else       *flg = PETSC_TRUE;
238:   return(0);
239: }

241: static PetscErrorCode ISUniqueLocal_Block(IS is,PetscBool *flg)
242: {
243:   IS_Block       *sub = (IS_Block*)is->data;
244:   PetscInt       n, bs, i, *idx, *idxcopy = NULL;
245:   PetscBool      sortedLocal;

249:   PetscLayoutGetLocalSize(is->map, &n);
250:   PetscLayoutGetBlockSize(is->map, &bs);
251:   n   /= bs;
252:   idx  = sub->idx;
253:   ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal);
254:   if (!sortedLocal) {
255:     PetscMalloc1(n, &idxcopy);
256:     PetscArraycpy(idxcopy, idx, n);
257:     PetscSortInt(n, idxcopy);
258:     idx = idxcopy;
259:   }
260:   for (i = 1; i < n; i++) if (idx[i] == idx[i - 1]) break;
261:   if (i < n) *flg = PETSC_FALSE;
262:   else       *flg = PETSC_TRUE;
263:   PetscFree(idxcopy);
264:   return(0);
265: }

267: static PetscErrorCode ISPermutationLocal_Block(IS is,PetscBool *flg)
268: {
269:   IS_Block       *sub = (IS_Block*)is->data;
270:   PetscInt       n, bs, i, *idx, *idxcopy = NULL;
271:   PetscBool      sortedLocal;

275:   PetscLayoutGetLocalSize(is->map, &n);
276:   PetscLayoutGetBlockSize(is->map, &bs);
277:   n   /= bs;
278:   idx  = sub->idx;
279:   ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,&sortedLocal);
280:   if (!sortedLocal) {
281:     PetscMalloc1(n, &idxcopy);
282:     PetscArraycpy(idxcopy, idx, n);
283:     PetscSortInt(n, idxcopy);
284:     idx = idxcopy;
285:   }
286:   for (i = 0; i < n; i++) if (idx[i] != i) break;
287:   if (i < n) *flg = PETSC_FALSE;
288:   else       *flg = PETSC_TRUE;
289:   PetscFree(idxcopy);
290:   return(0);
291: }

293: static PetscErrorCode ISIntervalLocal_Block(IS is,PetscBool *flg)
294: {
295:   IS_Block       *sub = (IS_Block*)is->data;
296:   PetscInt       n, bs, i, *idx;

300:   PetscLayoutGetLocalSize(is->map, &n);
301:   PetscLayoutGetBlockSize(is->map, &bs);
302:   n   /= bs;
303:   idx  = sub->idx;
304:   for (i = 1; i < n; i++) if (idx[i] != idx[i - 1] + 1) break;
305:   if (i < n) *flg = PETSC_FALSE;
306:   else       *flg = PETSC_TRUE;
307:   return(0);
308: }

310: static PetscErrorCode ISDuplicate_Block(IS is,IS *newIS)
311: {
313:   IS_Block       *sub = (IS_Block*)is->data;
314:   PetscInt        bs, n;

317:   PetscLayoutGetBlockSize(is->map, &bs);
318:   PetscLayoutGetLocalSize(is->map, &n);
319:   n   /= bs;
320:   ISCreateBlock(PetscObjectComm((PetscObject)is),bs,n,sub->idx,PETSC_COPY_VALUES,newIS);
321:   return(0);
322: }

324: static PetscErrorCode ISCopy_Block(IS is,IS isy)
325: {
326:   IS_Block       *is_block = (IS_Block*)is->data,*isy_block = (IS_Block*)isy->data;
327:   PetscInt       bs, n, N, bsy, ny, Ny;

331:   PetscLayoutGetBlockSize(is->map, &bs);
332:   PetscLayoutGetLocalSize(is->map, &n);
333:   PetscLayoutGetSize(is->map, &N);
334:   PetscLayoutGetBlockSize(isy->map, &bsy);
335:   PetscLayoutGetLocalSize(isy->map, &ny);
336:   PetscLayoutGetSize(isy->map, &Ny);
337:   if (n != ny || N != Ny || bs != bsy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Index sets incompatible");
338:   PetscArraycpy(isy_block->idx,is_block->idx,(n/bs));
339:   return(0);
340: }

342: static PetscErrorCode ISOnComm_Block(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
343: {
345:   IS_Block       *sub = (IS_Block*)is->data;
346:   PetscInt       bs, n;

349:   if (mode == PETSC_OWN_POINTER) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Cannot use PETSC_OWN_POINTER");
350:   PetscLayoutGetBlockSize(is->map, &bs);
351:   PetscLayoutGetLocalSize(is->map, &n);
352:   ISCreateBlock(comm,bs,n/bs,sub->idx,mode,newis);
353:   return(0);
354: }

356: static PetscErrorCode ISSetBlockSize_Block(IS is,PetscInt bs)
357: {

361:   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);
362:   PetscLayoutSetBlockSize(is->map, bs);
363:   return(0);
364: }

366: static PetscErrorCode ISToGeneral_Block(IS inis)
367: {
368:   IS_Block       *sub   = (IS_Block*)inis->data;
369:   PetscInt       bs,n;
370:   const PetscInt *idx;

374:   ISGetBlockSize(inis,&bs);
375:   ISGetLocalSize(inis,&n);
376:   ISGetIndices(inis,&idx);
377:   if (bs == 1) {
378:     PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER;
379:     sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/
380:     ISSetType(inis,ISGENERAL);
381:     ISGeneralSetIndices(inis,n,idx,mode);
382:   } else {
383:     ISSetType(inis,ISGENERAL);
384:     ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);
385:   }
386:   return(0);
387: }


390: static struct _ISOps myops = { ISGetIndices_Block,
391:                                ISRestoreIndices_Block,
392:                                ISInvertPermutation_Block,
393:                                ISSort_Block,
394:                                ISSortRemoveDups_Block,
395:                                ISSorted_Block,
396:                                ISDuplicate_Block,
397:                                ISDestroy_Block,
398:                                ISView_Block,
399:                                ISLoad_Default,
400:                                ISCopy_Block,
401:                                ISToGeneral_Block,
402:                                ISOnComm_Block,
403:                                ISSetBlockSize_Block,
404:                                0,
405:                                ISLocate_Block,
406:                                /* we can have specialized local routines for determining properties,
407:                                 * but unless the block size is the same on each process (which is not guaranteed at
408:                                 * the moment), then trying to do something specialized for global properties is too
409:                                 * complicated */
410:                                ISSortedLocal_Block,
411:                                NULL,
412:                                ISUniqueLocal_Block,
413:                                NULL,
414:                                ISPermutationLocal_Block,
415:                                NULL,
416:                                ISIntervalLocal_Block,
417:                                NULL};

419: /*@
420:    ISBlockSetIndices - The indices are relative to entries, not blocks.

422:    Collective on IS

424:    Input Parameters:
425: +  is - the index set
426: .  bs - number of elements in each block, one for each block and count of block not indices
427: .   n - the length of the index set (the number of blocks)
428: .  idx - the list of integers, these are by block, not by location
429: -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported


432:    Notes:
433:    When the communicator is not MPI_COMM_SELF, the operations on the
434:    index sets, IS, are NOT conceptually the same as MPI_Group operations.
435:    The index sets are then distributed sets of indices and thus certain operations
436:    on them are collective.

438:    Example:
439:    If you wish to index the values {0,1,4,5}, then use
440:    a block size of 2 and idx of {0,2}.

442:    Level: beginner

444: .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather()
445: @*/
446: PetscErrorCode  ISBlockSetIndices(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
447: {

451:   ISClearInfoCache(is,PETSC_FALSE);
452:   PetscUseMethod(is,"ISBlockSetIndices_C",(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode),(is,bs,n,idx,mode));
453:   return(0);
454: }

456: static PetscErrorCode  ISBlockSetIndices_Block(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
457: {
459:   PetscInt       i,min,max;
460:   IS_Block       *sub = (IS_Block*)is->data;
461:   PetscLayout    map;

464:   if (bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1");
465:   if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0");

468:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n*bs,is->map->N,bs,&map);
469:   PetscLayoutDestroy(&is->map);
470:   is->map = map;

472:   if (sub->allocated) {PetscFree(sub->idx);}
473:   if (mode == PETSC_COPY_VALUES) {
474:     PetscMalloc1(n,&sub->idx);
475:     PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));
476:     PetscArraycpy(sub->idx,idx,n);
477:     sub->allocated = PETSC_TRUE;
478:   } else if (mode == PETSC_OWN_POINTER) {
479:     sub->idx = (PetscInt*) idx;
480:     PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));
481:     sub->allocated = PETSC_TRUE;
482:   } else if (mode == PETSC_USE_POINTER) {
483:     sub->idx = (PetscInt*) idx;
484:     sub->allocated = PETSC_FALSE;
485:   }

487:   if (n) {
488:     min = max = idx[0];
489:     for (i=1; i<n; i++) {
490:       if (idx[i] < min) min = idx[i];
491:       if (idx[i] > max) max = idx[i];
492:     }
493:     is->min = bs*min;
494:     is->max = bs*max+bs-1;
495:   } else {
496:     is->min = PETSC_MAX_INT;
497:     is->max = PETSC_MIN_INT;
498:   }
499:   return(0);
500: }

502: /*@
503:    ISCreateBlock - Creates a data structure for an index set containing
504:    a list of integers. The indices are relative to entries, not blocks.

506:    Collective

508:    Input Parameters:
509: +  comm - the MPI communicator
510: .  bs - number of elements in each block
511: .  n - the length of the index set (the number of blocks)
512: .  idx - the list of integers, one for each block and count of block not indices
513: -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine

515:    Output Parameter:
516: .  is - the new index set

518:    Notes:
519:    When the communicator is not MPI_COMM_SELF, the operations on the
520:    index sets, IS, are NOT conceptually the same as MPI_Group operations.
521:    The index sets are then distributed sets of indices and thus certain operations
522:    on them are collective.

524:    Example:
525:    If you wish to index the values {0,1,6,7}, then use
526:    a block size of 2 and idx of {0,3}.

528:    Level: beginner

530: .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather()
531: @*/
532: PetscErrorCode  ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode,IS *is)
533: {

538:   if (bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"block size < 1");
539:   if (n < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"length < 0");

542:   ISCreate(comm,is);
543:   ISSetType(*is,ISBLOCK);
544:   ISBlockSetIndices(*is,bs,n,idx,mode);
545:   return(0);
546: }

548: static PetscErrorCode  ISBlockGetIndices_Block(IS is,const PetscInt *idx[])
549: {
550:   IS_Block *sub = (IS_Block*)is->data;

553:   *idx = sub->idx;
554:   return(0);
555: }

557: static PetscErrorCode  ISBlockRestoreIndices_Block(IS is,const PetscInt *idx[])
558: {
560:   return(0);
561: }

563: /*@C
564:    ISBlockGetIndices - Gets the indices associated with each block.

566:    Not Collective

568:    Input Parameter:
569: .  is - the index set

571:    Output Parameter:
572: .  idx - the integer indices, one for each block and count of block not indices

574:    Level: intermediate

576: .seealso: ISGetIndices(), ISBlockRestoreIndices()
577: @*/
578: PetscErrorCode  ISBlockGetIndices(IS is,const PetscInt *idx[])
579: {

583:   PetscUseMethod(is,"ISBlockGetIndices_C",(IS,const PetscInt*[]),(is,idx));
584:   return(0);
585: }

587: /*@C
588:    ISBlockRestoreIndices - Restores the indices associated with each block.

590:    Not Collective

592:    Input Parameter:
593: .  is - the index set

595:    Output Parameter:
596: .  idx - the integer indices

598:    Level: intermediate

600: .seealso: ISRestoreIndices(), ISBlockGetIndices()
601: @*/
602: PetscErrorCode  ISBlockRestoreIndices(IS is,const PetscInt *idx[])
603: {

607:   PetscUseMethod(is,"ISBlockRestoreIndices_C",(IS,const PetscInt*[]),(is,idx));
608:   return(0);
609: }

611: /*@
612:    ISBlockGetLocalSize - Returns the local number of blocks in the index set.

614:    Not Collective

616:    Input Parameter:
617: .  is - the index set

619:    Output Parameter:
620: .  size - the local number of blocks

622:    Level: intermediate


625: .seealso: ISGetBlockSize(), ISBlockGetSize(), ISGetSize(), ISCreateBlock()
626: @*/
627: PetscErrorCode  ISBlockGetLocalSize(IS is,PetscInt *size)
628: {

632:   PetscUseMethod(is,"ISBlockGetLocalSize_C",(IS,PetscInt*),(is,size));
633:   return(0);
634: }

636: static PetscErrorCode  ISBlockGetLocalSize_Block(IS is,PetscInt *size)
637: {
638:   PetscInt       bs, n;

642:   PetscLayoutGetBlockSize(is->map, &bs);
643:   PetscLayoutGetLocalSize(is->map, &n);
644:   *size = n/bs;
645:   return(0);
646: }

648: /*@
649:    ISBlockGetSize - Returns the global number of blocks in the index set.

651:    Not Collective

653:    Input Parameter:
654: .  is - the index set

656:    Output Parameter:
657: .  size - the global number of blocks

659:    Level: intermediate


662: .seealso: ISGetBlockSize(), ISBlockGetLocalSize(), ISGetSize(), ISCreateBlock()
663: @*/
664: PetscErrorCode  ISBlockGetSize(IS is,PetscInt *size)
665: {

669:   PetscUseMethod(is,"ISBlockGetSize_C",(IS,PetscInt*),(is,size));
670:   return(0);
671: }

673: static PetscErrorCode  ISBlockGetSize_Block(IS is,PetscInt *size)
674: {
675:   PetscInt       bs, N;

679:   PetscLayoutGetBlockSize(is->map, &bs);
680:   PetscLayoutGetSize(is->map, &N);
681:   *size = N/bs;
682:   return(0);
683: }

685: PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is)
686: {
688:   IS_Block       *sub;

691:   PetscNewLog(is,&sub);
692:   is->data = (void *) sub;
693:   PetscMemcpy(is->ops,&myops,sizeof(myops));
694:   PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",ISBlockSetIndices_Block);
695:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",ISBlockGetIndices_Block);
696:   PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",ISBlockRestoreIndices_Block);
697:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",ISBlockGetSize_Block);
698:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",ISBlockGetLocalSize_Block);
699:   return(0);
700: }