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 src/sys/ctable.h


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

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

 21:    Level: intermediate
 22: @*/
 23: int ISCompressIndicesGeneral(int n,int bs,int imax,const IS is_in[],IS is_out[])
 24: {
 25:   int                ierr,isz,len,i,j,*idx,ival,Nbs;
 26: #if defined (PETSC_USE_CTABLE)
 27:   PetscTable         gid1_lid1;
 28:   int                tt, gid1, *nidx;
 29:   PetscTablePosition tpos;
 30: #else
 31:   int                *nidx;
 32:   PetscBT            table;
 33: #endif

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

 94: int ISCompressIndicesSorted(int n,int bs,int imax,const IS is_in[],IS is_out[])
 95: {
 96:   int          ierr,i,j,k,val,len,*idx,*nidx,*idx_local;
 97:   PetscTruth   flg;
 98: #if defined (PETSC_USE_CTABLE)
 99:   int maxsz;
100: #else
101:   int Nbs=n/bs;
102: #endif
104:   for (i=0; i<imax; i++) {
105:     ISSorted(is_in[i],&flg);
106:     if (!flg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Indices are not sorted");
107:   }
108: #if defined (PETSC_USE_CTABLE)
109:   /* Now check max size */
110:   for (i=0,maxsz=0; i<imax; i++) {
111:     ISGetIndices(is_in[i],&idx);
112:     ISGetLocalSize(is_in[i],&len);
113:     if (len%bs !=0) SETERRQ(1,"Indices are not block ordered");
114:     len = len/bs; /* The reduced index size */
115:     if (len > maxsz) maxsz = len;
116:   }
117:   PetscMalloc((maxsz+1)*sizeof(int),&nidx);
118: #else
119:   PetscMalloc((Nbs+1)*sizeof(int),&nidx);
120: #endif
121:   /* Now check if the indices are in block order */
122:   for (i=0; i<imax; i++) {
123:     ISGetIndices(is_in[i],&idx);
124:     ISGetLocalSize(is_in[i],&len);
125:     if (len%bs !=0) SETERRQ(1,"Indices are not block ordered");

127:     len = len/bs; /* The reduced index size */
128:     idx_local = idx;
129:     for (j=0; j<len ; j++) {
130:       val = idx_local[0];
131:       if (val%bs != 0) SETERRQ(1,"Indices are not block ordered");
132:       for (k=0; k<bs; k++) {
133:         if (val+k != idx_local[k]) SETERRQ(1,"Indices are not block ordered");
134:       }
135:       nidx[j] = val/bs;
136:       idx_local +=bs;
137:     }
138:     ISRestoreIndices(is_in[i],&idx);
139:     ISCreateGeneral(PETSC_COMM_SELF,len,nidx,(is_out+i));
140:   }
141:   PetscFree(nidx);

143:   return(0);
144: }

148: int ISExpandIndicesGeneral(int n,int bs,int imax,const IS is_in[],IS is_out[])
149: {
150:   int          ierr,len,i,j,k,*idx,*nidx;
151: #if defined (PETSC_USE_CTABLE)
152:   int          maxsz;
153: #else
154:   int          Nbs;
155: #endif

158: #if defined (PETSC_USE_CTABLE)
159:   /* Now check max size */
160:   for (i=0,maxsz=0; i<imax; i++) {
161:     ISGetIndices(is_in[i],&idx);
162:     ISGetLocalSize(is_in[i],&len);
163:     if (len*bs > maxsz) maxsz = len*bs;
164:   }
165:   PetscMalloc((maxsz+1)*sizeof(int),&nidx);
166: #else
167:   Nbs = n/bs;
168:   PetscMalloc((Nbs*bs+1)*sizeof(int),&nidx);
169: #endif

171:   for (i=0; i<imax; i++) {
172:     ISGetIndices(is_in[i],&idx);
173:     ISGetLocalSize(is_in[i],&len);
174:     for (j=0; j<len ; ++j){
175:       for (k=0; k<bs; k++)
176:         nidx[j*bs+k] = idx[j]*bs+k;
177:     }
178:     ISRestoreIndices(is_in[i],&idx);
179:     ISCreateGeneral(PETSC_COMM_SELF,len*bs,nidx,is_out+i);
180:   }
181:   PetscFree(nidx);
182:   return(0);
183: }