Actual source code: block.c

petsc-3.12.2 2019-11-22
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;

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

 52:     *location = -1;
 53:     for (i = 0; i < numIdx; i++) {
 54:       if (idx[i] == bkey) {
 55:         *location = i;
 56:         break;
 57:       }
 58:     }
 59:   }
 60:   if (*location >= 0) {
 61:     *location = *location * bs + mkey;
 62:   }
 63:   return(0);
 64: }

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

 73:   PetscLayoutGetBlockSize(in->map, &bs);
 74:   PetscLayoutGetLocalSize(in->map, &n);
 75:   n   /= bs;
 76:   if (bs == 1) *idx = sub->idx;
 77:   else {
 78:     PetscMalloc1(bs*n,&jj);
 79:     *idx = jj;
 80:     k    = 0;
 81:     ii   = sub->idx;
 82:     for (i=0; i<n; i++)
 83:       for (j=0; j<bs; j++)
 84:         jj[k++] = bs*ii[i] + j;
 85:   }
 86:   return(0);
 87: }

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

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

105: static PetscErrorCode ISInvertPermutation_Block(IS is,PetscInt nlocal,IS *isout)
106: {
107:   IS_Block       *sub = (IS_Block*)is->data;
108:   PetscInt       i,*ii,bs,n,*idx = sub->idx;
109:   PetscMPIInt    size;

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

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

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

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

148:       PetscObjectGetName((PetscObject)is,&name);
149:       ISGetLocalSize(is,&n);
150:       ISGetIndices(is,&idx);
151:       ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idx,PETSC_USE_POINTER,&ist);
152:       PetscObjectSetName((PetscObject)ist,name);
153:       ISView(ist,viewer);
154:       ISDestroy(&ist);
155:       ISRestoreIndices(is,&idx);
156:     } else {
157:       PetscViewerASCIIPushSynchronized(viewer);
158:       if (is->isperm) {
159:         PetscViewerASCIISynchronizedPrintf(viewer,"Block Index set is permutation\n");
160:       }
161:       PetscViewerASCIISynchronizedPrintf(viewer,"Block size %D\n",bs);
162:       PetscViewerASCIISynchronizedPrintf(viewer,"Number of block indices in set %D\n",n);
163:       PetscViewerASCIISynchronizedPrintf(viewer,"The first indices of each block are\n");
164:       for (i=0; i<n; i++) {
165:         PetscViewerASCIISynchronizedPrintf(viewer,"Block %D Index %D\n",i,idx[i]);
166:       }
167:       PetscViewerFlush(viewer);
168:       PetscViewerASCIIPopSynchronized(viewer);
169:     }
170:   }
171:   return(0);
172: }

174: static PetscErrorCode ISSort_Block(IS is)
175: {
176:   IS_Block       *sub = (IS_Block*)is->data;
177:   PetscInt       bs, n;

181:   if (sub->sorted) return(0);
182:   PetscLayoutGetBlockSize(is->map, &bs);
183:   PetscLayoutGetLocalSize(is->map, &n);
184:   PetscSortInt(n/bs,sub->idx);
185:   sub->sorted = PETSC_TRUE;
186:   return(0);
187: }

189: static PetscErrorCode ISSortRemoveDups_Block(IS is)
190: {
191:   IS_Block       *sub = (IS_Block*)is->data;
192:   PetscInt       bs, n, nb;

196:   PetscLayoutGetBlockSize(is->map, &bs);
197:   PetscLayoutGetLocalSize(is->map, &n);
198:   nb   = n/bs;
199:   if (sub->sorted) {
200:     PetscSortedRemoveDupsInt(&nb,sub->idx);
201:   } else {
202:     PetscSortRemoveDupsInt(&nb,sub->idx);
203:   }
204:   PetscLayoutDestroy(&is->map);
205:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), nb*bs, PETSC_DECIDE, bs, &is->map);
206:   sub->sorted = PETSC_TRUE;
207:   return(0);
208: }

210: static PetscErrorCode ISSorted_Block(IS is,PetscBool  *flg)
211: {
212:   IS_Block *sub = (IS_Block*)is->data;

215:   *flg = sub->sorted;
216:   return(0);
217: }

