Actual source code: index.c

petsc-master 2019-12-16
Report Typos and Errors
  1: /*
  2:    Defines the abstract operations on index sets, i.e. the public interface.
  3: */
  4:  #include <petsc/private/isimpl.h>
  5:  #include <petscviewer.h>
  6:  #include <petscsf.h>

  8: /* Logging support */
  9: PetscClassId IS_CLASSID;
 10: /* TODO: Much more events are missing! */
 11: PetscLogEvent IS_Load;

 13: /*@
 14:    ISRenumber - Renumbers an index set (with multiplicities) in a contiguous way.

 16:    Collective on IS

 18:    Input Parmeters:
 19: +  subset - the index set
 20: -  subset_mult - the multiplcity of each entry in subset (optional, can be NULL)

 22:    Output Parameters:
 23: +  N - the maximum entry of the new IS
 24: -  subset_n - the new IS

 26:    Level: intermediate

 28: .seealso:
 29: @*/
 30: PetscErrorCode ISRenumber(IS subset, IS subset_mult, PetscInt *N, IS *subset_n)
 31: {
 32:   PetscSF        sf;
 33:   PetscLayout    map;
 34:   const PetscInt *idxs;
 35:   PetscInt       *leaf_data,*root_data,*gidxs;
 36:   PetscInt       N_n,n,i,lbounds[2],gbounds[2],Nl;
 37:   PetscInt       n_n,nlocals,start,first_index;
 38:   PetscMPIInt    commsize;
 39:   PetscBool      first_found;

 44:   if (subset_mult) {
 46:   }
 47:   if (!N && !subset_n) return(0);
 48:   ISGetLocalSize(subset,&n);
 49:   if (subset_mult) {
 50:     ISGetLocalSize(subset_mult,&i);
 51:     if (i != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Local subset and multiplicity sizes don't match! %d != %d",n,i);
 52:   }
 53:   /* create workspace layout for computing global indices of subset */
 54:   ISGetIndices(subset,&idxs);
 55:   lbounds[0] = lbounds[1] = 0;
 56:   for (i=0;i<n;i++) {
 57:     if (idxs[i] < lbounds[0]) lbounds[0] = idxs[i];
 58:     else if (idxs[i] > lbounds[1]) lbounds[1] = idxs[i];
 59:   }
 60:   lbounds[0] = -lbounds[0];
 61:   MPIU_Allreduce(lbounds,gbounds,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)subset));
 62:   gbounds[0] = -gbounds[0];
 63:   N_n  = gbounds[1] - gbounds[0] + 1;

 65:   PetscLayoutCreate(PetscObjectComm((PetscObject)subset),&map);
 66:   PetscLayoutSetBlockSize(map,1);
 67:   PetscLayoutSetSize(map,N_n);
 68:   PetscLayoutSetUp(map);
 69:   PetscLayoutGetLocalSize(map,&Nl);

 71:   /* create sf : leaf_data == multiplicity of indexes, root data == global index in layout */
 72:   PetscMalloc2(n,&leaf_data,Nl,&root_data);
 73:   if (subset_mult) {
 74:     const PetscInt* idxs_mult;

 76:     ISGetIndices(subset_mult,&idxs_mult);
 77:     PetscArraycpy(leaf_data,idxs_mult,n);
 78:     ISRestoreIndices(subset_mult,&idxs_mult);
 79:   } else {
 80:     for (i=0;i<n;i++) leaf_data[i] = 1;
 81:   }
 82:   /* local size of new subset */
 83:   n_n = 0;
 84:   for (i=0;i<n;i++) n_n += leaf_data[i];

 86:   /* global indexes in layout */
 87:   PetscMalloc1(n_n,&gidxs); /* allocating possibly extra space in gidxs which will be used later */
 88:   for (i=0;i<n;i++) gidxs[i] = idxs[i] - gbounds[0];
 89:   ISRestoreIndices(subset,&idxs);
 90:   PetscSFCreate(PetscObjectComm((PetscObject)subset),&sf);
 91:   PetscSFSetGraphLayout(sf,map,n,NULL,PETSC_COPY_VALUES,gidxs);
 92:   PetscLayoutDestroy(&map);

 94:   /* reduce from leaves to roots */
 95:   PetscArrayzero(root_data,Nl);
 96:   PetscSFReduceBegin(sf,MPIU_INT,leaf_data,root_data,MPI_MAX);
 97:   PetscSFReduceEnd(sf,MPIU_INT,leaf_data,root_data,MPI_MAX);

 99:   /* count indexes in local part of layout */
100:   nlocals = 0;
101:   first_index = -1;
102:   first_found = PETSC_FALSE;
103:   for (i=0;i<Nl;i++) {
104:     if (!first_found && root_data[i]) {
105:       first_found = PETSC_TRUE;
106:       first_index = i;
107:     }
108:     nlocals += root_data[i];
109:   }

111:   /* cumulative of number of indexes and size of subset without holes */
112: #if defined(PETSC_HAVE_MPI_EXSCAN)
113:   start = 0;
114:   MPI_Exscan(&nlocals,&start,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)subset));
115: #else
116:   MPI_Scan(&nlocals,&start,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)subset));
117:   start = start-nlocals;
118: #endif

120:   if (N) { /* compute total size of new subset if requested */
121:     *N   = start + nlocals;
122:     MPI_Comm_size(PetscObjectComm((PetscObject)subset),&commsize);
123:     MPI_Bcast(N,1,MPIU_INT,commsize-1,PetscObjectComm((PetscObject)subset));
124:   }

126:   if (!subset_n) {
127:     PetscFree(gidxs);
128:     PetscSFDestroy(&sf);
129:     PetscFree2(leaf_data,root_data);
130:     return(0);
131:   }

133:   /* adapt root data with cumulative */
134:   if (first_found) {
135:     PetscInt old_index;

137:     root_data[first_index] += start;
138:     old_index = first_index;
139:     for (i=first_index+1;i<Nl;i++) {
140:       if (root_data[i]) {
141:         root_data[i] += root_data[old_index];
142:         old_index = i;
143:       }
144:     }
145:   }

147:   /* from roots to leaves */
148:   PetscSFBcastBegin(sf,MPIU_INT,root_data,leaf_data);
149:   PetscSFBcastEnd(sf,MPIU_INT,root_data,leaf_data);
150:   PetscSFDestroy(&sf);

152:   /* create new IS with global indexes without holes */
153:   if (subset_mult) {
154:     const PetscInt* idxs_mult;
155:     PetscInt        cum;

157:     cum = 0;
158:     ISGetIndices(subset_mult,&idxs_mult);
159:     for (i=0;i<n;i++) {
160:       PetscInt j;
161:       for (j=0;j<idxs_mult[i];j++) gidxs[cum++] = leaf_data[i] - idxs_mult[i] + j;
162:     }
163:     ISRestoreIndices(subset_mult,&idxs_mult);
164:   } else {
165:     for (i=0;i<n;i++) {
166:       gidxs[i] = leaf_data[i]-1;
167:     }
168:   }
169:   ISCreateGeneral(PetscObjectComm((PetscObject)subset),n_n,gidxs,PETSC_OWN_POINTER,subset_n);
170:   PetscFree2(leaf_data,root_data);
171:   return(0);
172: }


175: /*@
176:    ISCreateSubIS - Create a sub index set from a global index set selecting some components.

178:    Collective on IS

180:    Input Parmeters:
181: +  is - the index set
182: -  comps - which components we will extract from is

184:    Output Parameters:
185: .  subis - the new sub index set

187:    Level: intermediate

189:    Example usage:
190:    We have an index set (is) living on 3 processes with the following values:
191:    | 4 9 0 | 2 6 7 | 10 11 1|
192:    and another index set (comps) used to indicate which components of is  we want to take,
193:    | 7 5  | 1 2 | 0 4|
194:    The output index set (subis) should look like:
195:    | 11 7 | 9 0 | 4 6|

197: .seealso: VecGetSubVector(), MatCreateSubMatrix()
198: @*/
199: PetscErrorCode ISCreateSubIS(IS is,IS comps,IS *subis)
200: {
201:   PetscSF         sf;
202:   const PetscInt  *is_indices,*comps_indices;
203:   PetscInt        *subis_indices,nroots,nleaves,*mine,i,lidx;
204:   PetscMPIInt     owner;
205:   PetscSFNode     *remote;
206:   PetscErrorCode  ierr;
207:   MPI_Comm        comm;


214:   PetscObjectGetComm((PetscObject)is, &comm);
215:   ISGetLocalSize(comps,&nleaves);
216:   ISGetLocalSize(is,&nroots);
217:   PetscMalloc1(nleaves,&remote);
218:   PetscMalloc1(nleaves,&mine);
219:   ISGetIndices(comps,&comps_indices);
220:   /*
221:    * Construct a PetscSF in which "is" data serves as roots and "subis" is leaves.
222:    * Root data are sent to leaves using PetscSFBcast().
223:    * */
224:   for (i=0; i<nleaves; i++) {
225:     mine[i] = i;
226:     /* Connect a remote root with the current leaf. The value on the remote root
227:      * will be received by the current local leaf.
228:      * */
229:     owner = -1;
230:     lidx =  -1;
231:     PetscLayoutFindOwnerIndex(is->map,comps_indices[i],&owner,&lidx);
232:     remote[i].rank = owner;
233:     remote[i].index = lidx;
234:   }
235:   ISRestoreIndices(comps,&comps_indices);
236:   PetscSFCreate(comm,&sf);
237:   PetscSFSetFromOptions(sf);\
238:   PetscSFSetGraph(sf,nroots,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);

240:   PetscMalloc1(nleaves,&subis_indices);
241:   ISGetIndices(is, &is_indices);
242:   PetscSFBcastBegin(sf,MPIU_INT,is_indices,subis_indices);
243:   PetscSFBcastEnd(sf,MPIU_INT,is_indices,subis_indices);
244:   ISRestoreIndices(is,&is_indices);
245:   PetscSFDestroy(&sf);
246:   ISCreateGeneral(comm,nleaves,subis_indices,PETSC_OWN_POINTER,subis);
247:   return(0);
248: }

