1 #include <../src/mat/impls/aij/seq/aij.h> 2 #include <../src/mat/impls/baij/seq/baij.h> 3 #include <../src/mat/impls/sell/seq/sell.h> 4 #include <petsc/private/isimpl.h> 5 6 /* 7 This routine is shared by SeqAIJ and SeqBAIJ matrices, 8 since it operators only on the nonzero structure of the elements or blocks. 9 */ 10 PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) 11 { 12 PetscInt bs, nis = iscoloring->n, m = mat->rmap->n; 13 PetscBool isBAIJ, isSELL; 14 15 PetscFunctionBegin; 16 /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */ 17 PetscCall(MatGetBlockSize(mat, &bs)); 18 PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQBAIJ, &isBAIJ)); 19 PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSELL, &isSELL)); 20 if (isBAIJ) { 21 c->brows = m; 22 c->bcols = 1; 23 } else { /* seqaij matrix */ 24 /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 25 PetscReal mem; 26 PetscInt nz, brows, bcols; 27 if (isSELL) { 28 Mat_SeqSELL *spA = (Mat_SeqSELL *)mat->data; 29 nz = spA->nz; 30 } else { 31 Mat_SeqAIJ *spA = (Mat_SeqAIJ *)mat->data; 32 nz = spA->nz; 33 } 34 35 bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 36 mem = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt); 37 bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar))); 38 brows = 1000 / bcols; 39 if (bcols > nis) bcols = nis; 40 if (brows == 0 || brows > m) brows = m; 41 c->brows = brows; 42 c->bcols = bcols; 43 } 44 45 c->M = mat->rmap->N / bs; /* set total rows, columns and local rows */ 46 c->N = mat->cmap->N / bs; 47 c->m = mat->rmap->N / bs; 48 c->rstart = 0; 49 c->ncolors = nis; 50 c->ctype = iscoloring->ctype; 51 PetscFunctionReturn(0); 52 } 53 54 /* 55 Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into Jacobian for improved cache performance 56 Input Parameters: 57 + mat - the matrix containing the nonzero structure of the Jacobian 58 . color - the coloring context 59 - nz - number of local non-zeros in mat 60 */ 61 PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat, MatFDColoring c, PetscInt nz) 62 { 63 PetscInt i, j, nrows, nbcols, brows = c->brows, bcols = c->bcols, mbs = c->m, nis = c->ncolors; 64 PetscInt *color_start, *row_start, *nrows_new, nz_new, row_end; 65 66 PetscFunctionBegin; 67 if (brows < 1 || brows > mbs) brows = mbs; 68 PetscCall(PetscMalloc2(bcols + 1, &color_start, bcols, &row_start)); 69 PetscCall(PetscCalloc1(nis, &nrows_new)); 70 PetscCall(PetscMalloc1(bcols * mat->rmap->n, &c->dy)); 71 72 nz_new = 0; 73 nbcols = 0; 74 color_start[bcols] = 0; 75 76 if (c->htype[0] == 'd') { /* ---- c->htype == 'ds', use MatEntry --------*/ 77 MatEntry *Jentry_new, *Jentry = c->matentry; 78 79 PetscCall(PetscMalloc1(nz, &Jentry_new)); 80 for (i = 0; i < nis; i += bcols) { /* loop over colors */ 81 if (i + bcols > nis) { 82 color_start[nis - i] = color_start[bcols]; 83 bcols = nis - i; 84 } 85 86 color_start[0] = color_start[bcols]; 87 for (j = 0; j < bcols; j++) { 88 color_start[j + 1] = c->nrows[i + j] + color_start[j]; 89 row_start[j] = 0; 90 } 91 92 row_end = brows; 93 if (row_end > mbs) row_end = mbs; 94 95 while (row_end <= mbs) { /* loop over block rows */ 96 for (j = 0; j < bcols; j++) { /* loop over block columns */ 97 nrows = c->nrows[i + j]; 98 nz = color_start[j]; 99 while (row_start[j] < nrows) { 100 if (Jentry[nz].row >= row_end) { 101 color_start[j] = nz; 102 break; 103 } else { /* copy Jentry[nz] to Jentry_new[nz_new] */ 104 Jentry_new[nz_new].row = Jentry[nz].row + j * mbs; /* index in dy-array */ 105 Jentry_new[nz_new].col = Jentry[nz].col; 106 Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 107 nz_new++; 108 nz++; 109 row_start[j]++; 110 } 111 } 112 } 113 if (row_end == mbs) break; 114 row_end += brows; 115 if (row_end > mbs) row_end = mbs; 116 } 117 nrows_new[nbcols++] = nz_new; 118 } 119 PetscCall(PetscFree(Jentry)); 120 c->matentry = Jentry_new; 121 } else { /* --------- c->htype == 'wp', use MatEntry2 ------------------*/ 122 MatEntry2 *Jentry2_new, *Jentry2 = c->matentry2; 123 124 PetscCall(PetscMalloc1(nz, &Jentry2_new)); 125 for (i = 0; i < nis; i += bcols) { /* loop over colors */ 126 if (i + bcols > nis) { 127 color_start[nis - i] = color_start[bcols]; 128 bcols = nis - i; 129 } 130 131 color_start[0] = color_start[bcols]; 132 for (j = 0; j < bcols; j++) { 133 color_start[j + 1] = c->nrows[i + j] + color_start[j]; 134 row_start[j] = 0; 135 } 136 137 row_end = brows; 138 if (row_end > mbs) row_end = mbs; 139 140 while (row_end <= mbs) { /* loop over block rows */ 141 for (j = 0; j < bcols; j++) { /* loop over block columns */ 142 nrows = c->nrows[i + j]; 143 nz = color_start[j]; 144 while (row_start[j] < nrows) { 145 if (Jentry2[nz].row >= row_end) { 146 color_start[j] = nz; 147 break; 148 } else { /* copy Jentry2[nz] to Jentry2_new[nz_new] */ 149 Jentry2_new[nz_new].row = Jentry2[nz].row + j * mbs; /* index in dy-array */ 150 Jentry2_new[nz_new].valaddr = Jentry2[nz].valaddr; 151 nz_new++; 152 nz++; 153 row_start[j]++; 154 } 155 } 156 } 157 if (row_end == mbs) break; 158 row_end += brows; 159 if (row_end > mbs) row_end = mbs; 160 } 161 nrows_new[nbcols++] = nz_new; 162 } 163 PetscCall(PetscFree(Jentry2)); 164 c->matentry2 = Jentry2_new; 165 } /* ---------------------------------------------*/ 166 167 PetscCall(PetscFree2(color_start, row_start)); 168 169 for (i = nbcols - 1; i > 0; i--) nrows_new[i] -= nrows_new[i - 1]; 170 PetscCall(PetscFree(c->nrows)); 171 c->nrows = nrows_new; 172 PetscFunctionReturn(0); 173 } 174 175 PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) 176 { 177 PetscInt i, n, nrows, mbs = c->m, j, k, m, ncols, col, nis = iscoloring->n, *rowhit, bs, bs2, *spidx, nz, tmp; 178 const PetscInt *is, *row, *ci, *cj; 179 PetscBool isBAIJ, isSELL; 180 const PetscScalar *A_val; 181 PetscScalar **valaddrhit; 182 MatEntry *Jentry; 183 MatEntry2 *Jentry2; 184 185 PetscFunctionBegin; 186 PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa)); 187 188 PetscCall(MatGetBlockSize(mat, &bs)); 189 PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQBAIJ, &isBAIJ)); 190 PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSELL, &isSELL)); 191 if (isBAIJ) { 192 Mat_SeqBAIJ *spA = (Mat_SeqBAIJ *)mat->data; 193 194 A_val = spA->a; 195 nz = spA->nz; 196 } else if (isSELL) { 197 Mat_SeqSELL *spA = (Mat_SeqSELL *)mat->data; 198 199 A_val = spA->val; 200 nz = spA->nz; 201 bs = 1; /* only bs=1 is supported for SeqSELL matrix */ 202 } else { 203 Mat_SeqAIJ *spA = (Mat_SeqAIJ *)mat->data; 204 205 A_val = spA->a; 206 nz = spA->nz; 207 bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 208 } 209 210 PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns)); 211 PetscCall(PetscMalloc1(nis, &c->nrows)); /* nrows is freed separately from ncolumns and columns */ 212 213 if (c->htype[0] == 'd') { 214 PetscCall(PetscMalloc1(nz, &Jentry)); 215 c->matentry = Jentry; 216 } else if (c->htype[0] == 'w') { 217 PetscCall(PetscMalloc1(nz, &Jentry2)); 218 c->matentry2 = Jentry2; 219 } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported"); 220 221 if (isBAIJ) { 222 PetscCall(MatGetColumnIJ_SeqBAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL)); 223 } else if (isSELL) { 224 PetscCall(MatGetColumnIJ_SeqSELL_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL)); 225 } else { 226 PetscCall(MatGetColumnIJ_SeqAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL)); 227 } 228 229 PetscCall(PetscCalloc1(c->m, &rowhit)); 230 PetscCall(PetscMalloc1(c->m, &valaddrhit)); 231 232 nz = 0; 233 for (i = 0; i < nis; i++) { /* loop over colors */ 234 PetscCall(ISGetLocalSize(c->isa[i], &n)); 235 PetscCall(ISGetIndices(c->isa[i], &is)); 236 237 c->ncolumns[i] = n; 238 c->columns[i] = (PetscInt *)is; 239 /* note: we know that c->isa is going to be around as long at the c->columns values */ 240 PetscCall(ISRestoreIndices(c->isa[i], &is)); 241 242 /* fast, crude version requires O(N*N) work */ 243 bs2 = bs * bs; 244 nrows = 0; 245 for (j = 0; j < n; j++) { /* loop over columns */ 246 col = is[j]; 247 tmp = ci[col]; 248 row = cj + tmp; 249 m = ci[col + 1] - tmp; 250 nrows += m; 251 for (k = 0; k < m; k++) { /* loop over columns marking them in rowhit */ 252 rowhit[*row] = col + 1; 253 valaddrhit[*row++] = (PetscScalar *)&A_val[bs2 * spidx[tmp + k]]; 254 } 255 } 256 c->nrows[i] = nrows; /* total num of rows for this color */ 257 258 if (c->htype[0] == 'd') { 259 for (j = 0; j < mbs; j++) { /* loop over rows */ 260 if (rowhit[j]) { 261 Jentry[nz].row = j; /* local row index */ 262 Jentry[nz].col = rowhit[j] - 1; /* local column index */ 263 Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 264 nz++; 265 rowhit[j] = 0.0; /* zero rowhit for reuse */ 266 } 267 } 268 } else { /* c->htype == 'wp' */ 269 for (j = 0; j < mbs; j++) { /* loop over rows */ 270 if (rowhit[j]) { 271 Jentry2[nz].row = j; /* local row index */ 272 Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 273 nz++; 274 rowhit[j] = 0.0; /* zero rowhit for reuse */ 275 } 276 } 277 } 278 } 279 280 if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 281 PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz)); 282 } 283 284 if (isBAIJ) { 285 PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL)); 286 PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy)); 287 } else if (isSELL) { 288 PetscCall(MatRestoreColumnIJ_SeqSELL_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL)); 289 } else { 290 PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL)); 291 } 292 PetscCall(PetscFree(rowhit)); 293 PetscCall(PetscFree(valaddrhit)); 294 PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa)); 295 296 PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->rmap->n, PETSC_DETERMINE, 0, NULL, &c->vscale)); 297 PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols)); 298 PetscFunctionReturn(0); 299 } 300