219: static PetscErrorCode ISDuplicate_Block(IS is,IS *newIS)
220: {
222:   IS_Block       *sub = (IS_Block*)is->data;
223:   PetscInt        bs, n;

226:   PetscLayoutGetBlockSize(is->map, &bs);
227:   PetscLayoutGetLocalSize(is->map, &n);
228:   n   /= bs;
229:   ISCreateBlock(PetscObjectComm((PetscObject)is),bs,n,sub->idx,PETSC_COPY_VALUES,newIS);
230:   return(0);
231: }

233: static PetscErrorCode ISIdentity_Block(IS is,PetscBool  *ident)
234: {
235:   IS_Block      *is_block = (IS_Block*)is->data;
236:   PetscInt       i,bs,n,*idx = is_block->idx;

240:   PetscLayoutGetBlockSize(is->map, &bs);
241:   PetscLayoutGetLocalSize(is->map, &n);
242:   n   /= bs;
243:   is->isidentity = PETSC_TRUE;
244:   *ident         = PETSC_TRUE;
245:   for (i=0; i<n; i++) {
246:     if (idx[i] != i) {
247:       is->isidentity = PETSC_FALSE;
248:       *ident         = PETSC_FALSE;
249:       return(0);
250:     }
251:   }
252:   return(0);
253: }

255: static PetscErrorCode ISCopy_Block(IS is,IS isy)
256: {
257:   IS_Block       *is_block = (IS_Block*)is->data,*isy_block = (IS_Block*)isy->data;
258:   PetscInt       bs, n, N, bsy, ny, Ny;

262:   PetscLayoutGetBlockSize(is->map, &bs);
263:   PetscLayoutGetLocalSize(is->map, &n);
264:   PetscLayoutGetSize(is->map, &N);
265:   PetscLayoutGetBlockSize(isy->map, &bsy);
266:   PetscLayoutGetLocalSize(isy->map, &ny);
267:   PetscLayoutGetSize(isy->map, &Ny);
268:   if (n != ny || N != Ny || bs != bsy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Index sets incompatible");
269:   isy_block->sorted = is_block->sorted;
270:   PetscArraycpy(isy_block->idx,is_block->idx,(n/bs));
271:   return(0);
272: }

274: static PetscErrorCode ISOnComm_Block(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
275: {
277:   IS_Block       *sub = (IS_Block*)is->data;
278:   PetscInt       bs, n;

281:   if (mode == PETSC_OWN_POINTER) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Cannot use PETSC_OWN_POINTER");
282:   PetscLayoutGetBlockSize(is->map, &bs);
283:   PetscLayoutGetLocalSize(is->map, &n);
284:   ISCreateBlock(comm,bs,n/bs,sub->idx,mode,newis);
285:   return(0);
286: }

288: static PetscErrorCode ISSetBlockSize_Block(IS is,PetscInt bs)
289: {

293:   PetscLayoutSetBlockSize(is->map, bs);
294:   return(0);
295: }

297: static PetscErrorCode ISToGeneral_Block(IS inis)
298: {
299:   IS_Block       *sub   = (IS_Block*)inis->data;
300:   PetscInt       bs,n;
301:   const PetscInt *idx;

305:   ISGetBlockSize(inis,&bs);
306:   ISGetLocalSize(inis,&n);
307:   ISGetIndices(inis,&idx);
308:   if (bs == 1) {
309:     PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER;
310:     sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/
311:     ISSetType(inis,ISGENERAL);
312:     ISGeneralSetIndices(inis,n,idx,mode);
313:   } else {
314:     ISSetType(inis,ISGENERAL);
315:     ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);
316:   }
317:   return(0);
318: }


321: static struct _ISOps myops = { ISGetIndices_Block,
322:                                ISRestoreIndices_Block,
323:                                ISInvertPermutation_Block,
324:                                ISSort_Block,
325:                                ISSortRemoveDups_Block,
326:                                ISSorted_Block,
327:                                ISDuplicate_Block,
328:                                ISDestroy_Block,
329:                                ISView_Block,
330:                                ISLoad_Default,
331:                                ISIdentity_Block,
332:                                ISCopy_Block,
333:                                ISToGeneral_Block,
334:                                ISOnComm_Block,
335:                                ISSetBlockSize_Block,
336:                                0,
337:                                ISLocate_Block};

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

342:    Collective on IS

344:    Input Parameters:
345: +  is - the index set
346: .  bs - number of elements in each block, one for each block and count of block not indices
347: .   n - the length of the index set (the number of blocks)
348: .  idx - the list of integers, these are by block, not by location
349: -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported


352:    Notes:
353:    When the communicator is not MPI_COMM_SELF, the operations on the
354:    index sets, IS, are NOT conceptually the same as MPI_Group operations.
355:    The index sets are then distributed sets of indices and thus certain operations
356:    on them are collective.

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

362:    Level: beginner

364: .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather()
365: @*/
366: PetscErrorCode  ISBlockSetIndices(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
367: {

371:   PetscUseMethod(is,"ISBlockSetIndices_C",(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode),(is,bs,n,idx,mode));
372:   return(0);
373: }

375: static PetscErrorCode  ISBlockSetIndices_Block(IS is,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode)
376: {
378:   PetscInt       i,min,max;
379:   IS_Block       *sub = (IS_Block*)is->data;
380:   PetscLayout    map;

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

387:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n*bs,is->map->N,bs,&map);
388:   PetscLayoutDestroy(&is->map);
389:   is->map = map;

391:   if (sub->allocated) {PetscFree(sub->idx);}
392:   if (mode == PETSC_COPY_VALUES) {
393:     PetscMalloc1(n,&sub->idx);
394:     PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));
395:     PetscArraycpy(sub->idx,idx,n);
396:     sub->allocated = PETSC_TRUE;
397:   } else if (mode == PETSC_OWN_POINTER) {
398:     sub->idx = (PetscInt*) idx;
399:     PetscLogObjectMemory((PetscObject)is,n*sizeof(PetscInt));
400:     sub->allocated = PETSC_TRUE;
401:   } else if (mode == PETSC_USE_POINTER) {
402:     sub->idx = (PetscInt*) idx;
403:     sub->allocated = PETSC_FALSE;
404:   }

