Actual source code: index.c
petsc-3.8.4 2018-03-24
2: /*
3: Defines the abstract operations on index sets, i.e. the public interface.
4: */
5: #include <petsc/private/isimpl.h>
6: #include <petscviewer.h>
7: #include <petscsf.h>
9: /* Logging support */
10: PetscClassId IS_CLASSID;
12: /*@
13: ISRenumber - Renumbers an index set (with multiplicities) in a contiguous way.
15: Collective on IS
17: Input Parmeters:
18: . subset - the index set
19: . subset_mult - the multiplcity of each entry in subset (optional, can be NULL)
21: Output Parameters:
22: . N - the maximum entry of the new IS
23: . subset_n - the new IS
25: Level: intermediate
27: .seealso:
28: @*/
29: PetscErrorCode ISRenumber(IS subset, IS subset_mult, PetscInt *N, IS *subset_n)
30: {
31: PetscSF sf;
32: PetscLayout map;
33: const PetscInt *idxs;
34: PetscInt *leaf_data,*root_data,*gidxs;
35: PetscInt N_n,n,i,lbounds[2],gbounds[2],Nl;
36: PetscInt n_n,nlocals,start,first_index;
37: PetscMPIInt commsize;
38: PetscBool first_found;
42: ISGetLocalSize(subset,&n);
43: if (subset_mult) {
45: ISGetLocalSize(subset,&i);
46: if (i != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Local subset and multiplicity sizes don't match! %d != %d",n,i);
47: }
48: /* create workspace layout for computing global indices of subset */
49: ISGetIndices(subset,&idxs);
50: lbounds[0] = lbounds[1] = 0;
51: for (i=0;i<n;i++) {
52: if (idxs[i] < lbounds[0]) lbounds[0] = idxs[i];
53: else if (idxs[i] > lbounds[1]) lbounds[1] = idxs[i];
54: }
55: lbounds[0] = -lbounds[0];
56: MPIU_Allreduce(lbounds,gbounds,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)subset));
57: gbounds[0] = -gbounds[0];
58: N_n= gbounds[1] - gbounds[0] + 1;
59: PetscLayoutCreate(PetscObjectComm((PetscObject)subset),&map);
60: PetscLayoutSetBlockSize(map,1);
61: PetscLayoutSetSize(map,N_n);
62: PetscLayoutSetUp(map);
63: PetscLayoutGetLocalSize(map,&Nl);
65: /* create sf : leaf_data == multiplicity of indexes, root data == global index in layout */
66: PetscMalloc2(n,&leaf_data,Nl,&root_data);
67: if (subset_mult) {
68: const PetscInt* idxs_mult;
70: ISGetIndices(subset_mult,&idxs_mult);
71: PetscMemcpy(leaf_data,idxs_mult,n*sizeof(PetscInt));
72: ISRestoreIndices(subset_mult,&idxs_mult);
73: } else {
74: for (i=0;i<n;i++) leaf_data[i] = 1;
75: }
76: /* local size of new subset */
77: n_n = 0;
78: for (i=0;i<n;i++) n_n += leaf_data[i];
80: /* global indexes in layout */
81: PetscMalloc1(n_n,&gidxs); /* allocating possibly extra space in gidxs which will be used later */
82: for (i=0;i<n;i++) gidxs[i] = idxs[i] - gbounds[0];
83: ISRestoreIndices(subset,&idxs);
84: PetscSFCreate(PetscObjectComm((PetscObject)subset),&sf);
85: PetscSFSetGraphLayout(sf,map,n,NULL,PETSC_COPY_VALUES,gidxs);
86: PetscLayoutDestroy(&map);
88: /* reduce from leaves to roots */
89: PetscMemzero(root_data,Nl*sizeof(PetscInt));
90: PetscSFReduceBegin(sf,MPIU_INT,leaf_data,root_data,MPI_MAX);
91: PetscSFReduceEnd(sf,MPIU_INT,leaf_data,root_data,MPI_MAX);
93: /* count indexes in local part of layout */
94: nlocals = 0;
95: first_index = -1;
96: first_found = PETSC_FALSE;
97: for (i=0;i<Nl;i++) {
98: if (!first_found && root_data[i]) {
99: first_found = PETSC_TRUE;
100: first_index = i;
101: }
102: nlocals += root_data[i];
103: }
105: /* cumulative of number of indexes and size of subset without holes */
106: #if defined(PETSC_HAVE_MPI_EXSCAN)
107: start = 0;
108: MPI_Exscan(&nlocals,&start,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)subset));
109: #else
110: MPI_Scan(&nlocals,&start,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)subset));
111: start = start-nlocals;
112: #endif
114: if (N) { /* compute total size of new subset if requested */
115: *N = start + nlocals;
116: MPI_Comm_size(PetscObjectComm((PetscObject)subset),&commsize);
117: MPI_Bcast(N,1,MPIU_INT,commsize-1,PetscObjectComm((PetscObject)subset));
118: }
120: /* adapt root data with cumulative */
121: if (first_found) {
122: PetscInt old_index;
124: root_data[first_index] += start;
125: old_index = first_index;
126: for (i=first_index+1;i<Nl;i++) {
127: if (root_data[i]) {
128: root_data[i] += root_data[old_index];
129: old_index = i;
130: }
131: }
132: }
134: /* from roots to leaves */
135: PetscSFBcastBegin(sf,MPIU_INT,root_data,leaf_data);
136: PetscSFBcastEnd(sf,MPIU_INT,root_data,leaf_data);
137: PetscSFDestroy(&sf);
139: /* create new IS with global indexes without holes */
140: if (subset_mult) {
141: const PetscInt* idxs_mult;
142: PetscInt cum;
144: cum = 0;
145: ISGetIndices(subset_mult,&idxs_mult);
146: for (i=0;i<n;i++) {
147: PetscInt j;
148: for (j=0;j<idxs_mult[i];j++) gidxs[cum++] = leaf_data[i] - idxs_mult[i] + j;
149: }
150: ISRestoreIndices(subset_mult,&idxs_mult);
151: } else {
152: for (i=0;i<n;i++) {
153: gidxs[i] = leaf_data[i]-1;
154: }
155: }
156: ISCreateGeneral(PetscObjectComm((PetscObject)subset),n_n,gidxs,PETSC_OWN_POINTER,subset_n);
157: PetscFree2(leaf_data,root_data);
158: return(0);
159: }
161: /*@
162: ISIdentity - Determines whether index set is the identity mapping.
164: Collective on IS
166: Input Parmeters:
167: . is - the index set
169: Output Parameters:
170: . ident - PETSC_TRUE if an identity, else PETSC_FALSE
172: Level: intermediate
174: Concepts: identity mapping
175: Concepts: index sets^is identity
177: .seealso: ISSetIdentity()
178: @*/
179: PetscErrorCode ISIdentity(IS is,PetscBool *ident)
180: {
186: *ident = is->isidentity;
187: if (*ident) return(0);
188: if (is->ops->identity) {
189: (*is->ops->identity)(is,ident);
190: }
191: return(0);
192: }
194: /*@
195: ISSetIdentity - Informs the index set that it is an identity.
197: Logically Collective on IS
199: Input Parmeters:
200: . is - the index set
202: Level: intermediate
204: Concepts: identity mapping
205: Concepts: index sets^is identity
207: .seealso: ISIdentity()
208: @*/
209: PetscErrorCode ISSetIdentity(IS is)
210: {
215: is->isidentity = PETSC_TRUE;
216: ISSetPermutation(is);
217: return(0);
218: }
220: /*@
221: ISContiguousLocal - Locates an index set with contiguous range within a global range, if possible
223: Not Collective
225: Input Parmeters:
226: + is - the index set
227: . gstart - global start
228: . gend - global end
230: Output Parameters:
231: + start - start of contiguous block, as an offset from gstart
232: - contig - PETSC_TRUE if the index set refers to contiguous entries on this process, else PETSC_FALSE
234: Level: developer
236: Concepts: index sets^is contiguous
238: .seealso: ISGetLocalSize(), VecGetOwnershipRange()
239: @*/
240: PetscErrorCode ISContiguousLocal(IS is,PetscInt gstart,PetscInt gend,PetscInt *start,PetscBool *contig)
241: {
248: if (is->ops->contiguous) {
249: (*is->ops->contiguous)(is,gstart,gend,start,contig);
250: } else {
251: *start = -1;
252: *contig = PETSC_FALSE;
253: }
254: return(0);
255: }
257: /*@
258: ISPermutation - PETSC_TRUE or PETSC_FALSE depending on whether the
259: index set has been declared to be a permutation.
261: Logically Collective on IS
263: Input Parmeters:
264: . is - the index set
266: Output Parameters:
267: . perm - PETSC_TRUE if a permutation, else PETSC_FALSE
269: Level: intermediate
271: Concepts: permutation
272: Concepts: index sets^is permutation
274: .seealso: ISSetPermutation()
275: @*/
276: PetscErrorCode ISPermutation(IS is,PetscBool *perm)
277: {
281: *perm = (PetscBool) is->isperm;
282: return(0);
283: }
285: /*@
286: ISSetPermutation - Informs the index set that it is a permutation.
288: Logically Collective on IS
290: Input Parmeters:
291: . is - the index set
293: Level: intermediate
295: Concepts: permutation
296: Concepts: index sets^permutation
298: The debug version of the libraries (./configure --with-debugging=1) checks if the
299: index set is actually a permutation. The optimized version just believes you.
301: .seealso: ISPermutation()
302: @*/
303: PetscErrorCode ISSetPermutation(IS is)
304: {
307: #if defined(PETSC_USE_DEBUG)
308: {
309: PetscMPIInt size;
312: MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
313: if (size == 1) {
314: PetscInt i,n,*idx;
315: const PetscInt *iidx;
317: ISGetSize(is,&n);
318: PetscMalloc1(n,&idx);
319: ISGetIndices(is,&iidx);
320: PetscMemcpy(idx,iidx,n*sizeof(PetscInt));
321: PetscSortInt(n,idx);
322: for (i=0; i<n; i++) {
323: if (idx[i] != i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Index set is not a permutation");
324: }
325: PetscFree(idx);
326: ISRestoreIndices(is,&iidx);
327: }
328: }
329: #endif
330: is->isperm = PETSC_TRUE;
331: return(0);
332: }
334: /*@
335: ISDestroy - Destroys an index set.
337: Collective on IS
339: Input Parameters:
340: . is - the index set
342: Level: beginner
344: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlocked()
345: @*/
346: PetscErrorCode ISDestroy(IS *is)
347: {
351: if (!*is) return(0);
353: if (--((PetscObject)(*is))->refct > 0) {*is = 0; return(0);}
354: if ((*is)->complement) {
355: PetscInt refcnt;
356: PetscObjectGetReference((PetscObject)((*is)->complement), &refcnt);
357: if (refcnt > 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlocal IS has not been restored");
358: ISDestroy(&(*is)->complement);
359: }
360: if ((*is)->ops->destroy) {
361: (*(*is)->ops->destroy)(*is);
362: }
363: PetscLayoutDestroy(&(*is)->map);
364: /* Destroy local representations of offproc data. */
365: PetscFree((*is)->total);
366: PetscFree((*is)->nonlocal);
367: PetscHeaderDestroy(is);
368: return(0);
369: }
371: /*@
372: ISInvertPermutation - Creates a new permutation that is the inverse of
373: a given permutation.
375: Collective on IS
377: Input Parameter:
378: + is - the index set
379: - nlocal - number of indices on this processor in result (ignored for 1 proccessor) or
380: use PETSC_DECIDE
382: Output Parameter:
383: . isout - the inverse permutation
385: Level: intermediate
387: Notes: For parallel index sets this does the complete parallel permutation, but the
388: code is not efficient for huge index sets (10,000,000 indices).
390: Concepts: inverse permutation
391: Concepts: permutation^inverse
392: Concepts: index sets^inverting
393: @*/
394: PetscErrorCode ISInvertPermutation(IS is,PetscInt nlocal,IS *isout)
395: {
401: if (!is->isperm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a permutation, must call ISSetPermutation() on the IS first");
402: if (is->isidentity) {
403: ISDuplicate(is,isout);
404: } else {
405: (*is->ops->invertpermutation)(is,nlocal,isout);
406: ISSetPermutation(*isout);
407: }
408: return(0);
409: }
411: /*@
412: ISGetSize - Returns the global length of an index set.
414: Not Collective
416: Input Parameter:
417: . is - the index set
419: Output Parameter:
420: . size - the global size
422: Level: beginner
424: Concepts: size^of index set
425: Concepts: index sets^size
427: @*/
428: PetscErrorCode ISGetSize(IS is,PetscInt *size)
429: {
435: (*is->ops->getsize)(is,size);
436: return(0);
437: }
439: /*@
440: ISGetLocalSize - Returns the local (processor) length of an index set.
442: Not Collective
444: Input Parameter:
445: . is - the index set
447: Output Parameter:
448: . size - the local size
450: Level: beginner
452: Concepts: size^of index set
453: Concepts: local size^of index set
454: Concepts: index sets^local size
456: @*/
457: PetscErrorCode ISGetLocalSize(IS is,PetscInt *size)
458: {
464: (*is->ops->getlocalsize)(is,size);
465: return(0);
466: }
468: /*@C
469: ISGetIndices - Returns a pointer to the indices. The user should call
470: ISRestoreIndices() after having looked at the indices. The user should
471: NOT change the indices.
473: Not Collective
475: Input Parameter:
476: . is - the index set
478: Output Parameter:
479: . ptr - the location to put the pointer to the indices
481: Fortran Note:
482: This routine has two different interfaces from Fortran; the first is not recommend, it does not require Fortran 90
483: $ IS is
484: $ integer is_array(1)
485: $ PetscOffset i_is
486: $ int ierr
487: $ call ISGetIndices(is,is_array,i_is,ierr)
488: $
489: $ Access first local entry in list
490: $ value = is_array(i_is + 1)
491: $
492: $ ...... other code
493: $ call ISRestoreIndices(is,is_array,i_is,ierr)
494: The second Fortran interface is recommended.
495: $ use petscisdef
496: $ PetscInt, pointer :: array(:)
497: $ PetscErrorCode ierr
498: $ IS i
499: $ call ISGetIndicesF90(i,array,ierr)
503: See the Fortran chapter of the users manual and
504: petsc/src/is/examples/[tutorials,tests] for details.
506: Level: intermediate
508: Concepts: index sets^getting indices
509: Concepts: indices of index set
511: .seealso: ISRestoreIndices(), ISGetIndicesF90()
512: @*/
513: PetscErrorCode ISGetIndices(IS is,const PetscInt *ptr[])
514: {
520: (*is->ops->getindices)(is,ptr);
521: return(0);
522: }
524: /*@C
525: ISGetMinMax - Gets the minimum and maximum values in an IS
527: Not Collective
529: Input Parameter:
530: . is - the index set
532: Output Parameter:
533: + min - the minimum value
534: - max - the maximum value
536: Level: intermediate
538: Notes: Empty index sets return min=PETSC_MAX_INT and max=PETSC_MIN_INT.
540: Concepts: index sets^getting indices
541: Concepts: indices of index set
543: .seealso: ISGetIndices(), ISRestoreIndices(), ISGetIndicesF90()
544: @*/
545: PetscErrorCode ISGetMinMax(IS is,PetscInt *min,PetscInt *max)
546: {
549: if (min) *min = is->min;
550: if (max) *max = is->max;
551: return(0);
552: }
554: /*@
555: ISLocate - determine the location of an index within the local component of an index set
557: Not Collective
559: Input Parameter:
560: + is - the index set
561: - key - the search key
563: Output Parameter:
564: . location - if >= 0, a location within the index set that is equal to the key, otherwise the key is not in the index set
566: Level: intermediate
567: @*/
568: PetscErrorCode ISLocate(IS is, PetscInt key, PetscInt *location)
569: {
573: if (is->ops->locate) {
574: (*is->ops->locate)(is,key,location);
575: } else {
576: PetscInt numIdx;
577: PetscBool sorted;
578: const PetscInt *idx;
580: ISGetLocalSize(is,&numIdx);
581: ISGetIndices(is,&idx);
582: ISSorted(is,&sorted);
583: if (sorted) {
584: PetscFindInt(key,numIdx,idx,location);
585: } else {
586: PetscInt i;
588: *location = -1;
589: for (i = 0; i < numIdx; i++) {
590: if (idx[i] == key) {
591: *location = i;
592: break;
593: }
594: }
595: }
596: ISRestoreIndices(is,&idx);
597: }
598: return(0);
599: }
601: /*@C
602: ISRestoreIndices - Restores an index set to a usable state after a call
603: to ISGetIndices().
605: Not Collective
607: Input Parameters:
608: + is - the index set
609: - ptr - the pointer obtained by ISGetIndices()
611: Fortran Note:
612: This routine is used differently from Fortran
613: $ IS is
614: $ integer is_array(1)
615: $ PetscOffset i_is
616: $ int ierr
617: $ call ISGetIndices(is,is_array,i_is,ierr)
618: $
619: $ Access first local entry in list
620: $ value = is_array(i_is + 1)
621: $
622: $ ...... other code
623: $ call ISRestoreIndices(is,is_array,i_is,ierr)
625: See the Fortran chapter of the users manual and
626: petsc/src/is/examples/[tutorials,tests] for details.
628: Level: intermediate
630: Note:
631: This routine zeros out ptr. This is to prevent accidental us of the array after it has been restored.
633: .seealso: ISGetIndices(), ISRestoreIndicesF90()
634: @*/
635: PetscErrorCode ISRestoreIndices(IS is,const PetscInt *ptr[])
636: {
642: if (is->ops->restoreindices) {
643: (*is->ops->restoreindices)(is,ptr);
644: }
645: return(0);
646: }
648: static PetscErrorCode ISGatherTotal_Private(IS is)
649: {
651: PetscInt i,n,N;
652: const PetscInt *lindices;
653: MPI_Comm comm;
654: PetscMPIInt rank,size,*sizes = NULL,*offsets = NULL,nn;
659: PetscObjectGetComm((PetscObject)is,&comm);
660: MPI_Comm_size(comm,&size);
661: MPI_Comm_rank(comm,&rank);
662: ISGetLocalSize(is,&n);
663: PetscMalloc2(size,&sizes,size,&offsets);
665: PetscMPIIntCast(n,&nn);
666: MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
667: offsets[0] = 0;
668: for (i=1; i<size; ++i) offsets[i] = offsets[i-1] + sizes[i-1];
669: N = offsets[size-1] + sizes[size-1];
671: PetscMalloc1(N,&(is->total));
672: ISGetIndices(is,&lindices);
673: MPI_Allgatherv((void*)lindices,nn,MPIU_INT,is->total,sizes,offsets,MPIU_INT,comm);
674: ISRestoreIndices(is,&lindices);
675: is->local_offset = offsets[rank];
676: PetscFree2(sizes,offsets);
677: return(0);
678: }
680: /*@C
681: ISGetTotalIndices - Retrieve an array containing all indices across the communicator.
683: Collective on IS
685: Input Parameter:
686: . is - the index set
688: Output Parameter:
689: . indices - total indices with rank 0 indices first, and so on; total array size is
690: the same as returned with ISGetSize().
692: Level: intermediate
694: Notes: this is potentially nonscalable, but depends on the size of the total index set
695: and the size of the communicator. This may be feasible for index sets defined on
696: subcommunicators, such that the set size does not grow with PETSC_WORLD_COMM.
697: Note also that there is no way to tell where the local part of the indices starts
698: (use ISGetIndices() and ISGetNonlocalIndices() to retrieve just the local and just
699: the nonlocal part (complement), respectively).
701: Concepts: index sets^getting nonlocal indices
702: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices(), ISGetSize()
703: @*/
704: PetscErrorCode ISGetTotalIndices(IS is, const PetscInt *indices[])
705: {
707: PetscMPIInt size;
712: MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
713: if (size == 1) {
714: (*is->ops->getindices)(is,indices);
715: } else {
716: if (!is->total) {
717: ISGatherTotal_Private(is);
718: }
719: *indices = is->total;
720: }
721: return(0);
722: }
724: /*@C
725: ISRestoreTotalIndices - Restore the index array obtained with ISGetTotalIndices().
727: Not Collective.
729: Input Parameter:
730: + is - the index set
731: - indices - index array; must be the array obtained with ISGetTotalIndices()
733: Level: intermediate
735: Concepts: index sets^getting nonlocal indices
736: Concepts: index sets^restoring nonlocal indices
737: .seealso: ISRestoreTotalIndices(), ISGetNonlocalIndices()
738: @*/
739: PetscErrorCode ISRestoreTotalIndices(IS is, const PetscInt *indices[])
740: {
742: PetscMPIInt size;
747: MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
748: if (size == 1) {
749: (*is->ops->restoreindices)(is,indices);
750: } else {
751: 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.");
752: }
753: return(0);
754: }
755: /*@C
756: ISGetNonlocalIndices - Retrieve an array of indices from remote processors
757: in this communicator.
759: Collective on IS
761: Input Parameter:
762: . is - the index set
764: Output Parameter:
765: . indices - indices with rank 0 indices first, and so on, omitting
766: the current rank. Total number of indices is the difference
767: total and local, obtained with ISGetSize() and ISGetLocalSize(),
768: respectively.
770: Level: intermediate
772: Notes: restore the indices using ISRestoreNonlocalIndices().
773: The same scalability considerations as those for ISGetTotalIndices
774: apply here.
776: Concepts: index sets^getting nonlocal indices
777: .seealso: ISGetTotalIndices(), ISRestoreNonlocalIndices(), ISGetSize(), ISGetLocalSize().
778: @*/
779: PetscErrorCode ISGetNonlocalIndices(IS is, const PetscInt *indices[])
780: {
782: PetscMPIInt size;
783: PetscInt n, N;
788: MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
789: if (size == 1) *indices = NULL;
790: else {
791: if (!is->total) {
792: ISGatherTotal_Private(is);
793: }
794: ISGetLocalSize(is,&n);
795: ISGetSize(is,&N);
796: PetscMalloc1(N-n, &(is->nonlocal));
797: PetscMemcpy(is->nonlocal, is->total, sizeof(PetscInt)*is->local_offset);
798: PetscMemcpy(is->nonlocal+is->local_offset, is->total+is->local_offset+n, sizeof(PetscInt)*(N - is->local_offset - n));
799: *indices = is->nonlocal;
800: }
801: return(0);
802: }
804: /*@C
805: ISRestoreTotalIndices - Restore the index array obtained with ISGetNonlocalIndices().
807: Not Collective.
809: Input Parameter:
810: + is - the index set
811: - indices - index array; must be the array obtained with ISGetNonlocalIndices()
813: Level: intermediate
815: Concepts: index sets^getting nonlocal indices
816: Concepts: index sets^restoring nonlocal indices
817: .seealso: ISGetTotalIndices(), ISGetNonlocalIndices(), ISRestoreTotalIndices()
818: @*/
819: PetscErrorCode ISRestoreNonlocalIndices(IS is, const PetscInt *indices[])
820: {
824: 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.");
825: return(0);
826: }
828: /*@
829: ISGetNonlocalIS - Gather all nonlocal indices for this IS and present
830: them as another sequential index set.
833: Collective on IS
835: Input Parameter:
836: . is - the index set
838: Output Parameter:
839: . complement - sequential IS with indices identical to the result of
840: ISGetNonlocalIndices()
842: Level: intermediate
844: Notes: complement represents the result of ISGetNonlocalIndices as an IS.
845: Therefore scalability issues similar to ISGetNonlocalIndices apply.
846: The resulting IS must be restored using ISRestoreNonlocalIS().
848: Concepts: index sets^getting nonlocal indices
849: .seealso: ISGetNonlocalIndices(), ISRestoreNonlocalIndices(), ISAllGather(), ISGetSize()
850: @*/
851: PetscErrorCode ISGetNonlocalIS(IS is, IS *complement)
852: {
858: /* Check if the complement exists already. */
859: if (is->complement) {
860: *complement = is->complement;
861: PetscObjectReference((PetscObject)(is->complement));
862: } else {
863: PetscInt N, n;
864: const PetscInt *idx;
865: ISGetSize(is, &N);
866: ISGetLocalSize(is,&n);
867: ISGetNonlocalIndices(is, &idx);
868: ISCreateGeneral(PETSC_COMM_SELF, N-n,idx, PETSC_USE_POINTER, &(is->complement));
869: PetscObjectReference((PetscObject)is->complement);
870: *complement = is->complement;
871: }
872: return(0);
873: }
876: /*@
877: ISRestoreNonlocalIS - Restore the IS obtained with ISGetNonlocalIS().
879: Not collective.
881: Input Parameter:
882: + is - the index set
883: - complement - index set of is's nonlocal indices
885: Level: intermediate
888: Concepts: index sets^getting nonlocal indices
889: Concepts: index sets^restoring nonlocal indices
890: .seealso: ISGetNonlocalIS(), ISGetNonlocalIndices(), ISRestoreNonlocalIndices()
891: @*/
892: PetscErrorCode ISRestoreNonlocalIS(IS is, IS *complement)
893: {
895: PetscInt refcnt;
900: if (*complement != is->complement) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Complement IS being restored was not obtained with ISGetNonlocalIS()");
901: PetscObjectGetReference((PetscObject)(is->complement), &refcnt);
902: if (refcnt <= 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate call to ISRestoreNonlocalIS() detected");
903: PetscObjectDereference((PetscObject)(is->complement));
904: return(0);
905: }
907: /*@C
908: ISView - Displays an index set.
910: Collective on IS
912: Input Parameters:
913: + is - the index set
914: - viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.
916: Level: intermediate
918: .seealso: PetscViewerASCIIOpen()
919: @*/
920: PetscErrorCode ISView(IS is,PetscViewer viewer)
921: {
926: if (!viewer) {
927: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)is),&viewer);
928: }
932: PetscObjectPrintClassNamePrefixType((PetscObject)is,viewer);
933: (*is->ops->view)(is,viewer);
934: return(0);
935: }
937: /*@
938: ISLoad - Loads a vector that has been stored in binary or HDF5 format with ISView().
940: Collective on PetscViewer
942: Input Parameters:
943: + is - the newly loaded vector, this needs to have been created with ISCreate() or some related function before a call to ISLoad().
944: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or HDF5 file viewer, obtained from PetscViewerHDF5Open()
946: Level: intermediate
948: Notes:
949: IF using HDF5, you must assign the IS the same name as was used in the IS
950: that was stored in the file using PetscObjectSetName(). Otherwise you will
951: get the error message: "Cannot H5DOpen2() with Vec name NAMEOFOBJECT"
953: Concepts: index set^loading from file
955: .seealso: PetscViewerBinaryOpen(), ISView(), MatLoad(), VecLoad()
956: @*/
957: PetscErrorCode ISLoad(IS is, PetscViewer viewer)
958: {
959: PetscBool isbinary, ishdf5;
965: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
966: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
967: if (!isbinary && !ishdf5) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
968: if (!((PetscObject)is)->type_name) {ISSetType(is, ISGENERAL);}
969: (*is->ops->load)(is, viewer);
970: return(0);
971: }
973: /*@
974: ISSort - Sorts the indices of an index set.
976: Collective on IS
978: Input Parameters:
979: . is - the index set
981: Level: intermediate
983: Concepts: index sets^sorting
984: Concepts: sorting^index set
986: .seealso: ISSortRemoveDups(), ISSorted()
987: @*/
988: PetscErrorCode ISSort(IS is)
989: {
994: (*is->ops->sort)(is);
995: return(0);
996: }
998: /*@
999: ISSortRemoveDups - Sorts the indices of an index set, removing duplicates.
1001: Collective on IS
1003: Input Parameters:
1004: . is - the index set
1006: Level: intermediate
1008: Concepts: index sets^sorting
1009: Concepts: sorting^index set
1011: .seealso: ISSort(), ISSorted()
1012: @*/
1013: PetscErrorCode ISSortRemoveDups(IS is)
1014: {
1019: (*is->ops->sortremovedups)(is);
1020: return(0);
1021: }
1023: /*@
1024: ISToGeneral - Converts an IS object of any type to ISGENERAL type
1026: Collective on IS
1028: Input Parameters:
1029: . is - the index set
1031: Level: intermediate
1033: Concepts: index sets^sorting
1034: Concepts: sorting^index set
1036: .seealso: ISSorted()
1037: @*/
1038: PetscErrorCode ISToGeneral(IS is)
1039: {
1044: if (is->ops->togeneral) {
1045: (*is->ops->togeneral)(is);
1046: } else SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_SUP,"Not written for this type %s",((PetscObject)is)->type_name);
1047: return(0);
1048: }
1050: /*@
1051: ISSorted - Checks the indices to determine whether they have been sorted.
1053: Collective on IS
1055: Input Parameter:
1056: . is - the index set
1058: Output Parameter:
1059: . flg - output flag, either PETSC_TRUE if the index set is sorted,
1060: or PETSC_FALSE otherwise.
1062: Notes: For parallel IS objects this only indicates if the local part of the IS
1063: is sorted. So some processors may return PETSC_TRUE while others may
1064: return PETSC_FALSE.
1066: Level: intermediate
1068: .seealso: ISSort(), ISSortRemoveDups()
1069: @*/
1070: PetscErrorCode ISSorted(IS is,PetscBool *flg)
1071: {
1077: (*is->ops->sorted)(is,flg);
1078: return(0);
1079: }
1081: /*@
1082: ISDuplicate - Creates a duplicate copy of an index set.
1084: Collective on IS
1086: Input Parmeters:
1087: . is - the index set
1089: Output Parameters:
1090: . isnew - the copy of the index set
1092: Level: beginner
1094: Concepts: index sets^duplicating
1096: .seealso: ISCreateGeneral(), ISCopy()
1097: @*/
1098: PetscErrorCode ISDuplicate(IS is,IS *newIS)
1099: {
1105: (*is->ops->duplicate)(is,newIS);
1106: (*newIS)->isidentity = is->isidentity;
1107: (*newIS)->isperm = is->isperm;
1108: return(0);
1109: }
1111: /*@
1112: ISCopy - Copies an index set.
1114: Collective on IS
1116: Input Parmeters:
1117: . is - the index set
1119: Output Parameters:
1120: . isy - the copy of the index set
1122: Level: beginner
1124: Concepts: index sets^copying
1126: .seealso: ISDuplicate()
1127: @*/
1128: PetscErrorCode ISCopy(IS is,IS isy)
1129: {
1136: if (is == isy) return(0);
1137: (*is->ops->copy)(is,isy);
1138: isy->isperm = is->isperm;
1139: isy->max = is->max;
1140: isy->min = is->min;
1141: isy->isidentity = is->isidentity;
1142: return(0);
1143: }
1145: /*@
1146: ISOnComm - Split a parallel IS on subcomms (usually self) or concatenate index sets on subcomms into a parallel index set
1148: Collective on IS and comm
1150: Input Arguments:
1151: + is - index set
1152: . comm - communicator for new index set
1153: - mode - copy semantics, PETSC_USE_POINTER for no-copy if possible, otherwise PETSC_COPY_VALUES
1155: Output Arguments:
1156: . newis - new IS on comm
1158: Level: advanced
1160: Notes:
1161: It is usually desirable to create a parallel IS and look at the local part when necessary.
1163: This function is useful if serial ISs must be created independently, or to view many
1164: logically independent serial ISs.
1166: The input IS must have the same type on every process.
1168: .seealso: ISSplit()
1169: @*/
1170: PetscErrorCode ISOnComm(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
1171: {
1173: PetscMPIInt match;
1178: MPI_Comm_compare(PetscObjectComm((PetscObject)is),comm,&match);
1179: if (mode != PETSC_COPY_VALUES && (match == MPI_IDENT || match == MPI_CONGRUENT)) {
1180: PetscObjectReference((PetscObject)is);
1181: *newis = is;
1182: } else {
1183: (*is->ops->oncomm)(is,comm,mode,newis);
1184: }
1185: return(0);
1186: }
1188: /*@
1189: ISSetBlockSize - informs an index set that it has a given block size
1191: Logicall Collective on IS
1193: Input Arguments:
1194: + is - index set
1195: - bs - block size
1197: Level: intermediate
1199: .seealso: ISGetBlockSize(), ISCreateBlock()
1200: @*/
1201: PetscErrorCode ISSetBlockSize(IS is,PetscInt bs)
1202: {
1208: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_OUTOFRANGE,"Block size %D, must be positive",bs);
1209: (*is->ops->setblocksize)(is,bs);
1210: return(0);
1211: }
1213: /*@
1214: ISGetBlockSize - Returns the number of elements in a block.
1216: Not Collective
1218: Input Parameter:
1219: . is - the index set
1221: Output Parameter:
1222: . size - the number of elements in a block
1224: Level: intermediate
1226: Concepts: IS^block size
1227: Concepts: index sets^block size
1229: .seealso: ISBlockGetSize(), ISGetSize(), ISCreateBlock(), ISSetBlockSize()
1230: @*/
1231: PetscErrorCode ISGetBlockSize(IS is,PetscInt *size)
1232: {
1236: PetscLayoutGetBlockSize(is->map, size);
1237: return(0);
1238: }
1240: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[])
1241: {
1243: PetscInt len,i;
1244: const PetscInt *ptr;
1247: ISGetSize(is,&len);
1248: ISGetIndices(is,&ptr);
1249: for (i=0; i<len; i++) idx[i] = ptr[i];
1250: ISRestoreIndices(is,&ptr);
1251: return(0);
1252: }
1254: /*MC
1255: ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1256: The users should call ISRestoreIndicesF90() after having looked at the
1257: indices. The user should NOT change the indices.
1259: Synopsis:
1260: ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
1262: Not collective
1264: Input Parameter:
1265: . x - index set
1267: Output Parameters:
1268: + xx_v - the Fortran90 pointer to the array
1269: - ierr - error code
1271: Example of Usage:
1272: .vb
1273: PetscScalar, pointer xx_v(:)
1274: ....
1275: call ISGetIndicesF90(x,xx_v,ierr)
1276: a = xx_v(3)
1277: call ISRestoreIndicesF90(x,xx_v,ierr)
1278: .ve
1280: Notes:
1281: Not yet supported for all F90 compilers.
1283: Level: intermediate
1285: .seealso: ISRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices()
1287: Concepts: index sets^getting indices in f90
1288: Concepts: indices of index set in f90
1290: M*/
1292: /*MC
1293: ISRestoreIndicesF90 - Restores an index set to a usable state after
1294: a call to ISGetIndicesF90().
1296: Synopsis:
1297: ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
1299: Not collective
1301: Input Parameters:
1302: . x - index set
1303: . xx_v - the Fortran90 pointer to the array
1305: Output Parameter:
1306: . ierr - error code
1309: Example of Usage:
1310: .vb
1311: PetscScalar, pointer xx_v(:)
1312: ....
1313: call ISGetIndicesF90(x,xx_v,ierr)
1314: a = xx_v(3)
1315: call ISRestoreIndicesF90(x,xx_v,ierr)
1316: .ve
1318: Notes:
1319: Not yet supported for all F90 compilers.
1321: Level: intermediate
1323: .seealso: ISGetIndicesF90(), ISGetIndices(), ISRestoreIndices()
1325: M*/
1327: /*MC
1328: ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1329: The users should call ISBlockRestoreIndicesF90() after having looked at the
1330: indices. The user should NOT change the indices.
1332: Synopsis:
1333: ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
1335: Not collective
1337: Input Parameter:
1338: . x - index set
1340: Output Parameters:
1341: + xx_v - the Fortran90 pointer to the array
1342: - ierr - error code
1343: Example of Usage:
1344: .vb
1345: PetscScalar, pointer xx_v(:)
1346: ....
1347: call ISBlockGetIndicesF90(x,xx_v,ierr)
1348: a = xx_v(3)
1349: call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1350: .ve
1352: Notes:
1353: Not yet supported for all F90 compilers
1355: Level: intermediate
1357: .seealso: ISBlockRestoreIndicesF90(), ISGetIndices(), ISRestoreIndices(),
1358: ISRestoreIndices()
1360: Concepts: index sets^getting block indices in f90
1361: Concepts: indices of index set in f90
1362: Concepts: block^ indices of index set in f90
1364: M*/
1366: /*MC
1367: ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
1368: a call to ISBlockGetIndicesF90().
1370: Synopsis:
1371: ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
1373: Not Collective
1375: Input Parameters:
1376: + x - index set
1377: - xx_v - the Fortran90 pointer to the array
1379: Output Parameter:
1380: . ierr - error code
1382: Example of Usage:
1383: .vb
1384: PetscScalar, pointer xx_v(:)
1385: ....
1386: call ISBlockGetIndicesF90(x,xx_v,ierr)
1387: a = xx_v(3)
1388: call ISBlockRestoreIndicesF90(x,xx_v,ierr)
1389: .ve
1391: Notes:
1392: Not yet supported for all F90 compilers
1394: Level: intermediate
1396: .seealso: ISBlockGetIndicesF90(), ISGetIndices(), ISRestoreIndices(), ISRestoreIndicesF90()
1398: M*/