Actual source code: isblock.c

  2: /* Routines to be used by MatIncreaseOverlap() for BAIJ and SBAIJ matrices */
  3: #include <petscis.h> 
  4: #include <petscbt.h>
  5: #include <petscctable.h>


 10: /*@C
 11:    ISCompressIndicesGeneral - convert the indices into block indices
 12:    Input Parameters:
 13: +  n - maximum possible length of the index set
 14: .  nkeys - expected number of keys when PETSC_USE_CTABLE
 15: .  bs - the size of block 
 16: .  imax - the number of index sets
 17: -  is_in - the non-blocked array of index sets 

 19:    Output Parameter:
 20: .  is_out - the blocked new index set

 22:    Level: intermediate

 24: .seealso: ISExpandIndicesGeneral()
 25: @*/
 26: PetscErrorCode  ISCompressIndicesGeneral(PetscInt n,PetscInt nkeys,PetscInt bs,PetscInt imax,const IS is_in[],IS is_out[])
 27: {
 28:   PetscErrorCode     ierr;
 29:   PetscInt           isz,len,i,j,ival,Nbs;
 30:   const PetscInt     *idx;
 31: #if defined (PETSC_USE_CTABLE)
 32:   PetscTable         gid1_lid1;
 33:   PetscInt           tt, gid1, *nidx,Nkbs;
 34:   PetscTablePosition tpos;
 35: #else
 36:   PetscInt           *nidx;
 37:   PetscBT            table;
 38: #endif

 41:   Nbs =n/bs;
 42: #if defined (PETSC_USE_CTABLE)
 43:   Nkbs = nkeys/bs;
 44:   PetscTableCreate(Nkbs,&gid1_lid1);
 45: #else
 46:   PetscMalloc(Nbs*sizeof(PetscInt),&nidx);
 47:   PetscBTCreate(Nbs,table);
 48: #endif
 49:   for (i=0; i<imax; i++) {
 50:     isz  = 0;
 51: #if defined (PETSC_USE_CTABLE)
 52:     PetscTableRemoveAll(gid1_lid1);
 53: #else
 54:     PetscBTMemzero(Nbs,table);
 55: #endif
 56:     ISGetIndices(is_in[i],&idx);
 57:     ISGetLocalSize(is_in[i],&len);
 58:     for (j=0; j<len ; j++) {
 59:       ival = idx[j]/bs; /* convert the indices into block indices */
 60: #if defined (PETSC_USE_CTABLE)
 61:       PetscTableFind(gid1_lid1,ival+1,&tt);
 62:       if (!tt) {
 63:         PetscTableAdd(gid1_lid1,ival+1,isz+1);
 64:         isz++;
 65:       }
 66: #else
 67:       if (ival>Nbs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
 68:       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
 69: #endif
 70:     }
 71:     ISRestoreIndices(is_in[i],&idx);
 72: 
 73: #if defined (PETSC_USE_CTABLE)
 74:     PetscMalloc(isz*sizeof(PetscInt),&nidx);
 75:     PetscTableGetHeadPosition(gid1_lid1,&tpos);
 76:     j = 0;
 77:     while (tpos) {
 78:       PetscTableGetNext(gid1_lid1,&tpos,&gid1,&tt);
 79:       if (tt-- > isz) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than array-dim"); }
 80:       nidx[tt] = gid1 - 1;
 81:       j++;
 82:     }
 83:     if (j != isz) { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"table error: jj != isz"); }
 84:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_OWN_POINTER,(is_out+i));
 85: #else
 86:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,(is_out+i));
 87: #endif
 88:   }
 89: #if defined (PETSC_USE_CTABLE)
 90:   PetscTableDestroy(&gid1_lid1);
 91: #else
 92:   PetscBTDestroy(table);
 93:   PetscFree(nidx);
 94: #endif
 95:   return(0);
 96: }

