xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
2525d23c0SHong Zhang #include <../src/mat/impls/baij/seq/baij.h>
3d4002b98SHong Zhang #include <../src/mat/impls/sell/seq/sell.h>
4af0996ceSBarry Smith #include <petsc/private/isimpl.h>
5525d23c0SHong Zhang 
6040ebd07SHong Zhang /*
7040ebd07SHong Zhang     This routine is shared by SeqAIJ and SeqBAIJ matrices,
8040ebd07SHong Zhang     since it operators only on the nonzero structure of the elements or blocks.
9040ebd07SHong Zhang */
10*9371c9d4SSatish Balay PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) {
11a8971b87SHong Zhang   PetscInt  bs, nis = iscoloring->n, m = mat->rmap->n;
12d4002b98SHong Zhang   PetscBool isBAIJ, isSELL;
1393dfae19SHong Zhang 
1493dfae19SHong Zhang   PetscFunctionBegin;
15531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */
169566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
179566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQBAIJ, &isBAIJ));
189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSELL, &isSELL));
19a8971b87SHong Zhang   if (isBAIJ) {
20a8971b87SHong Zhang     c->brows = m;
21a8971b87SHong Zhang     c->bcols = 1;
22531e53bdSHong Zhang   } else { /* seqaij matrix */
23531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */
24531e53bdSHong Zhang     PetscReal mem;
25a8971b87SHong Zhang     PetscInt  nz, brows, bcols;
26d4002b98SHong Zhang     if (isSELL) {
27d4002b98SHong Zhang       Mat_SeqSELL *spA = (Mat_SeqSELL *)mat->data;
28f4b713efSHong Zhang       nz               = spA->nz;
29f4b713efSHong Zhang     } else {
30f4b713efSHong Zhang       Mat_SeqAIJ *spA = (Mat_SeqAIJ *)mat->data;
31f4b713efSHong Zhang       nz              = spA->nz;
32f4b713efSHong Zhang     }
33531e53bdSHong Zhang 
3493dfae19SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
35531e53bdSHong Zhang     mem   = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt);
36531e53bdSHong Zhang     bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar)));
37531e53bdSHong Zhang     brows = 1000 / bcols;
38531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
39531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
40531e53bdSHong Zhang     c->brows = brows;
41531e53bdSHong Zhang     c->bcols = bcols;
4293dfae19SHong Zhang   }
43531e53bdSHong Zhang 
4493dfae19SHong Zhang   c->M       = mat->rmap->N / bs; /* set total rows, columns and local rows */
4593dfae19SHong Zhang   c->N       = mat->cmap->N / bs;
4693dfae19SHong Zhang   c->m       = mat->rmap->N / bs;
4793dfae19SHong Zhang   c->rstart  = 0;
4893dfae19SHong Zhang   c->ncolors = nis;
4991e7fa0fSBarry Smith   c->ctype   = iscoloring->ctype;
5093dfae19SHong Zhang   PetscFunctionReturn(0);
5193dfae19SHong Zhang }
5293dfae19SHong Zhang 
53b3e1f37bSHong Zhang /*
540df34763SHong Zhang  Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into Jacobian for improved cache performance
55b3e1f37bSHong Zhang    Input Parameters:
56b3e1f37bSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
57b3e1f37bSHong Zhang .  color - the coloring context
58b3e1f37bSHong Zhang -  nz - number of local non-zeros in mat
59b3e1f37bSHong Zhang */
60*9371c9d4SSatish Balay PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat, MatFDColoring c, PetscInt nz) {
61b3e1f37bSHong Zhang   PetscInt  i, j, nrows, nbcols, brows = c->brows, bcols = c->bcols, mbs = c->m, nis = c->ncolors;
62b3e1f37bSHong Zhang   PetscInt *color_start, *row_start, *nrows_new, nz_new, row_end;
6379c1e64dSHong Zhang 
6479c1e64dSHong Zhang   PetscFunctionBegin;
65a8971b87SHong Zhang   if (brows < 1 || brows > mbs) brows = mbs;
669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bcols + 1, &color_start, bcols, &row_start));
679566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nis, &nrows_new));
689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bcols * mat->rmap->n, &c->dy));
699566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)c, bcols * mat->rmap->n * sizeof(PetscScalar)));
706a509798SHong Zhang 
71e2c857f8SHong Zhang   nz_new             = 0;
72b3e1f37bSHong Zhang   nbcols             = 0;
73054951adSHong Zhang   color_start[bcols] = 0;
740df34763SHong Zhang 
750df34763SHong Zhang   if (c->htype[0] == 'd') { /* ----  c->htype == 'ds', use MatEntry --------*/
760df34763SHong Zhang     MatEntry *Jentry_new, *Jentry = c->matentry;
77e2cf4d64SStefano Zampini 
789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry_new));
79e2c857f8SHong Zhang     for (i = 0; i < nis; i += bcols) { /* loop over colors */
80054951adSHong Zhang       if (i + bcols > nis) {
81054951adSHong Zhang         color_start[nis - i] = color_start[bcols];
82054951adSHong Zhang         bcols                = nis - i;
83054951adSHong Zhang       }
84054951adSHong Zhang 
85054951adSHong Zhang       color_start[0] = color_start[bcols];
86054951adSHong Zhang       for (j = 0; j < bcols; j++) {
87054951adSHong Zhang         color_start[j + 1] = c->nrows[i + j] + color_start[j];
88054951adSHong Zhang         row_start[j]       = 0;
89054951adSHong Zhang       }
90e2c857f8SHong Zhang 
91e2c857f8SHong Zhang       row_end = brows;
92c8a9c622SHong Zhang       if (row_end > mbs) row_end = mbs;
93054951adSHong Zhang 
94c8a9c622SHong Zhang       while (row_end <= mbs) {        /* loop over block rows */
95e2c857f8SHong Zhang         for (j = 0; j < bcols; j++) { /* loop over block columns */
96e2c857f8SHong Zhang           nrows = c->nrows[i + j];
97054951adSHong Zhang           nz    = color_start[j];
98c8a9c622SHong Zhang           while (row_start[j] < nrows) {
99e2c857f8SHong Zhang             if (Jentry[nz].row >= row_end) {
100054951adSHong Zhang               color_start[j] = nz;
101e2c857f8SHong Zhang               break;
102c8a9c622SHong Zhang             } else {                                                 /* copy Jentry[nz] to Jentry_new[nz_new] */
103c8a9c622SHong Zhang               Jentry_new[nz_new].row     = Jentry[nz].row + j * mbs; /* index in dy-array */
104d880da65SHong Zhang               Jentry_new[nz_new].col     = Jentry[nz].col;
105e2c857f8SHong Zhang               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
106*9371c9d4SSatish Balay               nz_new++;
107*9371c9d4SSatish Balay               nz++;
108*9371c9d4SSatish Balay               row_start[j]++;
109e2c857f8SHong Zhang             }
110e2c857f8SHong Zhang           }
111e2c857f8SHong Zhang         }
112c8a9c622SHong Zhang         if (row_end == mbs) break;
113e2c857f8SHong Zhang         row_end += brows;
114c8a9c622SHong Zhang         if (row_end > mbs) row_end = mbs;
115e2c857f8SHong Zhang       }
1166a509798SHong Zhang       nrows_new[nbcols++] = nz_new;
117e2c857f8SHong Zhang     }
1189566063dSJacob Faibussowitsch     PetscCall(PetscFree(Jentry));
1190df34763SHong Zhang     c->matentry = Jentry_new;
1200df34763SHong Zhang   } else { /* ---------  c->htype == 'wp', use MatEntry2 ------------------*/
1210df34763SHong Zhang     MatEntry2 *Jentry2_new, *Jentry2 = c->matentry2;
122e2cf4d64SStefano Zampini 
1239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry2_new));
1240df34763SHong Zhang     for (i = 0; i < nis; i += bcols) { /* loop over colors */
1250df34763SHong Zhang       if (i + bcols > nis) {
1260df34763SHong Zhang         color_start[nis - i] = color_start[bcols];
1270df34763SHong Zhang         bcols                = nis - i;
1280df34763SHong Zhang       }
1290df34763SHong Zhang 
1300df34763SHong Zhang       color_start[0] = color_start[bcols];
1310df34763SHong Zhang       for (j = 0; j < bcols; j++) {
1320df34763SHong Zhang         color_start[j + 1] = c->nrows[i + j] + color_start[j];
1330df34763SHong Zhang         row_start[j]       = 0;
1340df34763SHong Zhang       }
1350df34763SHong Zhang 
1360df34763SHong Zhang       row_end = brows;
1370df34763SHong Zhang       if (row_end > mbs) row_end = mbs;
1380df34763SHong Zhang 
1390df34763SHong Zhang       while (row_end <= mbs) {        /* loop over block rows */
1400df34763SHong Zhang         for (j = 0; j < bcols; j++) { /* loop over block columns */
1410df34763SHong Zhang           nrows = c->nrows[i + j];
1420df34763SHong Zhang           nz    = color_start[j];
1430df34763SHong Zhang           while (row_start[j] < nrows) {
1440df34763SHong Zhang             if (Jentry2[nz].row >= row_end) {
1450df34763SHong Zhang               color_start[j] = nz;
1460df34763SHong Zhang               break;
1470df34763SHong Zhang             } else {                                                   /* copy Jentry2[nz] to Jentry2_new[nz_new] */
1480df34763SHong Zhang               Jentry2_new[nz_new].row     = Jentry2[nz].row + j * mbs; /* index in dy-array */
1490df34763SHong Zhang               Jentry2_new[nz_new].valaddr = Jentry2[nz].valaddr;
150*9371c9d4SSatish Balay               nz_new++;
151*9371c9d4SSatish Balay               nz++;
152*9371c9d4SSatish Balay               row_start[j]++;
1530df34763SHong Zhang             }
1540df34763SHong Zhang           }
1550df34763SHong Zhang         }
1560df34763SHong Zhang         if (row_end == mbs) break;
1570df34763SHong Zhang         row_end += brows;
1580df34763SHong Zhang         if (row_end > mbs) row_end = mbs;
1590df34763SHong Zhang       }
1600df34763SHong Zhang       nrows_new[nbcols++] = nz_new;
1610df34763SHong Zhang     }
1629566063dSJacob Faibussowitsch     PetscCall(PetscFree(Jentry2));
1630df34763SHong Zhang     c->matentry2 = Jentry2_new;
1640df34763SHong Zhang   } /* ---------------------------------------------*/
1650df34763SHong Zhang 
1669566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color_start, row_start));
167b3e1f37bSHong Zhang 
1686a509798SHong Zhang   for (i = nbcols - 1; i > 0; i--) nrows_new[i] -= nrows_new[i - 1];
1699566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->nrows));
170b3e1f37bSHong Zhang   c->nrows = nrows_new;
171a8971b87SHong Zhang   PetscFunctionReturn(0);
172e2c857f8SHong Zhang }
173a8971b87SHong Zhang 
174*9371c9d4SSatish Balay PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c) {
175071fcb05SBarry Smith   PetscInt           i, n, nrows, mbs = c->m, j, k, m, ncols, col, nis = iscoloring->n, *rowhit, bs, bs2, *spidx, nz, tmp;
176a8971b87SHong Zhang   const PetscInt    *is, *row, *ci, *cj;
177d4002b98SHong Zhang   PetscBool          isBAIJ, isSELL;
178071fcb05SBarry Smith   const PetscScalar *A_val;
179071fcb05SBarry Smith   PetscScalar      **valaddrhit;
180a8971b87SHong Zhang   MatEntry          *Jentry;
1810df34763SHong Zhang   MatEntry2         *Jentry2;
182a8971b87SHong Zhang 
183a8971b87SHong Zhang   PetscFunctionBegin;
1849566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa));
185a8971b87SHong Zhang 
1869566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
1879566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQBAIJ, &isBAIJ));
1889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSELL, &isSELL));
189a8971b87SHong Zhang   if (isBAIJ) {
190a8971b87SHong Zhang     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ *)mat->data;
191e2cf4d64SStefano Zampini 
192a8971b87SHong Zhang     A_val = spA->a;
193a8971b87SHong Zhang     nz    = spA->nz;
194d4002b98SHong Zhang   } else if (isSELL) {
195d4002b98SHong Zhang     Mat_SeqSELL *spA = (Mat_SeqSELL *)mat->data;
196e2cf4d64SStefano Zampini 
197f4b713efSHong Zhang     A_val = spA->val;
198f4b713efSHong Zhang     nz    = spA->nz;
199d4002b98SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqSELL matrix */
200a8971b87SHong Zhang   } else {
201a8971b87SHong Zhang     Mat_SeqAIJ *spA = (Mat_SeqAIJ *)mat->data;
202e2cf4d64SStefano Zampini 
203a8971b87SHong Zhang     A_val = spA->a;
204a8971b87SHong Zhang     nz    = spA->nz;
205a8971b87SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
206a8971b87SHong Zhang   }
207a8971b87SHong Zhang 
2089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns));
2099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nis, &c->nrows)); /* nrows is freeed separately from ncolumns and columns */
2109566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)c, 3 * nis * sizeof(PetscInt)));
211a8971b87SHong Zhang 
2120df34763SHong Zhang   if (c->htype[0] == 'd') {
2139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry));
2149566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)c, nz * sizeof(MatEntry)));
215a8971b87SHong Zhang     c->matentry = Jentry;
2160df34763SHong Zhang   } else if (c->htype[0] == 'w') {
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry2));
2189566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)c, nz * sizeof(MatEntry2)));
2190df34763SHong Zhang     c->matentry2 = Jentry2;
2200df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported");
221a8971b87SHong Zhang 
222a8971b87SHong Zhang   if (isBAIJ) {
2239566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqBAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
224d4002b98SHong Zhang   } else if (isSELL) {
2259566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqSELL_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
226a8971b87SHong Zhang   } else {
2279566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
228a8971b87SHong Zhang   }
229a8971b87SHong Zhang 
2309566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(c->m, &rowhit));
2319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->m, &valaddrhit));
232a8971b87SHong Zhang 
233a8971b87SHong Zhang   nz = 0;
234a8971b87SHong Zhang   for (i = 0; i < nis; i++) { /* loop over colors */
2359566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(c->isa[i], &n));
2369566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(c->isa[i], &is));
237a8971b87SHong Zhang 
238a8971b87SHong Zhang     c->ncolumns[i] = n;
239071fcb05SBarry Smith     c->columns[i]  = (PetscInt *)is;
240071fcb05SBarry Smith     /* note: we know that c->isa is going to be around as long at the c->columns values */
2419566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(c->isa[i], &is));
242a8971b87SHong Zhang 
243a8971b87SHong Zhang     /* fast, crude version requires O(N*N) work */
244a8971b87SHong Zhang     bs2   = bs * bs;
245a8971b87SHong Zhang     nrows = 0;
246a8971b87SHong Zhang     for (j = 0; j < n; j++) { /* loop over columns */
247a8971b87SHong Zhang       col = is[j];
248071fcb05SBarry Smith       tmp = ci[col];
249071fcb05SBarry Smith       row = cj + tmp;
250071fcb05SBarry Smith       m   = ci[col + 1] - tmp;
251a8971b87SHong Zhang       nrows += m;
252a8971b87SHong Zhang       for (k = 0; k < m; k++) { /* loop over columns marking them in rowhit */
253a8971b87SHong Zhang         rowhit[*row]       = col + 1;
254071fcb05SBarry Smith         valaddrhit[*row++] = (PetscScalar *)&A_val[bs2 * spidx[tmp + k]];
255a8971b87SHong Zhang       }
256a8971b87SHong Zhang     }
257a8971b87SHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
258a8971b87SHong Zhang 
2590df34763SHong Zhang     if (c->htype[0] == 'd') {
260a8971b87SHong Zhang       for (j = 0; j < mbs; j++) { /* loop over rows */
261a8971b87SHong Zhang         if (rowhit[j]) {
262a8971b87SHong Zhang           Jentry[nz].row     = j;             /* local row index */
263a8971b87SHong Zhang           Jentry[nz].col     = rowhit[j] - 1; /* local column index */
264a8971b87SHong Zhang           Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
265a8971b87SHong Zhang           nz++;
266a8971b87SHong Zhang           rowhit[j] = 0.0; /* zero rowhit for reuse */
267a8971b87SHong Zhang         }
268a8971b87SHong Zhang       }
2690df34763SHong Zhang     } else {                      /* c->htype == 'wp' */
2700df34763SHong Zhang       for (j = 0; j < mbs; j++) { /* loop over rows */
2710df34763SHong Zhang         if (rowhit[j]) {
2720df34763SHong Zhang           Jentry2[nz].row     = j;             /* local row index */
2730df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
2740df34763SHong Zhang           nz++;
2750df34763SHong Zhang           rowhit[j] = 0.0; /* zero rowhit for reuse */
2760df34763SHong Zhang         }
2770df34763SHong Zhang       }
2780df34763SHong Zhang     }
279a8971b87SHong Zhang   }
280a8971b87SHong Zhang 
281a8971b87SHong Zhang   if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
2829566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz));
283a8971b87SHong Zhang   }
28485740eacSHong Zhang 
285525d23c0SHong Zhang   if (isBAIJ) {
2869566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
2879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy));
2889566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)c, bs * mat->rmap->n * sizeof(PetscScalar)));
289d4002b98SHong Zhang   } else if (isSELL) {
2909566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqSELL_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
291525d23c0SHong Zhang   } else {
2929566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
293525d23c0SHong Zhang   }
2949566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowhit));
2959566063dSJacob Faibussowitsch   PetscCall(PetscFree(valaddrhit));
2969566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa));
297476e0d0aSHong Zhang 
2989566063dSJacob Faibussowitsch   PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->rmap->n, PETSC_DETERMINE, 0, NULL, &c->vscale));
2999566063dSJacob Faibussowitsch   PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols));
30079c1e64dSHong Zhang   PetscFunctionReturn(0);
30179c1e64dSHong Zhang }
302