406:   sub->sorted = PETSC_TRUE;
407:   for (i=1; i<n; i++) {
408:     if (idx[i] < idx[i-1]) {sub->sorted = PETSC_FALSE; break;}
409:   }
410:   if (n) {
411:     min = max = idx[0];
412:     for (i=1; i<n; i++) {
413:       if (idx[i] < min) min = idx[i];
414:       if (idx[i] > max) max = idx[i];
415:     }
416:     is->min = bs*min;
417:     is->max = bs*max+bs-1;
418:   } else {
419:     is->min = PETSC_MAX_INT;
420:     is->max = PETSC_MIN_INT;
421:   }
422:   is->isperm     = PETSC_FALSE;
423:   is->isidentity = PETSC_FALSE;
424:   return(0);
425: }

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

431:    Collective

433:    Input Parameters:
434: +  comm - the MPI communicator
435: .  bs - number of elements in each block
436: .  n - the length of the index set (the number of blocks)
437: .  idx - the list of integers, one for each block and count of block not indices
438: -  mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine

440:    Output Parameter:
441: .  is - the new index set

443:    Notes:
444:    When the communicator is not MPI_COMM_SELF, the operations on the
445:    index sets, IS, are NOT conceptually the same as MPI_Group operations.
446:    The index sets are then distributed sets of indices and thus certain operations
447:    on them are collective.

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

453:    Level: beginner

455: .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather()
456: @*/
457: PetscErrorCode  ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],PetscCopyMode mode,IS *is)
458: {

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

467:   ISCreate(comm,is);
468:   ISSetType(*is,ISBLOCK);
469:   ISBlockSetIndices(*is,bs,n,idx,mode);
470:   return(0);
471: }