100: PetscErrorCode  ISCompressIndicesSorted(PetscInt n,PetscInt bs,PetscInt imax,const IS is_in[],IS is_out[])
101: {
103:   PetscInt       i,j,k,val,len,*nidx,bbs;
104:   const PetscInt *idx,*idx_local;
105:   PetscBool      flg,isblock;
106: #if defined (PETSC_USE_CTABLE)
107:   PetscInt       maxsz;
108: #else
109:   PetscInt       Nbs=n/bs;
110: #endif

113:   for (i=0; i<imax; i++) {
114:     ISSorted(is_in[i],&flg);
115:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Indices are not sorted");
116:   }

118: #if defined (PETSC_USE_CTABLE)
119:   /* Now check max size */
120:   for (i=0,maxsz=0; i<imax; i++) {
121:     ISGetLocalSize(is_in[i],&len);
122:     if (len%bs !=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
123:     len = len/bs; /* The reduced index size */
124:     if (len > maxsz) maxsz = len;
125:   }
126:   PetscMalloc(maxsz*sizeof(PetscInt),&nidx);
127: #else
128:   PetscMalloc(Nbs*sizeof(PetscInt),&nidx);
129: #endif
130:   /* Now check if the indices are in block order */
131:   for (i=0; i<imax; i++) {
132:     ISGetLocalSize(is_in[i],&len);

134:     /* special case where IS is already block IS of the correct size */
135:     PetscTypeCompare((PetscObject)is_in[i],ISBLOCK,&isblock);
136:     if (isblock) {
137:       ISBlockGetLocalSize(is_in[i],&bbs);
138:       if (bs == bbs) {
139:         len = len/bs;
140:         ISBlockGetIndices(is_in[i],&idx);
141:         ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,(is_out+i));
142:         ISBlockRestoreIndices(is_in[i],&idx);
143:         continue;
144:       }
145:     }
146:     ISGetIndices(is_in[i],&idx);
147:     if (len%bs !=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");

149:     len = len/bs; /* The reduced index size */
150:     idx_local = idx;
151:     for (j=0; j<len ; j++) {
152:       val = idx_local[0];
153:       if (val%bs != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
154:       for (k=0; k<bs; k++) {
155:         if (val+k != idx_local[k]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
156:       }
157:       nidx[j] = val/bs;
158:       idx_local +=bs;
159:     }
160:     ISRestoreIndices(is_in[i],&idx);
161:     ISCreateGeneral(PETSC_COMM_SELF,len,nidx,PETSC_COPY_VALUES,(is_out+i));
162:   }
163:   PetscFree(nidx);

165:   return(0);
166: }

170: /*@C
171:    ISExpandIndicesGeneral - convert the indices into non-block indices
172:    Input Parameters:
173: +  n - the length of the index set   (not being used)
174: .  nkeys - expected number of keys when PETSC_USE_CTABLE (not being used)
175: .  bs - the size of block 
176: .  imax - the number of index sets
177: -  is_in - the blocked array of index sets 

179:    Output Parameter:
180: .  is_out - the non-blocked new index set

182:    Level: intermediate

184: .seealso: ISCompressIndicesGeneral()
185: @*/
186: PetscErrorCode  ISExpandIndicesGeneral(PetscInt n,PetscInt nkeys,PetscInt bs,PetscInt imax,const IS is_in[],IS is_out[])
187: {
189:   PetscInt       len,i,j,k,*nidx;
190:   const PetscInt *idx;
191:   PetscInt       maxsz;

194:   /* Check max size of is_in[] */
195:   maxsz=0;
196:   for (i=0; i<imax; i++) {
197:     ISGetLocalSize(is_in[i],&len);
198:     if (len > maxsz) maxsz = len;
199:   }
200:   PetscMalloc(maxsz*bs*sizeof(PetscInt),&nidx);

202:   for (i=0; i<imax; i++) {
203:     ISGetLocalSize(is_in[i],&len);
204:     ISGetIndices(is_in[i],&idx);
205:     for (j=0; j<len ; ++j){
206:       for (k=0; k<bs; k++) nidx[j*bs+k] = idx[j]*bs+k;
207:     }
208:     ISRestoreIndices(is_in[i],&idx);
209:     ISCreateGeneral(PETSC_COMM_SELF,len*bs,nidx,PETSC_COPY_VALUES,is_out+i);
210:   }
211:   PetscFree(nidx);
212:   return(0);
213: }