Actual source code: isdiff.c

 2:  #include petscis.h
 3:  #include petscbt.h

  7: /*@
  8:    ISDifference - Computes the difference between two index sets.

 10:    Collective on IS

 12:    Input Parameter:
 13: +  is1 - first index, to have items removed from it
 14: -  is2 - index values to be removed

 16:    Output Parameters:
 17: .  isout - is1 - is2

 19:    Notes:
 20:    Negative values are removed from the lists. is2 may have values
 21:    that are not in is1. This requires O(imax-imin) memory and O(imax-imin)
 22:    work, where imin and imax are the bounds on the indices in is1.

 24:    Level: intermediate

 26:    Concepts: index sets^difference
 27:    Concepts: IS^difference

 29: .seealso: ISDestroy(), ISView(), ISSum()

 31: @*/
 32: PetscErrorCode ISDifference(IS is1,IS is2,IS *isout)
 33: {
 35:   PetscInt      i,*i1,*i2,n1,n2,imin,imax,nout,*iout;
 36:   PetscBT  mask;
 37:   MPI_Comm comm;


 44:   ISGetIndices(is1,&i1);
 45:   ISGetLocalSize(is1,&n1);

 47:   /* Create a bit mask array to contain required values */
 48:   if (n1) {
 49:     imin = PETSC_MAX_INT;
 50:     imax = 0;
 51:     for (i=0; i<n1; i++) {
 52:       if (i1[i] < 0) continue;
 53:       imin = PetscMin(imin,i1[i]);
 54:       imax = PetscMax(imax,i1[i]);
 55:     }
 56:   } else {
 57:     imin = imax = 0;
 58:   }
 59:   PetscBTCreate(imax-imin,mask);
 60:   /* Put the values from is1 */
 61:   for (i=0; i<n1; i++) {
 62:     if (i1[i] < 0) continue;
 63:     PetscBTSet(mask,i1[i] - imin);
 64:   }
 65:   ISRestoreIndices(is1,&i1);
 66:   /* Remove the values from is2 */
 67:   ISGetIndices(is2,&i2);
 68:   ISGetLocalSize(is2,&n2);
 69:   for (i=0; i<n2; i++) {
 70:     if (i2[i] < imin || i2[i] > imax) continue;
 71:     PetscBTClear(mask,i2[i] - imin);
 72:   }
 73:   ISRestoreIndices(is2,&i2);
 74: 
 75:   /* Count the number in the difference */
 76:   nout = 0;
 77:   for (i=0; i<imax-imin+1; i++) {
 78:     if (PetscBTLookup(mask,i)) nout++;
 79:   }

 81:   /* create the new IS containing the difference */
 82:   PetscMalloc(nout*sizeof(PetscInt),&iout);
 83:   nout = 0;
 84:   for (i=0; i<imax-imin+1; i++) {
 85:     if (PetscBTLookup(mask,i)) iout[nout++] = i + imin;
 86:   }
 87:   PetscObjectGetComm((PetscObject)is1,&comm);
 88:   ISCreateGeneral(comm,nout,iout,isout);
 89:   PetscFree(iout);

 91:   PetscBTDestroy(mask);
 92:   return(0);
 93: }

 97: /*@
 98:    ISSum - Computes the sum (union) of two index sets in place.

100:    Only sequential version (at the moment)

102:    Input Parameter:
103: +  is1 - index set to be extended
104: -  is2 - index values to be added

106:    Notes:
107:    If n1 and n2 are the sizes of the sets, this takes O(n1+n2) time;
108:    if is2 is a subset of is1, is1 is left unchanged, otherwise is1
109:    is reallocated.
110:    Both index sets need to be sorted on input.

112:    Level: intermediate

114: .seealso: ISDestroy(), ISView(), ISDifference(), ISSum()

116:    Concepts: index sets^union
117:    Concepts: IS^union

119: @*/
120: PetscErrorCode ISSum(IS *is1,IS is2)
121: {
122:   MPI_Comm       comm;
123:   PetscTruth     f;
124:   PetscMPIInt    size;
125:   PetscInt       *i1,*i2,n1,n2,n3, p1,p2, *iout;

131:   PetscObjectGetComm((PetscObject)(*is1),&comm);
132:   MPI_Comm_size(comm,&size);
133:   if (size>1) SETERRQ(PETSC_ERR_SUP,"Currently only for uni-processor IS");

135:   ISSorted(*is1,&f);
136:   if (!f) SETERRQ(PETSC_ERR_ARG_INCOMP,"Arg 1 is not sorted");
137:   ISSorted(is2,&f);
138:   if (!f) SETERRQ(PETSC_ERR_ARG_INCOMP,"Arg 2 is not sorted");

140:   ISGetLocalSize(*is1,&n1);
141:   ISGetLocalSize(is2,&n2);
142:   if (!n2) return(0);
143:   ISGetIndices(*is1,&i1);
144:   ISGetIndices(is2,&i2);

146:   p1 = 0; p2 = 0; n3 = 0;
147:   do {
148:     if (p1==n1) { /* cleanup for is2 */ n3 += n2-p2; break;
149:     } else {
150:       while (p2<n2 && i2[p2]<i1[p1]) {n3++; p2++;}
151:       if (p2==n2) { /* cleanup for is1 */ n3 += n1-p1; break;
152:       } else {
153:         if (i2[p2]==i1[p1]) {n3++; p1++; p2++;}
154:       }
155:     }
156:     if (p2==n2) { /* cleanup for is1 */ n3 += n1-p1; break;
157:     } else {
158:       while (p1<n1 && i1[p1]<i2[p2]) {n3++; p1++;}
159:       if (p1==n1) { /* clean up for is2 */ n3 += n2-p2; break;
160:       } else {
161:         if (i1[p1]==i2[p2]) {n3++; p1++; p2++;}
162:       }
163:     }
164:   } while (p1<n1 || p2<n2);

166:   if (n3==n1) { /* no new elements to be added */
167:     ISRestoreIndices(*is1,&i1);
168:     ISRestoreIndices(is2,&i2);
169:     return(0);
170:   }
171:   PetscMalloc(n3*sizeof(PetscInt),&iout);

173:   p1 = 0; p2 = 0; n3 = 0;
174:   do {
175:     if (p1==n1) { /* cleanup for is2 */
176:       while (p2<n2) iout[n3++] = i2[p2++];
177:       break;
178:     } else {
179:       while (p2<n2 && i2[p2]<i1[p1]) iout[n3++] = i2[p2++];
180:       if (p2==n2) { /* cleanup for is1 */
181:         while (p1<n1) iout[n3++] = i1[p1++];
182:         break;
183:       } else {
184:         if (i2[p2]==i1[p1]) {iout[n3++] = i1[p1++]; p2++;}
185:       }
186:     }
187:     if (p2==n2) { /* cleanup for is1 */
188:       while (p1<n1) iout[n3++] = i1[p1++];
189:       break;
190:     } else {
191:       while (p1<n1 && i1[p1]<i2[p2]) iout[n3++] = i1[p1++];
192:       if (p1==n1) { /* clean up for is2 */
193:         while (p2<n2) iout[n3++] = i2[p2++];
194:         break;
195:       } else {
196:         if (i1[p1]==i2[p2]) {iout[n3++] = i1[p1++]; p2++;}
197:       }
198:     }
199:   } while (p1<n1 || p2<n2);

201:   ISRestoreIndices(*is1,&i1);
202:   ISRestoreIndices(is2,&i2);
203:   ISDestroy(*is1);
204:   ISCreateGeneral(PETSC_COMM_SELF,n3,iout,is1);
205:   PetscFree(iout);

207:   return(0);
208: }