473: static PetscErrorCode  ISBlockGetIndices_Block(IS is,const PetscInt *idx[])
474: {
475:   IS_Block *sub = (IS_Block*)is->data;

478:   *idx = sub->idx;
479:   return(0);
480: }

482: static PetscErrorCode  ISBlockRestoreIndices_Block(IS is,const PetscInt *idx[])
483: {
485:   return(0);
486: }

488: /*@C
489:    ISBlockGetIndices - Gets the indices associated with each block.

491:    Not Collective

493:    Input Parameter:
494: .  is - the index set

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

499:    Level: intermediate

501: .seealso: ISGetIndices(), ISBlockRestoreIndices()
502: @*/
503: PetscErrorCode  ISBlockGetIndices(IS is,const PetscInt *idx[])
504: {

508:   PetscUseMethod(is,"ISBlockGetIndices_C",(IS,const PetscInt*[]),(is,idx));
509:   return(0);
510: }

512: /*@C
513:    ISBlockRestoreIndices - Restores the indices associated with each block.

515:    Not Collective

517:    Input Parameter:
518: .  is - the index set

520:    Output Parameter:
521: .  idx - the integer indices

523:    Level: intermediate

525: .seealso: ISRestoreIndices(), ISBlockGetIndices()
526: @*/
527: PetscErrorCode  ISBlockRestoreIndices(IS is,const PetscInt *idx[])
528: {

532:   PetscUseMethod(is,"ISBlockRestoreIndices_C",(IS,const PetscInt*[]),(is,idx));
533:   return(0);
534: }

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

539:    Not Collective

541:    Input Parameter:
542: .  is - the index set

544:    Output Parameter:
545: .  size - the local number of blocks

547:    Level: intermediate


550: .seealso: ISGetBlockSize(), ISBlockGetSize(), ISGetSize(), ISCreateBlock()
551: @*/
552: PetscErrorCode  ISBlockGetLocalSize(IS is,PetscInt *size)
553: {

557:   PetscUseMethod(is,"ISBlockGetLocalSize_C",(IS,PetscInt*),(is,size));
558:   return(0);
559: }

561: static PetscErrorCode  ISBlockGetLocalSize_Block(IS is,PetscInt *size)
562: {
563:   PetscInt       bs, n;

567:   PetscLayoutGetBlockSize(is->map, &bs);
568:   PetscLayoutGetLocalSize(is->map, &n);
569:   *size = n/bs;
570:   return(0);
571: }

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

576:    Not Collective

578:    Input Parameter:
579: .  is - the index set

581:    Output Parameter:
582: .  size - the global number of blocks

584:    Level: intermediate


587: .seealso: ISGetBlockSize(), ISBlockGetLocalSize(), ISGetSize(), ISCreateBlock()
588: @*/
589: PetscErrorCode  ISBlockGetSize(IS is,PetscInt *size)
590: {

594:   PetscUseMethod(is,"ISBlockGetSize_C",(IS,PetscInt*),(is,size));
595:   return(0);
596: }

598: static PetscErrorCode  ISBlockGetSize_Block(IS is,PetscInt *size)
599: {
600:   PetscInt       bs, N;

604:   PetscLayoutGetBlockSize(is->map, &bs);
605:   PetscLayoutGetSize(is->map, &N);
606:   *size = N/bs;
607:   return(0);
608: }

610: PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is)
611: {
613:   IS_Block       *sub;

616:   PetscNewLog(is,&sub);
617:   is->data = (void *) sub;
618:   PetscMemcpy(is->ops,&myops,sizeof(myops));
619:   PetscObjectComposeFunction((PetscObject)is,"ISBlockSetIndices_C",ISBlockSetIndices_Block);
620:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetIndices_C",ISBlockGetIndices_Block);
621:   PetscObjectComposeFunction((PetscObject)is,"ISBlockRestoreIndices_C",ISBlockRestoreIndices_Block);
622:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetSize_C",ISBlockGetSize_Block);
623:   PetscObjectComposeFunction((PetscObject)is,"ISBlockGetLocalSize_C",ISBlockGetLocalSize_Block);
624:   return(0);
625: }