Actual source code: isblock.c

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

  6: /*@
  7:    ISCompressIndicesGeneral - convert the indices of an array of `IS` into an array of `ISGENERAL` of block indices

  9:    Input Parameters:
 10: +    n - maximum possible length of the index set
 11: .    nkeys - expected number of keys when using `PETSC_USE_CTABLE`
 12: .    bs - the size of block
 13: .    imax - the number of index sets
 14: -    is_in - the non-blocked array of index sets

 16:    Output Parameter:
 17: .    is_out - the blocked new index set, as `ISGENERAL`, not as `ISBLOCK`

 19:    Level: intermediate

 21: .seealso: [](sec_scatter), `IS`, `ISGENERAL`, `ISExpandIndicesGeneral()`
 22: @*/
 23: PetscErrorCode ISCompressIndicesGeneral(PetscInt n, PetscInt nkeys, PetscInt bs, PetscInt imax, const IS is_in[], IS is_out[])
 24: {
 25:   PetscInt        isz, len, i, j, ival;
 26:   const PetscInt *idx;
 27: #if defined(PETSC_USE_CTABLE)
 28:   PetscHMapI    gid1_lid1 = NULL;
 29:   PetscInt      tt, gid1, *nidx;
 30:   PetscHashIter tpos;
 31: #else
 32:   PetscInt *nidx;
 33:   PetscInt  Nbs;
 34:   PetscBT   table;
 35: #endif

 37:   PetscFunctionBegin;
 38: #if defined(PETSC_USE_CTABLE)
 39:   PetscCall(PetscHMapICreateWithSize(nkeys / bs, &gid1_lid1));
 40: #else
 41:   Nbs = n / bs;
 42:   PetscCall(PetscMalloc1(Nbs, &nidx));
 43:   PetscCall(PetscBTCreate(Nbs, &table));
 44: #endif
 45:   for (i = 0; i < imax; i++) {
 46:     isz = 0;
 47: #if defined(PETSC_USE_CTABLE)
 48:     PetscCall(PetscHMapIClear(gid1_lid1));
 49: #else
 50:     PetscCall(PetscBTMemzero(Nbs, table));
 51: #endif
 52:     PetscCall(ISGetIndices(is_in[i], &idx));
 53:     PetscCall(ISGetLocalSize(is_in[i], &len));
 54:     for (j = 0; j < len; j++) {
 55:       ival = idx[j] / bs; /* convert the indices into block indices */
 56: #if defined(PETSC_USE_CTABLE)
 57:       PetscCall(PetscHMapIGetWithDefault(gid1_lid1, ival + 1, 0, &tt));
 58:       if (!tt) {
 59:         PetscCall(PetscHMapISet(gid1_lid1, ival + 1, isz + 1));
 60:         isz++;
 61:       }
 62: #else
 63:       PetscCheck(ival <= Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index greater than mat-dim");
 64:       if (!PetscBTLookupSet(table, ival)) nidx[isz++] = ival;
 65: #endif
 66:     }
 67:     PetscCall(ISRestoreIndices(is_in[i], &idx));

 69: #if defined(PETSC_USE_CTABLE)
 70:     PetscCall(PetscMalloc1(isz, &nidx));
 71:     PetscHashIterBegin(gid1_lid1, tpos);
 72:     j = 0;
 73:     while (!PetscHashIterAtEnd(gid1_lid1, tpos)) {
 74:       PetscHashIterGetKey(gid1_lid1, tpos, gid1);
 75:       PetscHashIterGetVal(gid1_lid1, tpos, tt);
 76:       PetscCheck(tt-- <= isz, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index greater than array-dim");
 77:       nidx[tt] = gid1 - 1;
 78:       j++;
 79:       PetscHashIterNext(gid1_lid1, tpos);
 80:     }
 81:     PetscCheck(j == isz, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "table error: jj != isz");
 82:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_in[i]), isz, nidx, PETSC_OWN_POINTER, (is_out + i)));
 83: #else
 84:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_in[i]), isz, nidx, PETSC_COPY_VALUES, (is_out + i)));
 85: #endif
 86:   }
 87: #if defined(PETSC_USE_CTABLE)
 88:   PetscCall(PetscHMapIDestroy(&gid1_lid1));
 89: #else
 90:   PetscCall(PetscBTDestroy(&table));
 91:   PetscCall(PetscFree(nidx));
 92: #endif
 93:   PetscFunctionReturn(PETSC_SUCCESS);
 94: }

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

