Actual source code: block.c

petsc-3.12.4 2020-02-04
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:   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);
294:   PetscLayoutSetBlockSize(is->map, bs);
295:   return(0);
296: }

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

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


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

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

343:    Collective on IS

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


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

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

363:    Level: beginner

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

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

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

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

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

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

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

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

432:    Collective

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

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

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

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

454:    Level: beginner

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

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:   ISCreate(comm,is);
469:   ISSetType(*is,ISBLOCK);
470:   ISBlockSetIndices(*is,bs,n,idx,mode);
471:   return(0);
472: }

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

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

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

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

492:    Not Collective

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

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

500:    Level: intermediate

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

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

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

516:    Not Collective

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

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

524:    Level: intermediate

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

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

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

540:    Not Collective

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

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

548:    Level: intermediate


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

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

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

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

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

577:    Not Collective

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

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

585:    Level: intermediate


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

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

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

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

611: PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is)
612: {
614:   IS_Block       *sub;

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