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: PetscErrorCode ISCompressIndicesGeneral(PetscInt n,PetscInt bs,PetscInt imax,const IS is_in[],IS is_out[])
24: {
25: PetscErrorCode ierr;
26: PetscInt isz,len,i,j,*idx,ival,Nbs;
27: #if defined (PETSC_USE_CTABLE)
28: PetscTable gid1_lid1;
29: PetscInt tt, gid1, *nidx;
30: PetscTablePosition tpos;
31: #else
32: PetscInt *nidx;
33: PetscBT table;
34: #endif
37: Nbs =n/bs;
38: #if defined (PETSC_USE_CTABLE)
39: PetscTableCreate(Nbs,&gid1_lid1);
40: #else
41: PetscMalloc(Nbs*sizeof(PetscInt),&nidx);
42: PetscBTCreate(Nbs,table);
43: #endif
44: for (i=0; i<imax; i++) {
45: isz = 0;
46: #if defined (PETSC_USE_CTABLE)
47: PetscTableRemoveAll(gid1_lid1);
48: #else
49: PetscBTMemzero(Nbs,table);
50: #endif
51: ISGetIndices(is_in[i],&idx);
52: ISGetLocalSize(is_in[i],&len);
53: for (j=0; j<len ; j++) {
54: ival = idx[j]/bs; /* convert the indices into block indices */
55: #if defined (PETSC_USE_CTABLE)
56: PetscTableFind(gid1_lid1,ival+1,&tt);
57: if (!tt) {
58: PetscTableAdd(gid1_lid1,ival+1,isz+1);
59: isz++;
60: }
61: #else
62: if (ival>Nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
63: if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
64: #endif
65: }
66: ISRestoreIndices(is_in[i],&idx);
67: #if defined (PETSC_USE_CTABLE)
68: PetscMalloc(isz*sizeof(PetscInt),&nidx);
69: PetscTableGetHeadPosition(gid1_lid1,&tpos);
70: j = 0;
71: while (tpos) {
72: PetscTableGetNext(gid1_lid1,&tpos,&gid1,&tt);
73: if (tt-- > isz) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than array-dim"); }
74: nidx[tt] = gid1 - 1;
75: j++;
76: }
77: if (j != isz) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"table error: jj != isz"); }
78: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is_out+i));
79: PetscFree(nidx);
80: #else
81: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is_out+i));
82: #endif
83: }
84: #if defined (PETSC_USE_CTABLE)
85: PetscTableDelete(gid1_lid1);
86: #else
87: PetscBTDestroy(table);
88: PetscFree(nidx);
89: #endif
90: return(0);
91: }
95: PetscErrorCode ISCompressIndicesSorted(PetscInt n,PetscInt bs,PetscInt imax,const IS is_in[],IS is_out[])
96: {
98: PetscInt i,j,k,val,len,*idx,*nidx,*idx_local;
99: PetscTruth flg;
100: #if defined (PETSC_USE_CTABLE)
101: PetscInt maxsz;
102: #else
103: PetscInt Nbs=n/bs;
104: #endif
106: for (i=0; i<imax; i++) {
107: ISSorted(is_in[i],&flg);
108: if (!flg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Indices are not sorted");
109: }
110: #if defined (PETSC_USE_CTABLE)
111: /* Now check max size */
112: for (i=0,maxsz=0; i<imax; i++) {
113: ISGetLocalSize(is_in[i],&len);
114: if (len%bs !=0) SETERRQ(PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
115: len = len/bs; /* The reduced index size */
116: if (len > maxsz) maxsz = len;
117: }
118: PetscMalloc(maxsz*sizeof(PetscInt),&nidx);
119: #else
120: PetscMalloc(Nbs*sizeof(PetscInt),&nidx);
121: #endif
122: /* Now check if the indices are in block order */
123: for (i=0; i<imax; i++) {
124: ISGetIndices(is_in[i],&idx);
125: ISGetLocalSize(is_in[i],&len);
126: if (len%bs !=0) SETERRQ(PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
128: len = len/bs; /* The reduced index size */
129: idx_local = idx;
130: for (j=0; j<len ; j++) {
131: val = idx_local[0];
132: if (val%bs != 0) SETERRQ(PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
133: for (k=0; k<bs; k++) {
134: if (val+k != idx_local[k]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Indices are not block ordered");
135: }
136: nidx[j] = val/bs;
137: idx_local +=bs;
138: }
139: ISRestoreIndices(is_in[i],&idx);
140: ISCreateGeneral(PETSC_COMM_SELF,len,nidx,(is_out+i));
141: }
142: PetscFree(nidx);
144: return(0);
145: }
149: PetscErrorCode ISExpandIndicesGeneral(PetscInt n,PetscInt bs,PetscInt imax,const IS is_in[],IS is_out[])
150: {
152: PetscInt len,i,j,k,*idx,*nidx;
153: #if defined (PETSC_USE_CTABLE)
154: PetscInt maxsz;
155: #else
156: PetscInt Nbs;
157: #endif
160: #if defined (PETSC_USE_CTABLE)
161: /* Now check max size */
162: for (i=0,maxsz=0; i<imax; i++) {
163: ISGetIndices(is_in[i],&idx);
164: ISGetLocalSize(is_in[i],&len);
165: if (len*bs > maxsz) maxsz = len*bs;
166: }
167: PetscMalloc(maxsz*sizeof(PetscInt),&nidx);
168: #else
169: Nbs = n/bs;
170: PetscMalloc(Nbs*bs*sizeof(PetscInt),&nidx);
171: #endif
173: for (i=0; i<imax; i++) {
174: ISGetIndices(is_in[i],&idx);
175: ISGetLocalSize(is_in[i],&len);
176: for (j=0; j<len ; ++j){
177: for (k=0; k<bs; k++)
178: nidx[j*bs+k] = idx[j]*bs+k;
179: }
180: ISRestoreIndices(is_in[i],&idx);
181: ISCreateGeneral(PETSC_COMM_SELF,len*bs,nidx,is_out+i);
182: }
183: PetscFree(nidx);
184: return(0);
185: }