107:   PetscFunctionBegin;
108:   for (i = 0; i < imax; i++) {
109:     PetscCall(ISSorted(is_in[i], &flg));
110:     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Indices are not sorted");
111:   }

113: #if defined(PETSC_USE_CTABLE)
114:   /* Now check max size */
115:   for (i = 0, maxsz = 0; i < imax; i++) {
116:     PetscCall(ISGetLocalSize(is_in[i], &len));
117:     PetscCheck(len % bs == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Indices are not block ordered");
118:     len = len / bs; /* The reduced index size */
119:     if (len > maxsz) maxsz = len;
120:   }
121:   PetscCall(PetscMalloc1(maxsz, &nidx));
122: #else
123:   PetscCall(PetscMalloc1(Nbs, &nidx));
124: #endif
125:   /* Now check if the indices are in block order */
126:   for (i = 0; i < imax; i++) {
127:     PetscCall(ISGetLocalSize(is_in[i], &len));

129:     /* special case where IS is already block IS of the correct size */
130:     PetscCall(PetscObjectTypeCompare((PetscObject)is_in[i], ISBLOCK, &isblock));
131:     if (isblock) {
132:       PetscCall(ISBlockGetLocalSize(is_in[i], &bbs));
133:       if (bs == bbs) {
134:         len = len / bs;
135:         PetscCall(ISBlockGetIndices(is_in[i], &idx));
136:         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, len, idx, PETSC_COPY_VALUES, (is_out + i)));
137:         PetscCall(ISBlockRestoreIndices(is_in[i], &idx));
138:         continue;
139:       }
140:     }
141:     PetscCall(ISGetIndices(is_in[i], &idx));
142:     PetscCheck(len % bs == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Indices are not block ordered");

144:     len       = len / bs; /* The reduced index size */
145:     idx_local = idx;
146:     for (j = 0; j < len; j++) {
147:       val = idx_local[0];
148:       PetscCheck(val % bs == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Indices are not block ordered");
149:       for (k = 0; k < bs; k++) PetscCheck(val + k == idx_local[k], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Indices are not block ordered");
150:       nidx[j] = val / bs;
151:       idx_local += bs;
152:     }
153:     PetscCall(ISRestoreIndices(is_in[i], &idx));
154:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_in[i]), len, nidx, PETSC_COPY_VALUES, (is_out + i)));
155:   }
156:   PetscCall(PetscFree(nidx));
157:   PetscFunctionReturn(PETSC_SUCCESS);
158: }

160: /*@C
161:    ISExpandIndicesGeneral - convert the indices of an array `IS` into non-block indices in an array of `ISGENERAL`

163:    Input Parameters:
164: +    n - the length of the index set (not being used)
165: .    nkeys - expected number of keys when `PETSC_USE_CTABLE` is used
166: .    bs - the size of block
167: .    imax - the number of index sets
168: -    is_in - the blocked array of index sets

170:    Output Parameter:
171: .    is_out - the non-blocked new index set, as `ISGENERAL`

173:    Level: intermediate

175: .seealso: [](sec_scatter), `IS`, `ISGENERAL`, `ISCompressIndicesGeneral()`
176: @*/
177: PetscErrorCode ISExpandIndicesGeneral(PetscInt n, PetscInt nkeys, PetscInt bs, PetscInt imax, const IS is_in[], IS is_out[])
178: {
179:   PetscInt        len, i, j, k, *nidx;
180:   const PetscInt *idx;
181:   PetscInt        maxsz;

183:   PetscFunctionBegin;
184:   /* Check max size of is_in[] */
185:   maxsz = 0;
186:   for (i = 0; i < imax; i++) {
187:     PetscCall(ISGetLocalSize(is_in[i], &len));
188:     if (len > maxsz) maxsz = len;
189:   }
190:   PetscCall(PetscMalloc1(maxsz * bs, &nidx));

192:   for (i = 0; i < imax; i++) {
193:     PetscCall(ISGetLocalSize(is_in[i], &len));
194:     PetscCall(ISGetIndices(is_in[i], &idx));
195:     for (j = 0; j < len; ++j) {
196:       for (k = 0; k < bs; k++) nidx[j * bs + k] = idx[j] * bs + k;
197:     }
198:     PetscCall(ISRestoreIndices(is_in[i], &idx));
199:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_in[i]), len * bs, nidx, PETSC_COPY_VALUES, is_out + i));
200:   }
201:   PetscCall(PetscFree(nidx));
202:   PetscFunctionReturn(PETSC_SUCCESS);
203: }