250: /*@
251:    ISClearInfoCache - clear the cache of computed index set properties

253:    Not collective

255:    Input Parameters:
256: +  is - the index set
257: -  clear_permanent_local - whether to remove the permanent status of local properties

259:    NOTE: because all processes must agree on the global permanent status of a property,
260:    the permanent status can only be changed with ISSetInfo(), because this routine is not collective

262:    Level: developer

264: .seealso:  ISInfo, ISInfoType, ISSetInfo(), ISClearInfoCache()

266: @*/
267: PetscErrorCode ISClearInfoCache(IS is, PetscBool clear_permanent_local)
268: {
269:   PetscInt i, j;

274:   for (i = 0; i < IS_INFO_MAX; i++) {
275:     if (clear_permanent_local) is->info_permanent[IS_LOCAL][i] = PETSC_FALSE;
276:     for (j = 0; j < 2; j++) {
277:       if (!is->info_permanent[j][i]) is->info[j][i] = IS_INFO_UNKNOWN;
278:     }
279:   }
280:   return(0);
281: }

283: static PetscErrorCode ISSetInfo_Internal(IS is, ISInfo info, ISInfoType type, ISInfoBool ipermanent, PetscBool flg)
284: {
285:   ISInfoBool     iflg = flg ? IS_INFO_TRUE : IS_INFO_FALSE;
286:   PetscInt       itype = (type == IS_LOCAL) ? 0 : 1;
287:   PetscBool      permanent_set = (ipermanent == IS_INFO_UNKNOWN) ? PETSC_FALSE : PETSC_TRUE;
288:   PetscBool      permanent = (ipermanent == IS_INFO_TRUE) ? PETSC_TRUE : PETSC_FALSE;

291:   /* set this property */
292:   is->info[itype][(int)info] = iflg;
293:   if (permanent_set) is->info_permanent[itype][(int)info] = permanent;
294:   /* set implications */
295:   switch (info) {
296:   case IS_SORTED:
297:     if (flg && type == IS_GLOBAL) { /* an array that is globally sorted is also locally sorted */
298:       is->info[IS_LOCAL][(int)info] = IS_INFO_TRUE;
299:       /* global permanence implies local permanence */
300:       if (permanent_set && permanent) is->info_permanent[IS_LOCAL][(int)info] = PETSC_TRUE;
301:     }
302:     if (!flg) { /* if an array is not sorted, it cannot be an interval or the identity */
303:       is->info[itype][IS_INTERVAL] = IS_INFO_FALSE;
304:       is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
305:       if (permanent_set) {
306:         is->info_permanent[itype][IS_INTERVAL] = permanent;
307:         is->info_permanent[itype][IS_IDENTITY] = permanent;
308:       }
309:     }
310:     break;
311:   case IS_UNIQUE:
312:     if (flg && type == IS_GLOBAL) { /* an array that is globally unique is also locally unique */
313:       is->info[IS_LOCAL][(int)info] = IS_INFO_TRUE;
314:       /* global permanence implies local permanence */
315:       if (permanent_set && permanent) is->info_permanent[IS_LOCAL][(int)info] = PETSC_TRUE;
316:     }
317:     if (!flg) { /* if an array is not unique, it cannot be a permutation, and interval, or the identity */
318:       is->info[itype][IS_PERMUTATION] = IS_INFO_FALSE;
319:       is->info[itype][IS_INTERVAL]    = IS_INFO_FALSE;
320:       is->info[itype][IS_IDENTITY]    = IS_INFO_FALSE;
321:       if (permanent_set) {
322:         is->info_permanent[itype][IS_PERMUTATION] = permanent;
323:         is->info_permanent[itype][IS_INTERVAL]    = permanent;
324:         is->info_permanent[itype][IS_IDENTITY]    = permanent;
325:       }
326:     }
327:     break;
328:   case IS_PERMUTATION:
329:     if (flg) { /* an array that is a permutation is unique and is unique locally */
330:       is->info[itype][IS_UNIQUE] = IS_INFO_TRUE;
331:       is->info[IS_LOCAL][IS_UNIQUE] = IS_INFO_TRUE;
332:       if (permanent_set && permanent) {
333:         is->info_permanent[itype][IS_UNIQUE] = PETSC_TRUE;
334:         is->info_permanent[IS_LOCAL][IS_UNIQUE] = PETSC_TRUE;
335:       }
336:     } else { /* an array that is not a permutation cannot be the identity */
337:       is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
338:       if (permanent_set) is->info_permanent[itype][IS_IDENTITY] = permanent;
339:     }
340:     break;
341:   case IS_INTERVAL:
342:     if (flg) { /* an array that is an interval is sorted and unique */
343:       is->info[itype][IS_SORTED]         = IS_INFO_TRUE;
344:       is->info[IS_LOCAL][IS_SORTED]      = IS_INFO_TRUE;
345:       is->info[itype][IS_UNIQUE]         = IS_INFO_TRUE;
346:       is->info[IS_LOCAL][IS_UNIQUE]      = IS_INFO_TRUE;
347:       if (permanent_set && permanent) {
348:         is->info_permanent[itype][IS_SORTED]    = PETSC_TRUE;
349:         is->info_permanent[IS_LOCAL][IS_SORTED] = PETSC_TRUE;
350:         is->info_permanent[itype][IS_UNIQUE]    = PETSC_TRUE;
351:         is->info_permanent[IS_LOCAL][IS_UNIQUE] = PETSC_TRUE;
352:       }
353:     } else { /* an array that is not an interval cannot be the identity */
354:       is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
355:       if (permanent_set) is->info_permanent[itype][IS_IDENTITY] = permanent;
356:     }
357:     break;
358:   case IS_IDENTITY:
359:     if (flg) { /* an array that is the identity is sorted, unique, an interval, and a permutation */
360:       is->info[itype][IS_SORTED]         = IS_INFO_TRUE;
361:       is->info[IS_LOCAL][IS_SORTED]      = IS_INFO_TRUE;
362:       is->info[itype][IS_UNIQUE]         = IS_INFO_TRUE;
363:       is->info[IS_LOCAL][IS_UNIQUE]      = IS_INFO_TRUE;
364:       is->info[itype][IS_PERMUTATION]    = IS_INFO_TRUE;
365:       is->info[itype][IS_INTERVAL]       = IS_INFO_TRUE;
366:       is->info[IS_LOCAL][IS_INTERVAL]    = IS_INFO_TRUE;
367:       if (permanent_set && permanent) {
368:         is->info_permanent[itype][IS_SORTED]         = PETSC_TRUE;
369:         is->info_permanent[IS_LOCAL][IS_SORTED]      = PETSC_TRUE;
370:         is->info_permanent[itype][IS_UNIQUE]         = PETSC_TRUE;
371:         is->info_permanent[IS_LOCAL][IS_UNIQUE]      = PETSC_TRUE;
372:         is->info_permanent[itype][IS_PERMUTATION]    = PETSC_TRUE;
373:         is->info_permanent[itype][IS_INTERVAL]       = PETSC_TRUE;
374:         is->info_permanent[IS_LOCAL][IS_INTERVAL]    = PETSC_TRUE;
375:       }
376:     }
377:     break;
378:   default:
379:     if (type == IS_LOCAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown IS property");
380:     else SETERRQ(PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_OUTOFRANGE, "Unknown IS property");
381:   }
382:   return(0);
383: }

385: /*@
386:    ISSetInfo - Set known information about an index set.

388:    Logically Collective on IS if type is IS_GLOBAL

390:    Input Parameters:
391: +  is - the index set
392: .  info - describing a property of the index set, one of those listed below,
393: .  type - IS_LOCAL if the information describes the local portion of the index set,
394:           IS_GLOBAL if it describes the whole index set
395: .  permanent - PETSC_TRUE if it is known that the property will persist through changes to the index set, PETSC_FALSE otherwise
396:                If the user sets a property as permanently known, it will bypass computation of that property
397: -  flg - set the described property as true (PETSC_TRUE) or false (PETSC_FALSE)

399:   Info Describing IS Structure:
400: +    IS_SORTED - the [local part of the] index set is sorted in ascending order
401: .    IS_UNIQUE - each entry in the [local part of the] index set is unique
402: .    IS_PERMUTATION - the [local part of the] index set is a permutation of the integers {0, 1, ..., N-1}, where N is the size of the [local part of the] index set
403: .    IS_INTERVAL - the [local part of the] index set is equal to a contiguous range of integers {f, f + 1, ..., f + N-1}
404: -    IS_IDENTITY - the [local part of the] index set is equal to the integers {0, 1, ..., N-1}


407:    Notes:
408:    If type is IS_GLOBAL, all processes that share the index set must pass the same value in flg

410:    It is possible to set a property with ISSetInfo() that contradicts what would be previously computed with ISGetInfo()

412:    Level: advanced

414: .seealso:  ISInfo, ISInfoType, IS

416: @*/
417: PetscErrorCode ISSetInfo(IS is, ISInfo info, ISInfoType type, PetscBool permanent, PetscBool flg)
418: {
419:   MPI_Comm       comm, errcomm;
420:   PetscMPIInt    size;

426:   comm = PetscObjectComm((PetscObject)is);
427:   if (type == IS_GLOBAL) {
431:     errcomm = comm;
432:   } else {
433:     errcomm = PETSC_COMM_SELF;
434:   }

436:   if (((int) info) <= IS_INFO_MIN || ((int) info) >= IS_INFO_MAX) SETERRQ1(errcomm,PETSC_ERR_ARG_OUTOFRANGE,"Options %d is out of range",(int)info);

438:   MPI_Comm_size(comm, &size);
439:   /* do not use global values if size == 1: it makes it easier to keep the implications straight */
440:   if (size == 1) type = IS_LOCAL;
441:   ISSetInfo_Internal(is, info, type, permanent ? IS_INFO_TRUE : IS_INFO_FALSE, flg);
442:   return(0);
443: }

445: static PetscErrorCode ISGetInfo_Sorted(IS is, ISInfoType type, PetscBool *flg)
446: {
447:   MPI_Comm       comm;
448:   PetscMPIInt    size, rank;

452:   comm = PetscObjectComm((PetscObject)is);
453:   MPI_Comm_size(comm, &size);
454:   MPI_Comm_size(comm, &rank);
455:   if (type == IS_GLOBAL && is->ops->sortedglobal) {
456:     (*is->ops->sortedglobal)(is,flg);
457:   } else {
458:     PetscBool sortedLocal = PETSC_FALSE;

460:     /* determine if the array is locally sorted */
461:     if (type == IS_GLOBAL && size > 1) {
462:       /* call ISGetInfo so that a cached value will be used if possible */
463:       ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, &sortedLocal);
464:     } else if (is->ops->sortedlocal) {
465:       (*is->ops->sortedlocal)(is,&sortedLocal);
466:     } else {
467:       /* default: get the local indices and directly check */
468:       const PetscInt *idx;
469:       PetscInt n;

471:       ISGetIndices(is, &idx);
472:       ISGetLocalSize(is, &n);
473:       PetscSortedInt(n, idx, &sortedLocal);
474:       ISRestoreIndices(is, &idx);
475:     }

477:     if (type == IS_LOCAL || size == 1) {
478:       *flg = sortedLocal;
479:     } else {
480:       MPI_Allreduce(&sortedLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
481:       if (*flg) {
482:         PetscInt  n, min = PETSC_MAX_INT, max = PETSC_MIN_INT;
483:         PetscInt  maxprev;

485:         ISGetLocalSize(is, &n);
486:         if (n) {ISGetMinMax(is, &min, &max);}
487:         maxprev = PETSC_MIN_INT;
488:         MPI_Exscan(&max, &maxprev, 1, MPIU_INT, MPI_MAX, comm);
489:         if (rank && (maxprev > min)) sortedLocal = PETSC_FALSE;
490:         MPI_Allreduce(&sortedLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
491:       }
492:     }
493:   }
494:   return(0);
495: }

497: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[]);

499: static PetscErrorCode ISGetInfo_Unique(IS is, ISInfoType type, PetscBool *flg)
500: {
501:   MPI_Comm       comm;
502:   PetscMPIInt    size, rank;
503:   PetscInt       i;

507:   comm = PetscObjectComm((PetscObject)is);
508:   MPI_Comm_size(comm, &size);
509:   MPI_Comm_size(comm, &rank);
510:   if (type == IS_GLOBAL && is->ops->uniqueglobal) {
511:     (*is->ops->uniqueglobal)(is,flg);
512:   } else {
513:     PetscBool uniqueLocal;
514:     PetscInt  n = -1;
515:     PetscInt  *idx = NULL;

517:     /* determine if the array is locally unique */
518:     if (type == IS_GLOBAL && size > 1) {
519:       /* call ISGetInfo so that a cached value will be used if possible */
520:       ISGetInfo(is, IS_UNIQUE, IS_LOCAL, PETSC_TRUE, &uniqueLocal);
521:     } else if (is->ops->uniquelocal) {
522:       (*is->ops->uniquelocal)(is,&uniqueLocal);
523:     } else {
524:       /* default: get the local indices and directly check */
525:       uniqueLocal = PETSC_TRUE;
526:       ISGetLocalSize(is, &n);
527:       PetscMalloc1(n, &idx);
528:       ISGetIndicesCopy(is, idx);
529:       PetscSortInt(n, idx);
530:       for (i = 1; i < n; i++) if (idx[i] == idx[i-1]) break;
531:       if (i < n) uniqueLocal = PETSC_FALSE;
532:     }

534:     PetscFree(idx);
535:     if (type == IS_LOCAL || size == 1) {
536:       *flg = uniqueLocal;
537:     } else {
538:       MPI_Allreduce(&uniqueLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
539:       if (*flg) {
540:         PetscInt  min = PETSC_MAX_INT, max = PETSC_MIN_INT, maxprev;

542:         if (!idx) {
543:           ISGetLocalSize(is, &n);
544:           PetscMalloc1(n, &idx);
545:           ISGetIndicesCopy(is, idx);
546:         }
547:         PetscParallelSortInt(is->map, is->map, idx, idx);
548:         if (n) {
549:           min = idx[0];
550:           max = idx[n - 1];
551:         }
552:         for (i = 1; i < n; i++) if (idx[i] == idx[i-1]) break;
553:         if (i < n) uniqueLocal = PETSC_FALSE;
554:         maxprev = PETSC_MIN_INT;
555:         MPI_Exscan(&max, &maxprev, 1, MPIU_INT, MPI_MAX, comm);
556:         if (rank && (maxprev == min)) uniqueLocal = PETSC_FALSE;
557:         MPI_Allreduce(&uniqueLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
558:       }
559:     }
560:     PetscFree(idx);
561:   }
562:   return(0);
563: }

565: static PetscErrorCode ISGetInfo_Permutation(IS is, ISInfoType type, PetscBool *flg)
566: {
567:   MPI_Comm       comm;
568:   PetscMPIInt    size, rank;

572:   comm = PetscObjectComm((PetscObject)is);
573:   MPI_Comm_size(comm, &size);
574:   MPI_Comm_size(comm, &rank);
575:   if (type == IS_GLOBAL && is->ops->permglobal) {
576:     (*is->ops->permglobal)(is,flg);
577:   } else if (type == IS_LOCAL && is->ops->permlocal) {
578:     (*is->ops->permlocal)(is,flg);
579:   } else {
580:     PetscBool permLocal;
581:     PetscInt  n, i, rStart;
582:     PetscInt  *idx;

584:     ISGetLocalSize(is, &n);
585:     PetscMalloc1(n, &idx);
586:     ISGetIndicesCopy(is, idx);
587:     if (type == IS_GLOBAL) {
588:       PetscParallelSortInt(is->map, is->map, idx, idx);
589:       PetscLayoutGetRange(is->map, &rStart, NULL);
590:     } else {
591:       PetscSortInt(n, idx);
592:       rStart = 0;
593:     }
594:     permLocal = PETSC_TRUE;
595:     for (i = 0; i < n; i++) {
596:       if (idx[i] != rStart + i) break;
597:     }
598:     if (i < n) permLocal = PETSC_FALSE;
599:     if (type == IS_LOCAL || size == 1) {
600:       *flg = permLocal;
601:     } else {
602:       MPI_Allreduce(&permLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
603:     }
604:     PetscFree(idx);
605:   }
606:   return(0);
607: }

609: static PetscErrorCode ISGetInfo_Interval(IS is, ISInfoType type, PetscBool *flg)
610: {
611:   MPI_Comm       comm;
612:   PetscMPIInt    size, rank;
613:   PetscInt       i;

617:   comm = PetscObjectComm((PetscObject)is);
618:   MPI_Comm_size(comm, &size);
619:   MPI_Comm_size(comm, &rank);
620:   if (type == IS_GLOBAL && is->ops->intervalglobal) {
621:     (*is->ops->intervalglobal)(is,flg);
622:   } else {
623:     PetscBool intervalLocal;

625:     /* determine if the array is locally an interval */
626:     if (type == IS_GLOBAL && size > 1) {
627:       /* call ISGetInfo so that a cached value will be used if possible */
628:       ISGetInfo(is, IS_INTERVAL, IS_LOCAL, PETSC_TRUE, &intervalLocal);
629:     } else if (is->ops->intervallocal) {
630:       (*is->ops->intervallocal)(is,&intervalLocal);
631:     } else {
632:       PetscInt        n;
633:       const PetscInt  *idx;
634:       /* default: get the local indices and directly check */
635:       intervalLocal = PETSC_TRUE;
636:       ISGetLocalSize(is, &n);
637:       PetscMalloc1(n, &idx);
638:       ISGetIndices(is, &idx);
639:       for (i = 1; i < n; i++) if (idx[i] != idx[i-1] + 1) break;
640:       if (i < n) intervalLocal = PETSC_FALSE;
641:       ISRestoreIndices(is, &idx);
642:     }

644:     if (type == IS_LOCAL || size == 1) {
645:       *flg = intervalLocal;
646:     } else {
647:       MPI_Allreduce(&intervalLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
648:       if (*flg) {
649:         PetscInt  n, min = PETSC_MAX_INT, max = PETSC_MIN_INT;
650:         PetscInt  maxprev;

652:         ISGetLocalSize(is, &n);
653:         if (n) {ISGetMinMax(is, &min, &max);}
654:         maxprev = PETSC_MIN_INT;
655:         MPI_Exscan(&max, &maxprev, 1, MPIU_INT, MPI_MAX, comm);
656:         if (rank && n && (maxprev != min - 1)) intervalLocal = PETSC_FALSE;
657:         MPI_Allreduce(&intervalLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
658:       }
659:     }
660:   }
661:   return(0);
662: }

664: static PetscErrorCode ISGetInfo_Identity(IS is, ISInfoType type, PetscBool *flg)
665: {
666:   MPI_Comm       comm;
667:   PetscMPIInt    size, rank;

671:   comm = PetscObjectComm((PetscObject)is);
672:   MPI_Comm_size(comm, &size);
673:   MPI_Comm_size(comm, &rank);
674:   if (type == IS_GLOBAL && is->ops->intervalglobal) {
675:     PetscBool isinterval;

677:     (*is->ops->intervalglobal)(is,&isinterval);
678:     *flg = PETSC_FALSE;
679:     if (isinterval) {
680:       PetscInt  min;

682:       ISGetMinMax(is, &min, NULL);
683:       MPI_Bcast(&min, 1, MPIU_INT, 0, comm);
684:       if (min == 0) *flg = PETSC_TRUE;
685:     }
686:   } else if (type == IS_LOCAL && is->ops->intervallocal) {
687:     PetscBool isinterval;

689:     (*is->ops->intervallocal)(is,&isinterval);
690:     *flg = PETSC_FALSE;
691:     if (isinterval) {
692:       PetscInt  min;

694:       ISGetMinMax(is, &min, NULL);
695:       if (min == 0) *flg = PETSC_TRUE;
696:     }
697:   } else {
698:     PetscBool identLocal;
699:     PetscInt  n, i, rStart;
700:     const PetscInt *idx;

702:     ISGetLocalSize(is, &n);
703:     ISGetIndices(is, &idx);
704:     PetscLayoutGetRange(is->map, &rStart, NULL);
705:     identLocal = PETSC_TRUE;
706:     for (i = 0; i < n; i++) {
707:       if (idx[i] != rStart + i) break;
708:     }
709:     if (i < n) identLocal = PETSC_FALSE;
710:     if (type == IS_LOCAL || size == 1) {
711:       *flg = identLocal;
712:     } else {
713:       MPI_Allreduce(&identLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm);
714:     }
715:     ISRestoreIndices(is, &idx);
716:   }
717:   return(0);
718: }

720: /*@
721:    ISGetInfo - Determine whether an index set satisfies a given property

723:    Collective or logically collective on IS if the type is IS_GLOBAL (logically collective if the value of the property has been permanently set with ISSetInfo())

725:    Input Parameters:
726: +  is - the index set
727: .  info - describing a property of the index set, one of those listed in the documentation of ISSetInfo()
728: .  compute - if PETSC_FALSE, the property will not be computed if it is not already known and the property will be assumed to be false
729: -  type - whether the property is local (IS_LOCAL) or global (IS_GLOBAL)

731:    Output Parameter:
732: .  flg - wheter the property is true (PETSC_TRUE) or false (PETSC_FALSE)

734:    Note: ISGetInfo uses cached values when possible, which will be incorrect if ISSetInfo() has been called with incorrect information.  To clear cached values, use ISClearInfoCache().

736:    Level: advanced

738: .seealso:  ISInfo, ISInfoType, ISSetInfo(), ISClearInfoCache()

740: @*/
741: PetscErrorCode ISGetInfo(IS is, ISInfo info, ISInfoType type, PetscBool compute, PetscBool *flg)
742: {
743:   MPI_Comm       comm, errcomm;
744:   PetscMPIInt    rank, size;
745:   PetscInt       itype;
746:   PetscBool      hasprop;
747:   PetscBool      infer;

753:   comm = PetscObjectComm((PetscObject)is);
754:   if (type == IS_GLOBAL) {
756:     errcomm = comm;
757:   } else {
758:     errcomm = PETSC_COMM_SELF;
759:   }

761:   MPI_Comm_size(comm, &size);
762:   MPI_Comm_rank(comm, &rank);

764:   if (((int) info) <= IS_INFO_MIN || ((int) info) >= IS_INFO_MAX) SETERRQ1(errcomm,PETSC_ERR_ARG_OUTOFRANGE,"Options %d is out of range",(int)info);
765:   if (size == 1) type = IS_LOCAL;
766:   itype = (type == IS_LOCAL) ? 0 : 1;
767:   hasprop = PETSC_FALSE;
768:   infer = PETSC_FALSE;
769:   if (is->info_permanent[itype][(int)info]) {
770:     hasprop = (is->info[itype][(int)info] == IS_INFO_TRUE) ? PETSC_TRUE : PETSC_FALSE;
771:     infer = PETSC_TRUE;
772:   } else if ((itype == IS_LOCAL) && (is->info[IS_LOCAL][info] != IS_INFO_UNKNOWN)) {
773:     /* we can cache local properties as long as we clear them when the IS changes */
774:     /* NOTE: we only cache local values because there is no ISAssemblyBegin()/ISAssemblyEnd(),
775:      so we have no way of knowing when a cached value has been invalidated by changes on a different process */
776:     hasprop = (is->info[itype][(int)info] == IS_INFO_TRUE) ? PETSC_TRUE : PETSC_FALSE;
777:     infer = PETSC_TRUE;
778:   } else if (compute) {
779:     switch (info) {
780:     case IS_SORTED:
781:       ISGetInfo_Sorted(is, type, &hasprop);
782:       break;
783:     case IS_UNIQUE:
784:       ISGetInfo_Unique(is, type, &hasprop);
785:       break;
786:     case IS_PERMUTATION:
787:       ISGetInfo_Permutation(is, type, &hasprop);
788:       break;
789:     case IS_INTERVAL:
790:       ISGetInfo_Interval(is, type, &hasprop);
791:       break;
792:     case IS_IDENTITY:
793:       ISGetInfo_Identity(is, type, &hasprop);
794:       break;
795:     default:
796:       SETERRQ(errcomm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown IS property");
797:     }
798:     infer = PETSC_TRUE;
799:   }
800:   /* call ISSetInfo_Internal to keep all of the implications straight */
801:   if (infer) {ISSetInfo_Internal(is, info, type, IS_INFO_UNKNOWN, hasprop);}
802:   *flg = hasprop;
803:   return(0);
804: }

806: static PetscErrorCode ISCopyInfo(IS source, IS dest)
807: {

811:   PetscArraycpy(&dest->info[0], &source->info[0], 2);
812:   PetscArraycpy(&dest->info_permanent[0], &source->info_permanent[0], 2);
813:   return(0);
814: }

816: /*@
817:    ISIdentity - Determines whether index set is the identity mapping.

819:    Collective on IS

821:    Input Parmeters:
822: .  is - the index set

824:    Output Parameters:
825: .  ident - PETSC_TRUE if an identity, else PETSC_FALSE

827:    Level: intermediate

829:    Note: If ISSetIdentity() (or ISSetInfo() for a permanent property) has been called,
830:    ISIdentity() will return its answer without communication between processes, but
831:    otherwise the output ident will be computed from ISGetInfo(),
832:    which may require synchronization on the communicator of IS.  To avoid this computation,
833:    call ISGetInfo() directly with the compute flag set to PETSC_FALSE, and ident will be assumed false.

835: .seealso: ISSetIdentity(), ISGetInfo()
836: @*/
837: PetscErrorCode  ISIdentity(IS is,PetscBool  *ident)
838: {

844:   ISGetInfo(is,IS_IDENTITY,IS_GLOBAL,PETSC_TRUE,ident);
845:   return(0);
846: }

848: /*@
849:    ISSetIdentity - Informs the index set that it is an identity.

851:    Logically Collective on IS

853:    Input Parmeters:
854: .  is - the index set

856:    Level: intermediate

858:    Note: The IS will be considered the identity permanently, even if indices have been changes (for example, with
859:    ISGeneralSetIndices()).  It's a good idea to only set this property if the IS will not change in the future.
860:    To clear this property, use ISClearInfoCache().

862: .seealso: ISIdentity(), ISSetInfo(), ISClearInfoCache()
863: @*/
864: PetscErrorCode  ISSetIdentity(IS is)
865: {

870:   ISSetInfo(is,IS_IDENTITY,IS_GLOBAL,PETSC_TRUE,PETSC_TRUE);
871:   return(0);
872: }

874: /*@
875:    ISContiguousLocal - Locates an index set with contiguous range within a global range, if possible

877:    Not Collective

879:    Input Parmeters:
880: +  is - the index set
881: .  gstart - global start
882: -  gend - global end

884:    Output Parameters:
885: +  start - start of contiguous block, as an offset from gstart
886: -  contig - PETSC_TRUE if the index set refers to contiguous entries on this process, else PETSC_FALSE

888:    Level: developer

890: .seealso: ISGetLocalSize(), VecGetOwnershipRange()
891: @*/
892: PetscErrorCode  ISContiguousLocal(IS is,PetscInt gstart,PetscInt gend,PetscInt *start,PetscBool *contig)
893: {

900:   if (is->ops->contiguous) {
901:     (*is->ops->contiguous)(is,gstart,gend,start,contig);
902:   } else {
903:     *start  = -1;
904:     *contig = PETSC_FALSE;
905:   }
906:   return(0);
907: }

909: /*@
910:    ISPermutation - PETSC_TRUE or PETSC_FALSE depending on whether the
911:    index set has been declared to be a permutation.

913:    Logically Collective on IS

915:    Input Parmeters:
916: .  is - the index set

918:    Output Parameters:
919: .  perm - PETSC_TRUE if a permutation, else PETSC_FALSE

921:    Level: intermediate

923:    Note: If it is not alread known that the IS is a permutation (if ISSetPermutation()
924:    or ISSetInfo() has not been called), this routine will not attempt to compute
925:    whether the index set is a permutation and will assume perm is PETSC_FALSE.
926:    To compute the value when it is not already known, use ISGetInfo() with
927:    the compute flag set to PETSC_TRUE.

929: .seealso: ISSetPermutation(), ISGetInfo()
930: @*/
931: PetscErrorCode  ISPermutation(IS is,PetscBool  *perm)
932: {

938:   ISGetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_FALSE,perm);
939:   return(0);
940: }

942: /*@
943:    ISSetPermutation - Informs the index set that it is a permutation.

945:    Logically Collective on IS

947:    Input Parmeters:
948: .  is - the index set

950:    Level: intermediate


953:    The debug version of the libraries (./configure --with-debugging=1) checks if the
954:   index set is actually a permutation. The optimized version just believes you.

956:    Note: The IS will be considered a permutation permanently, even if indices have been changes (for example, with
957:    ISGeneralSetIndices()).  It's a good idea to only set this property if the IS will not change in the future.
958:    To clear this property, use ISClearInfoCache().

960: .seealso: ISPermutation(), ISSetInfo(), ISClearInfoCache().
961: @*/
962: PetscErrorCode  ISSetPermutation(IS is)
963: {

968: #if defined(PETSC_USE_DEBUG)
969:   {
970:     PetscMPIInt    size;

972:     MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
973:     if (size == 1) {
974:       PetscInt       i,n,*idx;
975:       const PetscInt *iidx;

977:       ISGetSize(is,&n);
978:       PetscMalloc1(n,&idx);
979:       ISGetIndices(is,&iidx);
980:       PetscArraycpy(idx,iidx,n);
981:       PetscSortInt(n,idx);
982:       for (i=0; i<n; i++) {
983:         if (idx[i] != i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Index set is not a permutation");
984:       }
985:       PetscFree(idx);
986:       ISRestoreIndices(is,&iidx);
987:     }
988:   }
989: #endif
990:   ISSetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_TRUE,PETSC_TRUE);
991:   return(0);
992: }

994: /*@
995:    ISDestroy - Destroys an index set.

997:    Collective on IS

999:    Input Parameters:
1000: .  is - the index set

1002:    Level: beginner

1004: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlocked()
1005: @*/
1006: PetscErrorCode  ISDestroy(IS *is)
1007: {

1011:   if (!*is) return(0);
1013:   if (--((PetscObject)(*is))->refct > 0) {*is = 0; return(0);}
1014:   if ((*is)->complement) {
1015:     PetscInt refcnt;
1016:     PetscObjectGetReference((PetscObject)((*is)->complement), &refcnt);
1017:     if (refcnt > 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlocal IS has not been restored");
1018:     ISDestroy(&(*is)->complement);
1019:   }
1020:   if ((*is)->ops->destroy) {
1021:     (*(*is)->ops->destroy)(*is);
1022:   }
1023:   PetscLayoutDestroy(&(*is)->map);
1024:   /* Destroy local representations of offproc data. */
1025:   PetscFree((*is)->total);
1026:   PetscFree((*is)->nonlocal);
1027:   PetscHeaderDestroy(is);
1028:   return(0);
1029: }

1031: /*@
1032:    ISInvertPermutation - Creates a new permutation that is the inverse of
1033:                          a given permutation.

1035:    Collective on IS

1037:    Input Parameter:
1038: +  is - the index set
1039: -  nlocal - number of indices on this processor in result (ignored for 1 proccessor) or
1040:             use PETSC_DECIDE

1042:    Output Parameter:
1043: .  isout - the inverse permutation

1045:    Level: intermediate

1047:    Notes:
1048:     For parallel index sets this does the complete parallel permutation, but the
1049:     code is not efficient for huge index sets (10,000,000 indices).

1051: @*/
1052: PetscErrorCode  ISInvertPermutation(IS is,PetscInt nlocal,IS *isout)
1053: {
1054:   PetscBool      isperm, isidentity, issame;

1060:   ISGetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_TRUE,&isperm);
1061:   if (!isperm) SETERRQ(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_WRONG,"Not a permutation");
1062:   ISGetInfo(is,IS_IDENTITY,IS_GLOBAL,PETSC_TRUE,&isidentity);
1063:   issame = PETSC_FALSE;
1064:   if (isidentity) {
1065:     PetscInt n;
1066:     PetscBool isallsame;

1068:     ISGetLocalSize(is, &n);
1069:     issame = (PetscBool) (n == nlocal);
1070:     MPI_Allreduce(&issame, &isallsame, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)is));
1071:     issame = isallsame;
1072:   }
1073:   if (issame) {
1074:     ISDuplicate(is,isout);
1075:   } else {
1076:     (*is->ops->invertpermutation)(is,nlocal,isout);
1077:     ISSetPermutation(*isout);
1078:   }
1079:   return(0);
1080: }

1082: /*@
1083:    ISGetSize - Returns the global length of an index set.

1085:    Not Collective

1087:    Input Parameter:
1088: .  is - the index set

1090:    Output Parameter:
1091: .  size - the global size

1093:    Level: beginner


1096: @*/
1097: PetscErrorCode  ISGetSize(IS is,PetscInt *size)
1098: {
1102:   *size = is->map->N;
1103:   return(0);
1104: }

1106: /*@
1107:    ISGetLocalSize - Returns the local (processor) length of an index set.

1109:    Not Collective

1111:    Input Parameter:
1112: .  is - the index set

1114:    Output Parameter:
1115: .  size - the local size

1117:    Level: beginner

1119: @*/
1120: PetscErrorCode  ISGetLocalSize(IS is,PetscInt *size)
1121: {
1125:   *size = is->map->n;
1126:   return(0);
1127: }

1129: /*@
1130:    ISGetLayout - get PetscLayout describing index set layout

1132:    Not Collective

1134:    Input Arguments:
1135: .  is - the index set

1137:    Output Arguments:
1138: .  map - the layout

1140:    Level: developer

1142: .seealso: ISGetSize(), ISGetLocalSize()
1143: @*/
1144: PetscErrorCode ISGetLayout(IS is,PetscLayout *map)
1145: {

1149:   *map = is->map;
1150:   return(0);
1151: }

1153: /*@C
1154:    ISGetIndices - Returns a pointer to the indices.  The user should call
1155:    ISRestoreIndices() after having looked at the indices.  The user should
1156:    NOT change the indices.

1158:    Not Collective

1160:    Input Parameter:
1161: .  is - the index set

1163:    Output Parameter:
1164: .  ptr - the location to put the pointer to the indices

1166:    Fortran Note:
1167:    This routine has two different interfaces from Fortran; the first is not recommend, it does not require Fortran 90
1168: $    IS          is
1169: $    integer     is_array(1)
1170: $    PetscOffset i_is
1171: $    int         ierr
1172: $       call ISGetIndices(is,is_array,i_is,ierr)
1173: $
1174: $   Access first local entry in list
1175: $      value = is_array(i_is + 1)
1176: $
1177: $      ...... other code
1178: $       call ISRestoreIndices(is,is_array,i_is,ierr)
1179:    The second Fortran interface is recommended.
1180: $          use petscisdef
1181: $          PetscInt, pointer :: array(:)
1182: $          PetscErrorCode  ierr
1183: $          IS       i
1184: $          call ISGetIndicesF90(i,array,ierr)



1188:    See the Fortran chapter of the users manual and
1189:    petsc/src/is/examples/[tutorials,tests] for details.

1191:    Level: intermediate


1194: .seealso: ISRestoreIndices(), ISGetIndicesF90()
1195: @*/
1196: PetscErrorCode  ISGetIndices(IS is,const PetscInt *ptr[])
1197: {

1203:   (*is->ops->getindices)(is,ptr);
1204:   return(0);
1205: }

1207: /*@C
1208:    ISGetMinMax - Gets the minimum and maximum values in an IS

1210:    Not Collective

1212:    Input Parameter:
1213: .  is - the index set

1215:    Output Parameter:
1216: +   min - the minimum value
1217: -   max - the maximum value

1219:    Level: intermediate

1221:    Notes:
1222:     Empty index sets return min=PETSC_MAX_INT and max=PETSC_MIN_INT.
1223:     In parallel, it returns the min and max of the local portion of the IS


1226: .seealso: ISGetIndices(), ISRestoreIndices(), ISGetIndicesF90()
1227: @*/
1228: PetscErrorCode  ISGetMinMax(IS is,PetscInt *min,PetscInt *max)
1229: {
1232:   if (min) *min = is->min;
1233:   if (max) *max = is->max;
1234:   return(0);
1235: }

1237: /*@
1238:   ISLocate - determine the location of an index within the local component of an index set

1240:   Not Collective

1242:   Input Parameter:
1243: + is - the index set
1244: - key - the search key

1246:   Output Parameter:
1247: . location - if >= 0, a location within the index set that is equal to the key, otherwise the key is not in the index set

1249:   Level: intermediate
1250: @*/
1251: PetscErrorCode ISLocate(IS is, PetscInt key, PetscInt *location)
1252: {

1256:   if (is->ops->locate) {
1257:     (*is->ops->locate)(is,key,location);
1258:   } else {
1259:     PetscInt       numIdx;
1260:     PetscBool      sorted;
1261:     const PetscInt *idx;

1263:     ISGetLocalSize(is,&numIdx);
1264:     ISGetIndices(is,&idx);
1265:     ISSorted(is,&sorted);
1266:     if (sorted) {
1267:       PetscFindInt(key,numIdx,idx,location);
1268:     } else {
1269:       PetscInt i;

1271:       *location = -1;
1272:       for (i = 0; i < numIdx; i++) {
1273:         if (idx[i] == key) {
1274:           *location = i;
1275:           break;
1276:         }
1277:       }
1278:     }
1279:     ISRestoreIndices(is,&idx);
1280:   }
1281:   return(0);
1282: }

1284: /*@C
1285:    ISRestoreIndices - Restores an index set to a usable state after a call
1286:                       to ISGetIndices().

1288:    Not Collective

1290:    Input Parameters:
1291: +  is - the index set
1292: -  ptr - the pointer obtained by ISGetIndices()

1294:    Fortran Note:
1295:    This routine is used differently from Fortran
1296: $    IS          is
1297: $    integer     is_array(1)
1298: $    PetscOffset i_is
1299: $    int         ierr
1300: $       call ISGetIndices(is,is_array,i_is,ierr)
1301: $
1302: $   Access first local entry in list
1303: $      value = is_array(i_is + 1)
1304: $
1305: $      ...... other code
1306: $       call ISRestoreIndices(is,is_array,i_is,ierr)

1308:    See the Fortran chapter of the users manual and
1309:    petsc/src/is/examples/[tutorials,tests] for details.

1311:    Level: intermediate

1313:    Note:
1314:    This routine zeros out ptr. This is to prevent accidental us of the array after it has been restored.

1316: .seealso: ISGetIndices(), ISRestoreIndicesF90()
1317: @*/
1318: PetscErrorCode  ISRestoreIndices(IS is,const PetscInt *ptr[])
1319: {

1325:   if (is->ops->restoreindices) {
1326:     (*is->ops->restoreindices)(is,ptr);
1327:   }
1328:   return(0);
1329: }

1331: static PetscErrorCode ISGatherTotal_Private(IS is)
1332: {
1334:   PetscInt       i,n,N;
1335:   const PetscInt *lindices;
1336:   MPI_Comm       comm;
1337:   PetscMPIInt    rank,size,*sizes = NULL,*offsets = NULL,nn;


1342:   PetscObjectGetComm((PetscObject)is,&comm);
1343:   MPI_Comm_size(comm,&size);
1344:   MPI_Comm_rank(comm,&rank);
1345:   ISGetLocalSize(is,&n);
1346:   PetscMalloc2(size,&sizes,size,&offsets);

1348:   PetscMPIIntCast(n,&nn);
1349:   MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
1350:   offsets[0] = 0;
1351:   for (i=1; i<size; ++i) offsets[i] = offsets[i-1] + sizes[i-1];
1352:   N = offsets[size-1] + sizes[size-1];

1354:   PetscMalloc1(N,&(is->total));
1355:   ISGetIndices(is,&lindices);
1356:   MPI_Allgatherv((void*)lindices,nn,MPIU_INT,is->total,sizes,offsets,MPIU_INT,comm);
1357:   ISRestoreIndices(is,&lindices);
1358:   is->local_offset = offsets[rank];
1359:   PetscFree2(sizes,offsets);
1360:   return(0);
1361: }

1363: /*@C
1364:    ISGetTotalIndices - Retrieve an array containing all indices across the communicator.

1366:    Collective on IS

1368:    Input Parameter:
1369: .  is - the index set

1371:    Output Parameter:
1372: .  indices - total indices with rank 0 indices first, and so on; total array size is
1373:              the same as returned with ISGetSize().

1375:    Level: intermediate

1377:    Notes:
1378:     this is potentially nonscalable, but depends on the size of the total index set
1379:      and the size of the communicator. This may be feasible for index sets defined on
1380:      subcommunicators, such that the set size does not grow with PETSC_WORLD_COMM.
1381:      Note also that there is no way to tell where the local part of the indices starts
1382:      (use ISGetIndices() and ISGetNonlocalIndices() to retrieve just the local and just
1383:       the nonlocal part (complement), respectively).

1385: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices(), ISGetSize()
1386: @*/
1387: PetscErrorCode ISGetTotalIndices(IS is, const PetscInt *indices[])
1388: {
1390:   PetscMPIInt    size;

1395:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
1396:   if (size == 1) {
1397:     (*is->ops->getindices)(is,indices);
1398:   } else {
1399:     if (!is->total) {
1400:       ISGatherTotal_Private(is);
1401:     }
1402:     *indices = is->total;
1403:   }
1404:   return(0);
1405: }

1407: /*@C
1408:    ISRestoreTotalIndices - Restore the index array obtained with ISGetTotalIndices().

1410:    Not Collective.

1412:    Input Parameter:
1413: +  is - the index set
1414: -  indices - index array; must be the array obtained with ISGetTotalIndices()

1416:    Level: intermediate

1418: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices()
1419: @*/
1420: PetscErrorCode  ISRestoreTotalIndices(IS is, const PetscInt *indices[])
1421: {
1423:   PetscMPIInt    size;

1428:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
1429:   if (size == 1) {
1430:     (*is->ops->restoreindices)(is,indices);
1431:   } else {
1432:     if (is->total != *indices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
1433:   }
1434:   return(0);
1435: }
1436: /*@C
1437:    ISGetNonlocalIndices - Retrieve an array of indices from remote processors
1438:                        in this communicator.

1440:    Collective on IS

1442:    Input Parameter:
1443: .  is - the index set

1445:    Output Parameter:
1446: .  indices - indices with rank 0 indices first, and so on,  omitting
1447:              the current rank.  Total number of indices is the difference
1448:              total and local, obtained with ISGetSize() and ISGetLocalSize(),
1449:              respectively.

1451:    Level: intermediate

1453:    Notes:
1454:     restore the indices using ISRestoreNonlocalIndices().
1455:           The same scalability considerations as those for ISGetTotalIndices
1456:           apply here.

1458: .seealso: ISGetTotalIndices(), ISRestoreNonlocalIndices(), ISGetSize(), ISGetLocalSize().
1459: @*/
1460: PetscErrorCode  ISGetNonlocalIndices(IS is, const PetscInt *indices[])
1461: {
1463:   PetscMPIInt    size;
1464:   PetscInt       n, N;

1469:   MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
1470:   if (size == 1) *indices = NULL;
1471:   else {
1472:     if (!is->total) {
1473:       ISGatherTotal_Private(is);
1474:     }
1475:     ISGetLocalSize(is,&n);
1476:     ISGetSize(is,&N);
1477:     PetscMalloc1(N-n, &(is->nonlocal));
1478:     PetscArraycpy(is->nonlocal, is->total, is->local_offset);
1479:     PetscArraycpy(is->nonlocal+is->local_offset, is->total+is->local_offset+n,N - is->local_offset - n);
1480:     *indices = is->nonlocal;
1481:   }
1482:   return(0);
1483: }

1485: /*@C
1486:    ISRestoreTotalIndices - Restore the index array obtained with ISGetNonlocalIndices().

1488:    Not Collective.

1490:    Input Parameter:
1491: +  is - the index set
1492: -  indices - index array; must be the array obtained with ISGetNonlocalIndices()

1494:    Level: intermediate

1496: .seealso: ISGetTotalIndices(), ISGetNonlocalIndices(), ISRestoreTotalIndices()
1497: @*/
1498: PetscErrorCode  ISRestoreNonlocalIndices(IS is, const PetscInt *indices[])
1499: {
1503:   if (is->nonlocal != *indices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
1504:   return(0);
1505: }

1507: /*@
1508:    ISGetNonlocalIS - Gather all nonlocal indices for this IS and present
1509:                      them as another sequential index set.


1512:    Collective on IS

1514:    Input Parameter:
1515: .  is - the index set

1517:    Output Parameter:
1518: .  complement - sequential IS with indices identical to the result of
1519:                 ISGetNonlocalIndices()

1521:    Level: intermediate

1523:    Notes:
1524:     complement represents the result of ISGetNonlocalIndices as an IS.
1525:           Therefore scalability issues similar to ISGetNonlocalIndices apply.
1526:           The resulting IS must be restored using ISRestoreNonlocalIS().

1528: .seealso: ISGetNonlocalIndices(), ISRestoreNonlocalIndices(),  ISAllGather(), ISGetSize()
1529: @*/
1530: PetscErrorCode  ISGetNonlocalIS(IS is, IS *complement)
1531: {

1537:   /* Check if the complement exists already. */
1538:   if (is->complement) {
1539:     *complement = is->complement;
1540:     PetscObjectReference((PetscObject)(is->complement));
1541:   } else {
1542:     PetscInt       N, n;
1543:     const PetscInt *idx;
1544:     ISGetSize(is, &N);
1545:     ISGetLocalSize(is,&n);
1546:     ISGetNonlocalIndices(is, &idx);
1547:     ISCreateGeneral(PETSC_COMM_SELF, N-n,idx, PETSC_USE_POINTER, &(is->complement));
1548:     PetscObjectReference((PetscObject)is->complement);
1549:     *complement = is->complement;
1550:   }
1551:   return(0);
1552: }


1555: /*@
1556:    ISRestoreNonlocalIS - Restore the IS obtained with ISGetNonlocalIS().

1558:    Not collective.

1560:    Input Parameter:
1561: +  is         - the index set
1562: -  complement - index set of is's nonlocal indices

1564:    Level: intermediate


1567: .seealso: ISGetNonlocalIS(), ISGetNonlocalIndices(), ISRestoreNonlocalIndices()
1568: @*/
1569: PetscErrorCode  ISRestoreNonlocalIS(IS is, IS *complement)
1570: {
1572:   PetscInt       refcnt;

1577:   if (*complement != is->complement) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Complement IS being restored was not obtained with ISGetNonlocalIS()");
1578:   PetscObjectGetReference((PetscObject)(is->complement), &refcnt);
1579:   if (refcnt <= 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate call to ISRestoreNonlocalIS() detected");
1580:   PetscObjectDereference((PetscObject)(is->complement));
1581:   return(0);
1582: }

1584: /*@C
1585:    ISViewFromOptions - View from Options

1587:    Collective on IS

1589:    Input Parameters:
1590: +  A - the index set
1591: .  obj - Optional object
1592: -  name - command line option

1594:    Level: intermediate
1595: .seealso:  IS, ISView, PetscObjectViewFromOptions(), ISCreate()
1596: @*/
1597: PetscErrorCode  ISViewFromOptions(IS A,PetscObject obj,const char name[])
1598: {

1603:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
1604:   return(0);
1605: }

1607: /*@C
1608:    ISView - Displays an index set.

1610:    Collective on IS

1612:    Input Parameters:
1613: +  is - the index set
1614: -  viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.

1616:    Level: intermediate

1618: .seealso: PetscViewerASCIIOpen()
1619: @*/
1620: PetscErrorCode  ISView(IS is,PetscViewer viewer)
1621: {

1626:   if (!viewer) {
1627:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)is),&viewer);
1628:   }

1632:   PetscObjectPrintClassNamePrefixType((PetscObject)is,viewer);
1633:   (*is->ops->view)(is,viewer);
1634:   return(0);
1635: }

1637: /*@
1638:   ISLoad - Loads a vector that has been stored in binary or HDF5 format with ISView().

1640:   Collective on PetscViewer

1642:   Input Parameters:
1643: + is - the newly loaded vector, this needs to have been created with ISCreate() or some related function before a call to ISLoad().
1644: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or HDF5 file viewer, obtained from PetscViewerHDF5Open()

1646:   Level: intermediate

1648:   Notes:
1649:   IF using HDF5, you must assign the IS the same name as was used in the IS
1650:   that was stored in the file using PetscObjectSetName(). Otherwise you will
1651:   get the error message: "Cannot H5DOpen2() with Vec name NAMEOFOBJECT"

1653: .seealso: PetscViewerBinaryOpen(), ISView(), MatLoad(), VecLoad()
1654: @*/
1655: PetscErrorCode ISLoad(IS is, PetscViewer viewer)
1656: {
1657:   PetscBool      isbinary, ishdf5;

1663:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
1664:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1665:   if (!isbinary && !ishdf5) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1666:   if (!((PetscObject)is)->type_name) {ISSetType(is, ISGENERAL);}
1667:   PetscLogEventBegin(IS_Load,viewer,0,0,0);
1668:   (*is->ops->load)(is, viewer);
1669:   PetscLogEventEnd(IS_Load,viewer,0,0,0);
1670:   return(0);
1671: }

1673: /*@
1674:    ISSort - Sorts the indices of an index set.

1676:    Collective on IS

1678:    Input Parameters:
1679: .  is - the index set

1681:    Level: intermediate


1684: .seealso: ISSortRemoveDups(), ISSorted()
1685: @*/
1686: PetscErrorCode  ISSort(IS is)
1687: {

1692:   (*is->ops->sort)(is);
1693:   ISSetInfo(is,IS_SORTED,IS_LOCAL,is->info_permanent[IS_LOCAL][IS_SORTED],PETSC_TRUE);
1694:   return(0);
1695: }

1697: /*@
1698:   ISSortRemoveDups - Sorts the indices of an index set, removing duplicates.

1700:   Collective on IS

1702:   Input Parameters:
1703: . is - the index set

1705:   Level: intermediate


1708: .seealso: ISSort(), ISSorted()
1709: @*/
1710: PetscErrorCode ISSortRemoveDups(IS is)
1711: {

1716:   ISClearInfoCache(is,PETSC_FALSE);
1717:   (*is->ops->sortremovedups)(is);
1718:   ISSetInfo(is,IS_SORTED,IS_LOCAL,is->info_permanent[IS_LOCAL][IS_SORTED],PETSC_TRUE);
1719:   ISSetInfo(is,IS_UNIQUE,IS_LOCAL,is->info_permanent[IS_LOCAL][IS_UNIQUE],PETSC_TRUE);
1720:   return(0);
1721: }

1723: /*@
1724:    ISToGeneral - Converts an IS object of any type to ISGENERAL type

1726:    Collective on IS

1728:    Input Parameters:
1729: .  is - the index set

1731:    Level: intermediate


1734: .seealso: ISSorted()
1735: @*/
1736: PetscErrorCode  ISToGeneral(IS is)
1737: {

1742:   if (is->ops->togeneral) {
1743:     (*is->ops->togeneral)(is);
1744:   } else SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_SUP,"Not written for this type %s",((PetscObject)is)->type_name);
1745:   return(0);
1746: }

1748: /*@
1749:    ISSorted - Checks the indices to determine whether they have been sorted.

1751:    Collective on IS

1753:    Input Parameter:
1754: .  is - the index set

1756:    Output Parameter:
1757: .  flg - output flag, either PETSC_TRUE if the index set is sorted,
1758:          or PETSC_FALSE otherwise.

1760:    Notes:
1761:     For parallel IS objects this only indicates if the local part of the IS
1762:           is sorted. So some processors may return PETSC_TRUE while others may
1763:           return PETSC_FALSE.

1765:    Level: intermediate

1767: .seealso: ISSort(), ISSortRemoveDups()
1768: @*/
1769: PetscErrorCode  ISSorted(IS is,PetscBool  *flg)
1770: {

1776:   ISGetInfo(is,IS_SORTED,IS_LOCAL,PETSC_TRUE,flg);
1777:   return(0);
1778: }

1780: /*@
1781:    ISDuplicate - Creates a duplicate copy of an index set.

1783:    Collective on IS

1785:    Input Parmeters:
1786: .  is - the index set

1788:    Output Parameters:
1789: .  isnew - the copy of the index set

1791:    Level: beginner

1793: .seealso: ISCreateGeneral(), ISCopy()
1794: @*/
1795: PetscErrorCode  ISDuplicate(IS is,IS *newIS)
1796: {

1802:   (*is->ops->duplicate)(is,newIS);
1803:   ISCopyInfo(is,*newIS);
1804:   return(0);
1805: }

1807: /*@
1808:    ISCopy - Copies an index set.

1810:    Collective on IS

1812:    Input Parmeters:
1813: .  is - the index set

1815:    Output Parameters:
1816: .  isy - the copy of the index set

1818:    Level: beginner

1820: .seealso: ISDuplicate()
1821: @*/
1822: PetscErrorCode  ISCopy(IS is,IS isy)
1823: {

1830:   if (is == isy) return(0);
1831:   ISCopyInfo(is,isy);
1832:   isy->max        = is->max;
1833:   isy->min        = is->min;
1834:   (*is->ops->copy)(is,isy);
1835:   return(0);
1836: }

1838: /*@
1839:    ISOnComm - Split a parallel IS on subcomms (usually self) or concatenate index sets on subcomms into a parallel index set

1841:    Collective on IS

1843:    Input Arguments:
1844: + is - index set
1845: . comm - communicator for new index set
1846: - mode - copy semantics, PETSC_USE_POINTER for no-copy if possible, otherwise PETSC_COPY_VALUES

1848:    Output Arguments:
1849: . newis - new IS on comm

1851:    Level: advanced

1853:    Notes:
1854:    It is usually desirable to create a parallel IS and look at the local part when necessary.

1856:    This function is useful if serial ISs must be created independently, or to view many
1857:    logically independent serial ISs.

1859:    The input IS must have the same type on every process.

1861: .seealso: ISSplit()
1862: @*/
1863: PetscErrorCode  ISOnComm(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
1864: {
1866:   PetscMPIInt    match;

1871:   MPI_Comm_compare(PetscObjectComm((PetscObject)is),comm,&match);
1872:   if (mode != PETSC_COPY_VALUES && (match == MPI_IDENT || match == MPI_CONGRUENT)) {
1873:     PetscObjectReference((PetscObject)is);
1874:     *newis = is;
1875:   } else {
1876:     (*is->ops->oncomm)(is,comm,mode,newis);
1877:   }
1878:   return(0);
1879: }

1881: /*@
1882:    ISSetBlockSize - informs an index set that it has a given block size

1884:    Logicall Collective on IS

1886:    Input Arguments:
1887: + is - index set
1888: - bs - block size

1890:    Level: intermediate
1891:    
1892:    Notes: 
1893:    This is much like the block size for Vecs. It indicates that one can think of the indices as 
1894:    being in a collection of equal size blocks. For ISBlock() these collections of blocks are all contiquous
1895:    within a block but this is not the case for other IS.
1896:    ISBlockGetIndices() only works for ISBlock IS, not others.

1898: .seealso: ISGetBlockSize(), ISCreateBlock(), ISBlockGetIndices(), 
1899: @*/
1900: PetscErrorCode  ISSetBlockSize(IS is,PetscInt bs)
1901: {

1907:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_OUTOFRANGE,"Block size %D, must be positive",bs);
1908:   (*is->ops->setblocksize)(is,bs);
1909:   return(0);
1910: }

1912: /*@
1913:    ISGetBlockSize - Returns the number of elements in a block.

1915:    Not Collective

1917:    Input Parameter:
1918: .  is - the index set

1920:    Output Parameter:
1921: .  size - the number of elements in a block

1923:    Level: intermediate

1925: Notes: 
1926:    This is much like the block size for Vecs. It indicates that one can think of the indices as 
1927:    being in a collection of equal size blocks. For ISBlock() these collections of blocks are all contiquous
1928:    within a block but this is not the case for other IS.
1929:    ISBlockGetIndices() only works for ISBlock IS, not others.

1931: .seealso: ISBlockGetSize(), ISGetSize(), ISCreateBlock(), ISSetBlockSize()
1932: @*/
1933: PetscErrorCode  ISGetBlockSize(IS is,PetscInt *size)
1934: {

1938:   PetscLayoutGetBlockSize(is->map, size);
1939:   return(0);
1940: }

1942: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[])
1943: {
1945:   PetscInt       len,i;
1946:   const PetscInt *ptr;

1949:   ISGetLocalSize(is,&len);
1950:   ISGetIndices(is,&ptr);
1951:   for (i=0; i<len; i++) idx[i] = ptr[i];
1952:   ISRestoreIndices(is,&ptr);
1953:   return(0);
1954: }

1956: /*MC
1957:     ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1958:     The users should call ISRestoreIndicesF90() after having looked at the
1959:     indices.  The user should NOT change the indices.

1961:     Synopsis:
1962:     ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1964:     Not collective

1966:     Input Parameter:
1967: .   x - index set

1969:     Output Parameters:
1970: +   xx_v - the Fortran90 pointer to the array
1971: -   ierr - error code

1973:     Example of Usage:
1974: .vb
1975:     PetscInt, pointer xx_v(:)
1976:     ....
1977:     call ISGetIndicesF90(x,xx_v,ierr)
1978:     a = xx_v(3)
1979:     call ISRestoreIndicesF90(x,xx_v,ierr)
1980: .ve

1982:     Level: intermediate

1984: .seealso:  ISRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices()


1987: M*/

1989: /*MC
1990:     ISRestoreIndicesF90 - Restores an index set to a usable state after
1991:     a call to ISGetIndicesF90().

1993:     Synopsis:
1994:     ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

1996:     Not collective

1998:     Input Parameters:
1999: +   x - index set
2000: -   xx_v - the Fortran90 pointer to the array

2002:     Output Parameter:
2003: .   ierr - error code


2006:     Example of Usage:
2007: .vb
2008:     PetscInt, pointer xx_v(:)
2009:     ....
2010:     call ISGetIndicesF90(x,xx_v,ierr)
2011:     a = xx_v(3)
2012:     call ISRestoreIndicesF90(x,xx_v,ierr)
2013: .ve

2015:     Level: intermediate

2017: .seealso:  ISGetIndicesF90(), ISGetIndices(), ISRestoreIndices()

2019: M*/

2021: /*MC
2022:     ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
2023:     The users should call ISBlockRestoreIndicesF90() after having looked at the
2024:     indices.  The user should NOT change the indices.

2026:     Synopsis:
2027:     ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

2029:     Not collective

2031:     Input Parameter:
2032: .   x - index set

2034:     Output Parameters:
2035: +   xx_v - the Fortran90 pointer to the array
2036: -   ierr - error code
2037:     Example of Usage:
2038: .vb
2039:     PetscInt, pointer xx_v(:)
2040:     ....
2041:     call ISBlockGetIndicesF90(x,xx_v,ierr)
2042:     a = xx_v(3)
2043:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
2044: .ve

2046:     Level: intermediate

2048: .seealso:  ISBlockRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices(),
2049:            ISRestoreIndices()

2051: M*/

2053: /*MC
2054:     ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
2055:     a call to ISBlockGetIndicesF90().

2057:     Synopsis:
2058:     ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)

2060:     Not Collective

2062:     Input Parameters:
2063: +   x - index set
2064: -   xx_v - the Fortran90 pointer to the array

2066:     Output Parameter:
2067: .   ierr - error code

2069:     Example of Usage:
2070: .vb
2071:     PetscInt, pointer xx_v(:)
2072:     ....
2073:     call ISBlockGetIndicesF90(x,xx_v,ierr)
2074:     a = xx_v(3)
2075:     call ISBlockRestoreIndicesF90(x,xx_v,ierr)
2076: .ve

2078:     Notes:
2079:     Not yet supported for all F90 compilers

2081:     Level: intermediate

2083: .seealso:  ISBlockGetIndicesF90(), ISGetIndices(), ISRestoreIndices(), ISRestoreIndicesF90()

2085: M*/