xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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 */
10d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c)
11d71ae5a4SJacob Faibussowitsch {
12a8971b87SHong Zhang   PetscInt  bs, nis = iscoloring->n, m = mat->rmap->n;
13d4002b98SHong Zhang   PetscBool isBAIJ, isSELL;
1493dfae19SHong Zhang 
1593dfae19SHong Zhang   PetscFunctionBegin;
16531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */
179566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
189566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQBAIJ, &isBAIJ));
199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSELL, &isSELL));
20a8971b87SHong Zhang   if (isBAIJ) {
21a8971b87SHong Zhang     c->brows = m;
22a8971b87SHong Zhang     c->bcols = 1;
23531e53bdSHong Zhang   } else { /* seqaij matrix */
24531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */
25531e53bdSHong Zhang     PetscReal mem;
26a8971b87SHong Zhang     PetscInt  nz, brows, bcols;
27d4002b98SHong Zhang     if (isSELL) {
28d4002b98SHong Zhang       Mat_SeqSELL *spA = (Mat_SeqSELL *)mat->data;
29f4b713efSHong Zhang       nz               = spA->nz;
30f4b713efSHong Zhang     } else {
31f4b713efSHong Zhang       Mat_SeqAIJ *spA = (Mat_SeqAIJ *)mat->data;
32f4b713efSHong Zhang       nz              = spA->nz;
33f4b713efSHong Zhang     }
34531e53bdSHong Zhang 
3593dfae19SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
36531e53bdSHong Zhang     mem   = nz * (sizeof(PetscScalar) + sizeof(PetscInt)) + 3 * m * sizeof(PetscInt);
37531e53bdSHong Zhang     bcols = (PetscInt)(0.5 * mem / (m * sizeof(PetscScalar)));
38531e53bdSHong Zhang     brows = 1000 / bcols;
39531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
40531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
41531e53bdSHong Zhang     c->brows = brows;
42531e53bdSHong Zhang     c->bcols = bcols;
4393dfae19SHong Zhang   }
44531e53bdSHong Zhang 
4593dfae19SHong Zhang   c->M       = mat->rmap->N / bs; /* set total rows, columns and local rows */
4693dfae19SHong Zhang   c->N       = mat->cmap->N / bs;
4793dfae19SHong Zhang   c->m       = mat->rmap->N / bs;
4893dfae19SHong Zhang   c->rstart  = 0;
4993dfae19SHong Zhang   c->ncolors = nis;
5091e7fa0fSBarry Smith   c->ctype   = iscoloring->ctype;
5193dfae19SHong Zhang   PetscFunctionReturn(0);
5293dfae19SHong Zhang }
5393dfae19SHong Zhang 
54b3e1f37bSHong Zhang /*
550df34763SHong Zhang  Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into Jacobian for improved cache performance
56b3e1f37bSHong Zhang    Input Parameters:
57b3e1f37bSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
58b3e1f37bSHong Zhang .  color - the coloring context
59b3e1f37bSHong Zhang -  nz - number of local non-zeros in mat
60b3e1f37bSHong Zhang */
61d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat, MatFDColoring c, PetscInt nz)
62d71ae5a4SJacob Faibussowitsch {
63b3e1f37bSHong Zhang   PetscInt  i, j, nrows, nbcols, brows = c->brows, bcols = c->bcols, mbs = c->m, nis = c->ncolors;
64b3e1f37bSHong Zhang   PetscInt *color_start, *row_start, *nrows_new, nz_new, row_end;
6579c1e64dSHong Zhang 
6679c1e64dSHong Zhang   PetscFunctionBegin;
67a8971b87SHong Zhang   if (brows < 1 || brows > mbs) brows = mbs;
689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bcols + 1, &color_start, bcols, &row_start));
699566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nis, &nrows_new));
709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bcols * mat->rmap->n, &c->dy));
716a509798SHong Zhang 
72e2c857f8SHong Zhang   nz_new             = 0;
73b3e1f37bSHong Zhang   nbcols             = 0;
74054951adSHong Zhang   color_start[bcols] = 0;
750df34763SHong Zhang 
760df34763SHong Zhang   if (c->htype[0] == 'd') { /* ----  c->htype == 'ds', use MatEntry --------*/
770df34763SHong Zhang     MatEntry *Jentry_new, *Jentry = c->matentry;
78e2cf4d64SStefano Zampini 
799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry_new));
80e2c857f8SHong Zhang     for (i = 0; i < nis; i += bcols) { /* loop over colors */
81054951adSHong Zhang       if (i + bcols > nis) {
82054951adSHong Zhang         color_start[nis - i] = color_start[bcols];
83054951adSHong Zhang         bcols                = nis - i;
84054951adSHong Zhang       }
85054951adSHong Zhang 
86054951adSHong Zhang       color_start[0] = color_start[bcols];
87054951adSHong Zhang       for (j = 0; j < bcols; j++) {
88054951adSHong Zhang         color_start[j + 1] = c->nrows[i + j] + color_start[j];
89054951adSHong Zhang         row_start[j]       = 0;
90054951adSHong Zhang       }
91e2c857f8SHong Zhang 
92e2c857f8SHong Zhang       row_end = brows;
93c8a9c622SHong Zhang       if (row_end > mbs) row_end = mbs;
94054951adSHong Zhang 
95c8a9c622SHong Zhang       while (row_end <= mbs) {        /* loop over block rows */
96e2c857f8SHong Zhang         for (j = 0; j < bcols; j++) { /* loop over block columns */
97e2c857f8SHong Zhang           nrows = c->nrows[i + j];
98054951adSHong Zhang           nz    = color_start[j];
99c8a9c622SHong Zhang           while (row_start[j] < nrows) {
100e2c857f8SHong Zhang             if (Jentry[nz].row >= row_end) {
101054951adSHong Zhang               color_start[j] = nz;
102e2c857f8SHong Zhang               break;
103c8a9c622SHong Zhang             } else {                                                 /* copy Jentry[nz] to Jentry_new[nz_new] */
104c8a9c622SHong Zhang               Jentry_new[nz_new].row     = Jentry[nz].row + j * mbs; /* index in dy-array */
105d880da65SHong Zhang               Jentry_new[nz_new].col     = Jentry[nz].col;
106e2c857f8SHong Zhang               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
1079371c9d4SSatish Balay               nz_new++;
1089371c9d4SSatish Balay               nz++;
1099371c9d4SSatish Balay               row_start[j]++;
110e2c857f8SHong Zhang             }
111e2c857f8SHong Zhang           }
112e2c857f8SHong Zhang         }
113c8a9c622SHong Zhang         if (row_end == mbs) break;
114e2c857f8SHong Zhang         row_end += brows;
115c8a9c622SHong Zhang         if (row_end > mbs) row_end = mbs;
116e2c857f8SHong Zhang       }
1176a509798SHong Zhang       nrows_new[nbcols++] = nz_new;
118e2c857f8SHong Zhang     }
1199566063dSJacob Faibussowitsch     PetscCall(PetscFree(Jentry));
1200df34763SHong Zhang     c->matentry = Jentry_new;
1210df34763SHong Zhang   } else { /* ---------  c->htype == 'wp', use MatEntry2 ------------------*/
1220df34763SHong Zhang     MatEntry2 *Jentry2_new, *Jentry2 = c->matentry2;
123e2cf4d64SStefano Zampini 
1249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry2_new));
1250df34763SHong Zhang     for (i = 0; i < nis; i += bcols) { /* loop over colors */
1260df34763SHong Zhang       if (i + bcols > nis) {
1270df34763SHong Zhang         color_start[nis - i] = color_start[bcols];
1280df34763SHong Zhang         bcols                = nis - i;
1290df34763SHong Zhang       }
1300df34763SHong Zhang 
1310df34763SHong Zhang       color_start[0] = color_start[bcols];
1320df34763SHong Zhang       for (j = 0; j < bcols; j++) {
1330df34763SHong Zhang         color_start[j + 1] = c->nrows[i + j] + color_start[j];
1340df34763SHong Zhang         row_start[j]       = 0;
1350df34763SHong Zhang       }
1360df34763SHong Zhang 
1370df34763SHong Zhang       row_end = brows;
1380df34763SHong Zhang       if (row_end > mbs) row_end = mbs;
1390df34763SHong Zhang 
1400df34763SHong Zhang       while (row_end <= mbs) {        /* loop over block rows */
1410df34763SHong Zhang         for (j = 0; j < bcols; j++) { /* loop over block columns */
1420df34763SHong Zhang           nrows = c->nrows[i + j];
1430df34763SHong Zhang           nz    = color_start[j];
1440df34763SHong Zhang           while (row_start[j] < nrows) {
1450df34763SHong Zhang             if (Jentry2[nz].row >= row_end) {
1460df34763SHong Zhang               color_start[j] = nz;
1470df34763SHong Zhang               break;
1480df34763SHong Zhang             } else {                                                   /* copy Jentry2[nz] to Jentry2_new[nz_new] */
1490df34763SHong Zhang               Jentry2_new[nz_new].row     = Jentry2[nz].row + j * mbs; /* index in dy-array */
1500df34763SHong Zhang               Jentry2_new[nz_new].valaddr = Jentry2[nz].valaddr;
1519371c9d4SSatish Balay               nz_new++;
1529371c9d4SSatish Balay               nz++;
1539371c9d4SSatish Balay               row_start[j]++;
1540df34763SHong Zhang             }
1550df34763SHong Zhang           }
1560df34763SHong Zhang         }
1570df34763SHong Zhang         if (row_end == mbs) break;
1580df34763SHong Zhang         row_end += brows;
1590df34763SHong Zhang         if (row_end > mbs) row_end = mbs;
1600df34763SHong Zhang       }
1610df34763SHong Zhang       nrows_new[nbcols++] = nz_new;
1620df34763SHong Zhang     }
1639566063dSJacob Faibussowitsch     PetscCall(PetscFree(Jentry2));
1640df34763SHong Zhang     c->matentry2 = Jentry2_new;
1650df34763SHong Zhang   } /* ---------------------------------------------*/
1660df34763SHong Zhang 
1679566063dSJacob Faibussowitsch   PetscCall(PetscFree2(color_start, row_start));
168b3e1f37bSHong Zhang 
1696a509798SHong Zhang   for (i = nbcols - 1; i > 0; i--) nrows_new[i] -= nrows_new[i - 1];
1709566063dSJacob Faibussowitsch   PetscCall(PetscFree(c->nrows));
171b3e1f37bSHong Zhang   c->nrows = nrows_new;
172a8971b87SHong Zhang   PetscFunctionReturn(0);
173e2c857f8SHong Zhang }
174a8971b87SHong Zhang 
175d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat, ISColoring iscoloring, MatFDColoring c)
176d71ae5a4SJacob Faibussowitsch {
177071fcb05SBarry Smith   PetscInt           i, n, nrows, mbs = c->m, j, k, m, ncols, col, nis = iscoloring->n, *rowhit, bs, bs2, *spidx, nz, tmp;
178a8971b87SHong Zhang   const PetscInt    *is, *row, *ci, *cj;
179d4002b98SHong Zhang   PetscBool          isBAIJ, isSELL;
180071fcb05SBarry Smith   const PetscScalar *A_val;
181071fcb05SBarry Smith   PetscScalar      **valaddrhit;
182a8971b87SHong Zhang   MatEntry          *Jentry;
1830df34763SHong Zhang   MatEntry2         *Jentry2;
184a8971b87SHong Zhang 
185a8971b87SHong Zhang   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(iscoloring, PETSC_OWN_POINTER, PETSC_IGNORE, &c->isa));
187a8971b87SHong Zhang 
1889566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &bs));
1899566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)mat, MATSEQBAIJ, &isBAIJ));
1909566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATSEQSELL, &isSELL));
191a8971b87SHong Zhang   if (isBAIJ) {
192a8971b87SHong Zhang     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ *)mat->data;
193e2cf4d64SStefano Zampini 
194a8971b87SHong Zhang     A_val = spA->a;
195a8971b87SHong Zhang     nz    = spA->nz;
196d4002b98SHong Zhang   } else if (isSELL) {
197d4002b98SHong Zhang     Mat_SeqSELL *spA = (Mat_SeqSELL *)mat->data;
198e2cf4d64SStefano Zampini 
199f4b713efSHong Zhang     A_val = spA->val;
200f4b713efSHong Zhang     nz    = spA->nz;
201d4002b98SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqSELL matrix */
202a8971b87SHong Zhang   } else {
203a8971b87SHong Zhang     Mat_SeqAIJ *spA = (Mat_SeqAIJ *)mat->data;
204e2cf4d64SStefano Zampini 
205a8971b87SHong Zhang     A_val = spA->a;
206a8971b87SHong Zhang     nz    = spA->nz;
207a8971b87SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
208a8971b87SHong Zhang   }
209a8971b87SHong Zhang 
2109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nis, &c->ncolumns, nis, &c->columns));
211*d5b43468SJose E. Roman   PetscCall(PetscMalloc1(nis, &c->nrows)); /* nrows is freed separately from ncolumns and columns */
212a8971b87SHong Zhang 
2130df34763SHong Zhang   if (c->htype[0] == 'd') {
2149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry));
215a8971b87SHong Zhang     c->matentry = Jentry;
2160df34763SHong Zhang   } else if (c->htype[0] == 'w') {
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nz, &Jentry2));
2180df34763SHong Zhang     c->matentry2 = Jentry2;
2190df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "htype is not supported");
220a8971b87SHong Zhang 
221a8971b87SHong Zhang   if (isBAIJ) {
2229566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqBAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
223d4002b98SHong Zhang   } else if (isSELL) {
2249566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqSELL_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
225a8971b87SHong Zhang   } else {
2269566063dSJacob Faibussowitsch     PetscCall(MatGetColumnIJ_SeqAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
227a8971b87SHong Zhang   }
228a8971b87SHong Zhang 
2299566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(c->m, &rowhit));
2309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(c->m, &valaddrhit));
231a8971b87SHong Zhang 
232a8971b87SHong Zhang   nz = 0;
233a8971b87SHong Zhang   for (i = 0; i < nis; i++) { /* loop over colors */
2349566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(c->isa[i], &n));
2359566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(c->isa[i], &is));
236a8971b87SHong Zhang 
237a8971b87SHong Zhang     c->ncolumns[i] = n;
238071fcb05SBarry Smith     c->columns[i]  = (PetscInt *)is;
239071fcb05SBarry Smith     /* note: we know that c->isa is going to be around as long at the c->columns values */
2409566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(c->isa[i], &is));
241a8971b87SHong Zhang 
242a8971b87SHong Zhang     /* fast, crude version requires O(N*N) work */
243a8971b87SHong Zhang     bs2   = bs * bs;
244a8971b87SHong Zhang     nrows = 0;
245a8971b87SHong Zhang     for (j = 0; j < n; j++) { /* loop over columns */
246a8971b87SHong Zhang       col = is[j];
247071fcb05SBarry Smith       tmp = ci[col];
248071fcb05SBarry Smith       row = cj + tmp;
249071fcb05SBarry Smith       m   = ci[col + 1] - tmp;
250a8971b87SHong Zhang       nrows += m;
251a8971b87SHong Zhang       for (k = 0; k < m; k++) { /* loop over columns marking them in rowhit */
252a8971b87SHong Zhang         rowhit[*row]       = col + 1;
253071fcb05SBarry Smith         valaddrhit[*row++] = (PetscScalar *)&A_val[bs2 * spidx[tmp + k]];
254a8971b87SHong Zhang       }
255a8971b87SHong Zhang     }
256a8971b87SHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
257a8971b87SHong Zhang 
2580df34763SHong Zhang     if (c->htype[0] == 'd') {
259a8971b87SHong Zhang       for (j = 0; j < mbs; j++) { /* loop over rows */
260a8971b87SHong Zhang         if (rowhit[j]) {
261a8971b87SHong Zhang           Jentry[nz].row     = j;             /* local row index */
262a8971b87SHong Zhang           Jentry[nz].col     = rowhit[j] - 1; /* local column index */
263a8971b87SHong Zhang           Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
264a8971b87SHong Zhang           nz++;
265a8971b87SHong Zhang           rowhit[j] = 0.0; /* zero rowhit for reuse */
266a8971b87SHong Zhang         }
267a8971b87SHong Zhang       }
2680df34763SHong Zhang     } else {                      /* c->htype == 'wp' */
2690df34763SHong Zhang       for (j = 0; j < mbs; j++) { /* loop over rows */
2700df34763SHong Zhang         if (rowhit[j]) {
2710df34763SHong Zhang           Jentry2[nz].row     = j;             /* local row index */
2720df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
2730df34763SHong Zhang           nz++;
2740df34763SHong Zhang           rowhit[j] = 0.0; /* zero rowhit for reuse */
2750df34763SHong Zhang         }
2760df34763SHong Zhang       }
2770df34763SHong Zhang     }
278a8971b87SHong Zhang   }
279a8971b87SHong Zhang 
280a8971b87SHong Zhang   if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
2819566063dSJacob Faibussowitsch     PetscCall(MatFDColoringSetUpBlocked_AIJ_Private(mat, c, nz));
282a8971b87SHong Zhang   }
28385740eacSHong Zhang 
284525d23c0SHong Zhang   if (isBAIJ) {
2859566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqBAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
2869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * mat->rmap->n, &c->dy));
287d4002b98SHong Zhang   } else if (isSELL) {
2889566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqSELL_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
289525d23c0SHong Zhang   } else {
2909566063dSJacob Faibussowitsch     PetscCall(MatRestoreColumnIJ_SeqAIJ_Color(mat, 0, PETSC_FALSE, PETSC_FALSE, &ncols, &ci, &cj, &spidx, NULL));
291525d23c0SHong Zhang   }
2929566063dSJacob Faibussowitsch   PetscCall(PetscFree(rowhit));
2939566063dSJacob Faibussowitsch   PetscCall(PetscFree(valaddrhit));
2949566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(iscoloring, PETSC_OWN_POINTER, &c->isa));
295476e0d0aSHong Zhang 
2969566063dSJacob Faibussowitsch   PetscCall(VecCreateGhost(PetscObjectComm((PetscObject)mat), mat->rmap->n, PETSC_DETERMINE, 0, NULL, &c->vscale));
2979566063dSJacob Faibussowitsch   PetscCall(PetscInfo(c, "ncolors %" PetscInt_FMT ", brows %" PetscInt_FMT " and bcols %" PetscInt_FMT " are used.\n", c->ncolors, c->brows, c->bcols));
29879c1e64dSHong Zhang   PetscFunctionReturn(0);
29979c1e64dSHong Zhang }
300