xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 4cc2b5b5fe4ffd09e5956b56d7cdc4f43e324103)
1421e10b8SBarry Smith /*
2421e10b8SBarry Smith    This provides a matrix that consists of Mats
3421e10b8SBarry Smith */
4421e10b8SBarry Smith 
5af0996ceSBarry Smith #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
6c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> /* use the common AIJ data-structure */
7421e10b8SBarry Smith 
8421e10b8SBarry Smith typedef struct {
9421e10b8SBarry Smith   SEQAIJHEADER(Mat);
10421e10b8SBarry Smith   SEQBAIJHEADER;
11ccb205f8SBarry Smith   Mat *diags;
12421e10b8SBarry Smith 
136e63c7a1SBarry Smith   Vec left, right, middle, workb; /* dummy vectors to perform local parts of product */
14421e10b8SBarry Smith } Mat_BlockMat;
15421e10b8SBarry Smith 
16d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
17d71ae5a4SJacob Faibussowitsch {
18290bbb0aSBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat *)A->data;
19290bbb0aSBarry Smith   PetscScalar       *x;
20a2ea699eSBarry Smith   const Mat         *v;
21290bbb0aSBarry Smith   const PetscScalar *b;
22d0f46423SBarry Smith   PetscInt           n = A->cmap->n, i, mbs = n / A->rmap->bs, j, bs = A->rmap->bs;
23290bbb0aSBarry Smith   const PetscInt    *idx;
24290bbb0aSBarry Smith   IS                 row, col;
25290bbb0aSBarry Smith   MatFactorInfo      info;
266e63c7a1SBarry Smith   Vec                left = a->left, right = a->right, middle = a->middle;
27290bbb0aSBarry Smith   Mat               *diag;
28290bbb0aSBarry Smith 
29290bbb0aSBarry Smith   PetscFunctionBegin;
30290bbb0aSBarry Smith   its = its * lits;
3108401ef6SPierre Jolivet   PetscCheck(its > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive", its, lits);
32aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
3308401ef6SPierre Jolivet   PetscCheck(omega == 1.0, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for omega not equal to 1.0");
3428b400f6SJacob Faibussowitsch   PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for fshift");
357a46b595SBarry Smith   PetscCheck(!((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot do backward sweep without forward sweep");
36290bbb0aSBarry Smith 
37290bbb0aSBarry Smith   if (!a->diags) {
389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs, &a->diags));
399566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
40290bbb0aSBarry Smith     for (i = 0; i < mbs; i++) {
419566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND, &row, &col));
429566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(a->diags[i], a->a[a->diag[i]], row, &info));
439566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(a->diags[i], a->a[a->diag[i]], &info));
449566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&row));
459566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&col));
46290bbb0aSBarry Smith     }
479566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(bb, &a->workb));
48290bbb0aSBarry Smith   }
49290bbb0aSBarry Smith   diag = a->diags;
50290bbb0aSBarry Smith 
519566063dSJacob Faibussowitsch   PetscCall(VecSet(xx, 0.0));
529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
53dd8e379bSPierre Jolivet   /* copy right-hand side because it must be modified during iteration */
549566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, a->workb));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(a->workb, &b));
56290bbb0aSBarry Smith 
5741f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
58290bbb0aSBarry Smith   while (its--) {
59290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
60290bbb0aSBarry Smith       for (i = 0; i < mbs; i++) {
616e63c7a1SBarry Smith         n   = a->i[i + 1] - a->i[i] - 1;
626e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
636e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
64290bbb0aSBarry Smith 
659566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
66290bbb0aSBarry Smith         for (j = 0; j < n; j++) {
679566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(right, x + idx[j] * bs));
689566063dSJacob Faibussowitsch           PetscCall(MatMultAdd(v[j], right, left, left));
699566063dSJacob Faibussowitsch           PetscCall(VecResetArray(right));
70290bbb0aSBarry Smith         }
719566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
729566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
739566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
74290bbb0aSBarry Smith 
759566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
769566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
776e63c7a1SBarry Smith 
78dd8e379bSPierre Jolivet         /* now adjust right-hand side, see MatSOR_SeqSBAIJ */
796e63c7a1SBarry Smith         for (j = 0; j < n; j++) {
809566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(v[j], right, left));
819566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(middle, b + idx[j] * bs));
829566063dSJacob Faibussowitsch           PetscCall(VecAXPY(middle, -1.0, left));
839566063dSJacob Faibussowitsch           PetscCall(VecResetArray(middle));
846e63c7a1SBarry Smith         }
859566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
86290bbb0aSBarry Smith       }
87290bbb0aSBarry Smith     }
88290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
89290bbb0aSBarry Smith       for (i = mbs - 1; i >= 0; i--) {
906e63c7a1SBarry Smith         n   = a->i[i + 1] - a->i[i] - 1;
916e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
926e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
93290bbb0aSBarry Smith 
949566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
95290bbb0aSBarry Smith         for (j = 0; j < n; j++) {
969566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(right, x + idx[j] * bs));
979566063dSJacob Faibussowitsch           PetscCall(MatMultAdd(v[j], right, left, left));
989566063dSJacob Faibussowitsch           PetscCall(VecResetArray(right));
99290bbb0aSBarry Smith         }
1009566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
1019566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
1029566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
103290bbb0aSBarry Smith 
1049566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
1059566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
1069566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
107290bbb0aSBarry Smith       }
108290bbb0aSBarry Smith     }
109290bbb0aSBarry Smith   }
1109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(a->workb, &b));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113290bbb0aSBarry Smith }
114290bbb0aSBarry Smith 
115d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSOR_BlockMat(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
116d71ae5a4SJacob Faibussowitsch {
117ccb205f8SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat *)A->data;
118ccb205f8SBarry Smith   PetscScalar       *x;
119a2ea699eSBarry Smith   const Mat         *v;
120ccb205f8SBarry Smith   const PetscScalar *b;
121d0f46423SBarry Smith   PetscInt           n = A->cmap->n, i, mbs = n / A->rmap->bs, j, bs = A->rmap->bs;
122ccb205f8SBarry Smith   const PetscInt    *idx;
123ccb205f8SBarry Smith   IS                 row, col;
124ccb205f8SBarry Smith   MatFactorInfo      info;
125ccb205f8SBarry Smith   Vec                left = a->left, right = a->right;
126ccb205f8SBarry Smith   Mat               *diag;
127ccb205f8SBarry Smith 
128ccb205f8SBarry Smith   PetscFunctionBegin;
129ccb205f8SBarry Smith   its = its * lits;
13008401ef6SPierre Jolivet   PetscCheck(its > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive", its, lits);
131aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
13208401ef6SPierre Jolivet   PetscCheck(omega == 1.0, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for omega not equal to 1.0");
13328b400f6SJacob Faibussowitsch   PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for fshift");
134ccb205f8SBarry Smith 
135ccb205f8SBarry Smith   if (!a->diags) {
1369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs, &a->diags));
1379566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
138ccb205f8SBarry Smith     for (i = 0; i < mbs; i++) {
1399566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND, &row, &col));
1409566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(a->diags[i], a->a[a->diag[i]], row, col, &info));
1419566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(a->diags[i], a->a[a->diag[i]], &info));
1429566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&row));
1439566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&col));
144ccb205f8SBarry Smith     }
145ccb205f8SBarry Smith   }
146ccb205f8SBarry Smith   diag = a->diags;
147ccb205f8SBarry Smith 
1489566063dSJacob Faibussowitsch   PetscCall(VecSet(xx, 0.0));
1499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
151ccb205f8SBarry Smith 
15241f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
153ccb205f8SBarry Smith   while (its--) {
154ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
155ccb205f8SBarry Smith       for (i = 0; i < mbs; i++) {
156ccb205f8SBarry Smith         n   = a->i[i + 1] - a->i[i];
157ccb205f8SBarry Smith         idx = a->j + a->i[i];
158ccb205f8SBarry Smith         v   = a->a + a->i[i];
159ccb205f8SBarry Smith 
1609566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
161ccb205f8SBarry Smith         for (j = 0; j < n; j++) {
162ccb205f8SBarry Smith           if (idx[j] != i) {
1639566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(right, x + idx[j] * bs));
1649566063dSJacob Faibussowitsch             PetscCall(MatMultAdd(v[j], right, left, left));
1659566063dSJacob Faibussowitsch             PetscCall(VecResetArray(right));
166ccb205f8SBarry Smith           }
167ccb205f8SBarry Smith         }
1689566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
1699566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
1709566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
171ccb205f8SBarry Smith 
1729566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
1739566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
1749566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
175ccb205f8SBarry Smith       }
176ccb205f8SBarry Smith     }
177ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
178ccb205f8SBarry Smith       for (i = mbs - 1; i >= 0; i--) {
179ccb205f8SBarry Smith         n   = a->i[i + 1] - a->i[i];
180ccb205f8SBarry Smith         idx = a->j + a->i[i];
181ccb205f8SBarry Smith         v   = a->a + a->i[i];
182ccb205f8SBarry Smith 
1839566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
184ccb205f8SBarry Smith         for (j = 0; j < n; j++) {
185ccb205f8SBarry Smith           if (idx[j] != i) {
1869566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(right, x + idx[j] * bs));
1879566063dSJacob Faibussowitsch             PetscCall(MatMultAdd(v[j], right, left, left));
1889566063dSJacob Faibussowitsch             PetscCall(VecResetArray(right));
189ccb205f8SBarry Smith           }
190ccb205f8SBarry Smith         }
1919566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
1929566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
1939566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
194ccb205f8SBarry Smith 
1959566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
1969566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
1979566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
198ccb205f8SBarry Smith       }
199ccb205f8SBarry Smith     }
200ccb205f8SBarry Smith   }
2019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204ccb205f8SBarry Smith }
205ccb205f8SBarry Smith 
206d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValues_BlockMat(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
207d71ae5a4SJacob Faibussowitsch {
208421e10b8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
209421e10b8SBarry Smith   PetscInt     *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N, lastcol = -1;
210421e10b8SBarry Smith   PetscInt     *imax = a->imax, *ai = a->i, *ailen = a->ilen;
211d0f46423SBarry Smith   PetscInt     *aj = a->j, nonew = a->nonew, bs = A->rmap->bs, brow, bcol;
212421e10b8SBarry Smith   PetscInt      ridx, cidx;
213ace3abfcSBarry Smith   PetscBool     roworiented = a->roworiented;
214421e10b8SBarry Smith   MatScalar     value;
215421e10b8SBarry Smith   Mat          *ap, *aa = a->a;
216421e10b8SBarry Smith 
217421e10b8SBarry Smith   PetscFunctionBegin;
218421e10b8SBarry Smith   for (k = 0; k < m; k++) { /* loop over added rows */
219421e10b8SBarry Smith     row  = im[k];
220421e10b8SBarry Smith     brow = row / bs;
221421e10b8SBarry Smith     if (row < 0) continue;
2226bdcaf15SBarry Smith     PetscCheck(row < A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->N - 1);
223421e10b8SBarry Smith     rp   = aj + ai[brow];
224421e10b8SBarry Smith     ap   = aa + ai[brow];
225421e10b8SBarry Smith     rmax = imax[brow];
226421e10b8SBarry Smith     nrow = ailen[brow];
227421e10b8SBarry Smith     low  = 0;
228421e10b8SBarry Smith     high = nrow;
229421e10b8SBarry Smith     for (l = 0; l < n; l++) { /* loop over added columns */
230421e10b8SBarry Smith       if (in[l] < 0) continue;
2316bdcaf15SBarry Smith       PetscCheck(in[l] < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[l], A->cmap->n - 1);
2329371c9d4SSatish Balay       col  = in[l];
2339371c9d4SSatish Balay       bcol = col / bs;
234b94d7dedSBarry Smith       if (A->symmetric == PETSC_BOOL3_TRUE && brow > bcol) continue;
2359371c9d4SSatish Balay       ridx = row % bs;
2369371c9d4SSatish Balay       cidx = col % bs;
2372205254eSKarl Rupp       if (roworiented) value = v[l + k * n];
2382205254eSKarl Rupp       else value = v[k + l * m];
2392205254eSKarl Rupp 
2402205254eSKarl Rupp       if (col <= lastcol) low = 0;
2412205254eSKarl Rupp       else high = nrow;
242421e10b8SBarry Smith       lastcol = col;
243421e10b8SBarry Smith       while (high - low > 7) {
244421e10b8SBarry Smith         t = (low + high) / 2;
245421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
246421e10b8SBarry Smith         else low = t;
247421e10b8SBarry Smith       }
248421e10b8SBarry Smith       for (i = low; i < high; i++) {
249421e10b8SBarry Smith         if (rp[i] > bcol) break;
2502205254eSKarl Rupp         if (rp[i] == bcol) goto noinsert1;
251421e10b8SBarry Smith       }
252421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
25308401ef6SPierre Jolivet       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
254fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A, a->mbs, 1, nrow, brow, bcol, rmax, aa, ai, aj, rp, ap, imax, nonew, Mat);
2559371c9d4SSatish Balay       N = nrow++ - 1;
2569371c9d4SSatish Balay       high++;
257421e10b8SBarry Smith       /* shift up all the later entries in this row */
258421e10b8SBarry Smith       for (ii = N; ii >= i; ii--) {
259421e10b8SBarry Smith         rp[ii + 1] = rp[ii];
260421e10b8SBarry Smith         ap[ii + 1] = ap[ii];
261421e10b8SBarry Smith       }
262f4259b30SLisandro Dalcin       if (N >= i) ap[i] = NULL;
263421e10b8SBarry Smith       rp[i] = bcol;
264421e10b8SBarry Smith       a->nz++;
265421e10b8SBarry Smith     noinsert1:;
26648a46eb9SPierre Jolivet       if (!*(ap + i)) PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, NULL, ap + i));
2679566063dSJacob Faibussowitsch       PetscCall(MatSetValues(ap[i], 1, &ridx, 1, &cidx, &value, is));
268421e10b8SBarry Smith       low = i;
269421e10b8SBarry Smith     }
270421e10b8SBarry Smith     ailen[brow] = nrow;
271421e10b8SBarry Smith   }
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273421e10b8SBarry Smith }
274421e10b8SBarry Smith 
275d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
276d71ae5a4SJacob Faibussowitsch {
277a5973382SShri Abhyankar   Mat                tmpA;
278a5973382SShri Abhyankar   PetscInt           i, j, m, n, bs = 1, ncols, *lens, currentcol, mbs, **ii, *ilens, nextcol, *llens, cnt = 0;
279a5973382SShri Abhyankar   const PetscInt    *cols;
280a5973382SShri Abhyankar   const PetscScalar *values;
281ace3abfcSBarry Smith   PetscBool          flg = PETSC_FALSE, notdone;
282a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
283a5973382SShri Abhyankar   Mat_BlockMat      *amat;
284a5973382SShri Abhyankar 
285a5973382SShri Abhyankar   PetscFunctionBegin;
286c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
2879566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
2889566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &tmpA));
2899566063dSJacob Faibussowitsch   PetscCall(MatSetType(tmpA, MATSEQAIJ));
2909566063dSJacob Faibussowitsch   PetscCall(MatLoad_SeqAIJ(tmpA, viewer));
291a5973382SShri Abhyankar 
2929566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(tmpA, &m, &n));
293d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, NULL, "Options for loading BlockMat matrix 1", "Mat");
2949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-matload_block_size", "Set the blocksize used to store the matrix", "MatLoad", bs, &bs, NULL));
2959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-matload_symmetric", "Store the matrix as symmetric", "MatLoad", flg, &flg, NULL));
296d0609cedSBarry Smith   PetscOptionsEnd();
297a5973382SShri Abhyankar 
298a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
299a5973382SShri Abhyankar   a   = (Mat_SeqAIJ *)tmpA->data;
300a5973382SShri Abhyankar   mbs = m / bs;
3019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &lens, bs, &ii, bs, &ilens));
3029566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lens, mbs));
303a5973382SShri Abhyankar 
304a5973382SShri Abhyankar   for (i = 0; i < mbs; i++) {
305a5973382SShri Abhyankar     for (j = 0; j < bs; j++) {
306a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i * bs + j];
307a5973382SShri Abhyankar       ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j];
308a5973382SShri Abhyankar     }
309a5973382SShri Abhyankar 
310a5973382SShri Abhyankar     currentcol = -1;
311a5973382SShri Abhyankar     while (PETSC_TRUE) {
312a5973382SShri Abhyankar       notdone = PETSC_FALSE;
313a5973382SShri Abhyankar       nextcol = 1000000000;
314a5973382SShri Abhyankar       for (j = 0; j < bs; j++) {
315f4f49eeaSPierre Jolivet         while (ilens[j] > 0 && ii[j][0] / bs <= currentcol) {
316a5973382SShri Abhyankar           ii[j]++;
317a5973382SShri Abhyankar           ilens[j]--;
318a5973382SShri Abhyankar         }
319a5973382SShri Abhyankar         if (ilens[j] > 0) {
320a5973382SShri Abhyankar           notdone = PETSC_TRUE;
321a5973382SShri Abhyankar           nextcol = PetscMin(nextcol, ii[j][0] / bs);
322a5973382SShri Abhyankar         }
323a5973382SShri Abhyankar       }
324a5973382SShri Abhyankar       if (!notdone) break;
325a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
326a5973382SShri Abhyankar       currentcol = nextcol;
327a5973382SShri Abhyankar     }
328a5973382SShri Abhyankar   }
329a5973382SShri Abhyankar 
33048a46eb9SPierre Jolivet   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) PetscCall(MatSetSizes(newmat, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
3319566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(newmat, bs, 0, lens));
3321baa6e33SBarry Smith   if (flg) PetscCall(MatSetOption(newmat, MAT_SYMMETRIC, PETSC_TRUE));
333a5973382SShri Abhyankar   amat = (Mat_BlockMat *)(newmat)->data;
334a5973382SShri Abhyankar 
335a5973382SShri Abhyankar   /* preallocate the submatrices */
3369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &llens));
337a5973382SShri Abhyankar   for (i = 0; i < mbs; i++) { /* loops for block rows */
338a5973382SShri Abhyankar     for (j = 0; j < bs; j++) {
339a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i * bs + j];
340a5973382SShri Abhyankar       ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j];
341a5973382SShri Abhyankar     }
342a5973382SShri Abhyankar 
343a5973382SShri Abhyankar     currentcol = 1000000000;
344a5973382SShri Abhyankar     for (j = 0; j < bs; j++) { /* loop over rows in block finding first nonzero block */
345ad540459SPierre Jolivet       if (ilens[j] > 0) currentcol = PetscMin(currentcol, ii[j][0] / bs);
346a5973382SShri Abhyankar     }
347a5973382SShri Abhyankar 
348a5973382SShri Abhyankar     while (PETSC_TRUE) { /* loops over blocks in block row */
349a5973382SShri Abhyankar       notdone = PETSC_FALSE;
350a5973382SShri Abhyankar       nextcol = 1000000000;
3519566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(llens, bs));
352a5973382SShri Abhyankar       for (j = 0; j < bs; j++) {                              /* loop over rows in block */
353f4f49eeaSPierre Jolivet         while (ilens[j] > 0 && ii[j][0] / bs <= currentcol) { /* loop over columns in row */
354a5973382SShri Abhyankar           ii[j]++;
355a5973382SShri Abhyankar           ilens[j]--;
356a5973382SShri Abhyankar           llens[j]++;
357a5973382SShri Abhyankar         }
358a5973382SShri Abhyankar         if (ilens[j] > 0) {
359a5973382SShri Abhyankar           notdone = PETSC_TRUE;
360a5973382SShri Abhyankar           nextcol = PetscMin(nextcol, ii[j][0] / bs);
361a5973382SShri Abhyankar         }
362a5973382SShri Abhyankar       }
36308401ef6SPierre Jolivet       PetscCheck(cnt < amat->maxnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of blocks found greater than expected %" PetscInt_FMT, cnt);
364a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
365a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
3669566063dSJacob Faibussowitsch         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, llens, amat->a + cnt++));
367a5973382SShri Abhyankar       }
368a5973382SShri Abhyankar 
369a5973382SShri Abhyankar       if (!notdone) break;
370a5973382SShri Abhyankar       currentcol = nextcol;
371a5973382SShri Abhyankar     }
372a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
373a5973382SShri Abhyankar   }
374a5973382SShri Abhyankar 
3759566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lens, ii, ilens));
3769566063dSJacob Faibussowitsch   PetscCall(PetscFree(llens));
377a5973382SShri Abhyankar 
378a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
379a5973382SShri Abhyankar   for (i = 0; i < m; i++) {
3809566063dSJacob Faibussowitsch     PetscCall(MatGetRow(tmpA, i, &ncols, &cols, &values));
3819566063dSJacob Faibussowitsch     PetscCall(MatSetValues(newmat, 1, &i, ncols, cols, values, INSERT_VALUES));
3829566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(tmpA, i, &ncols, &cols, &values));
383a5973382SShri Abhyankar   }
3849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(newmat, MAT_FINAL_ASSEMBLY));
3859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(newmat, MAT_FINAL_ASSEMBLY));
3863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
387a5973382SShri Abhyankar }
388a5973382SShri Abhyankar 
389d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_BlockMat(Mat A, PetscViewer viewer)
390d71ae5a4SJacob Faibussowitsch {
39164075487SBarry Smith   Mat_BlockMat     *a = (Mat_BlockMat *)A->data;
392ccb205f8SBarry Smith   const char       *name;
393ccb205f8SBarry Smith   PetscViewerFormat format;
394ccb205f8SBarry Smith 
395ccb205f8SBarry Smith   PetscFunctionBegin;
3969566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)A, &name));
3979566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
398ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
3999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Nonzero block matrices = %" PetscInt_FMT " \n", a->nz));
400b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(PetscViewerASCIIPrintf(viewer, "Only upper triangular part of symmetric matrix is stored\n"));
401ccb205f8SBarry Smith   }
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
403ccb205f8SBarry Smith }
404421e10b8SBarry Smith 
405d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_BlockMat(Mat mat)
406d71ae5a4SJacob Faibussowitsch {
407421e10b8SBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)mat->data;
408ccb205f8SBarry Smith   PetscInt      i;
409421e10b8SBarry Smith 
410421e10b8SBarry Smith   PetscFunctionBegin;
4119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->right));
4129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->left));
4139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->middle));
4149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->workb));
415ccb205f8SBarry Smith   if (bmat->diags) {
41648a46eb9SPierre Jolivet     for (i = 0; i < mat->rmap->n / mat->rmap->bs; i++) PetscCall(MatDestroy(&bmat->diags[i]));
417ccb205f8SBarry Smith   }
418ccb205f8SBarry Smith   if (bmat->a) {
41948a46eb9SPierre Jolivet     for (i = 0; i < bmat->nz; i++) PetscCall(MatDestroy(&bmat->a[i]));
420ccb205f8SBarry Smith   }
4219566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(mat, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i));
4229566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424421e10b8SBarry Smith }
425421e10b8SBarry Smith 
426d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_BlockMat(Mat A, Vec x, Vec y)
427d71ae5a4SJacob Faibussowitsch {
428421e10b8SBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
429421e10b8SBarry Smith   PetscScalar  *xx, *yy;
430d0f46423SBarry Smith   PetscInt     *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j;
431421e10b8SBarry Smith   Mat          *aa;
432421e10b8SBarry Smith 
433421e10b8SBarry Smith   PetscFunctionBegin;
434421e10b8SBarry Smith   /*
435421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
436421e10b8SBarry Smith   */
4379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
438ccb205f8SBarry Smith 
4399566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4409566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
441421e10b8SBarry Smith   aj = bmat->j;
442421e10b8SBarry Smith   aa = bmat->a;
443421e10b8SBarry Smith   ii = bmat->i;
444421e10b8SBarry Smith   for (i = 0; i < m; i++) {
445421e10b8SBarry Smith     jrow = ii[i];
4469566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left, yy + bs * i));
447421e10b8SBarry Smith     n = ii[i + 1] - jrow;
448421e10b8SBarry Smith     for (j = 0; j < n; j++) {
4499566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow]));
4509566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4519566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
452421e10b8SBarry Smith       jrow++;
453421e10b8SBarry Smith     }
4549566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
455421e10b8SBarry Smith   }
4569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
4579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
459421e10b8SBarry Smith }
460421e10b8SBarry Smith 
46166976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_BlockMat_Symmetric(Mat A, Vec x, Vec y)
462d71ae5a4SJacob Faibussowitsch {
463290bbb0aSBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
464290bbb0aSBarry Smith   PetscScalar  *xx, *yy;
465d0f46423SBarry Smith   PetscInt     *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j;
466290bbb0aSBarry Smith   Mat          *aa;
467290bbb0aSBarry Smith 
468290bbb0aSBarry Smith   PetscFunctionBegin;
469290bbb0aSBarry Smith   /*
470290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
471290bbb0aSBarry Smith   */
4729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
473290bbb0aSBarry Smith 
4749566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
476290bbb0aSBarry Smith   aj = bmat->j;
477290bbb0aSBarry Smith   aa = bmat->a;
478290bbb0aSBarry Smith   ii = bmat->i;
479290bbb0aSBarry Smith   for (i = 0; i < m; i++) {
480290bbb0aSBarry Smith     jrow = ii[i];
481290bbb0aSBarry Smith     n    = ii[i + 1] - jrow;
4829566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left, yy + bs * i));
4839566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->middle, xx + bs * i));
484290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
485290bbb0aSBarry Smith     if (aj[jrow] == i) {
4869566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow]));
4879566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4889566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
489290bbb0aSBarry Smith       jrow++;
490290bbb0aSBarry Smith       n--;
491290bbb0aSBarry Smith     }
492290bbb0aSBarry Smith     for (j = 0; j < n; j++) {
4939566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow])); /* upper triangular part */
4949566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4959566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
496290bbb0aSBarry Smith 
4979566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, yy + bs * aj[jrow])); /* lower triangular part */
4989566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(aa[jrow], bmat->middle, bmat->right, bmat->right));
4999566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
500290bbb0aSBarry Smith       jrow++;
501290bbb0aSBarry Smith     }
5029566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
5039566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->middle));
504290bbb0aSBarry Smith   }
5059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
5069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
5073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
508290bbb0aSBarry Smith }
509290bbb0aSBarry Smith 
510d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_BlockMat(Mat A, Vec x, Vec y, Vec z)
511d71ae5a4SJacob Faibussowitsch {
512421e10b8SBarry Smith   PetscFunctionBegin;
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
514421e10b8SBarry Smith }
515421e10b8SBarry Smith 
516d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_BlockMat(Mat A, Vec x, Vec y)
517d71ae5a4SJacob Faibussowitsch {
518421e10b8SBarry Smith   PetscFunctionBegin;
5193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
520421e10b8SBarry Smith }
521421e10b8SBarry Smith 
522d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A, Vec x, Vec y, Vec z)
523d71ae5a4SJacob Faibussowitsch {
524421e10b8SBarry Smith   PetscFunctionBegin;
5253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
526421e10b8SBarry Smith }
527421e10b8SBarry Smith 
528ccb205f8SBarry Smith /*
529ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
530ccb205f8SBarry Smith */
531d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
532d71ae5a4SJacob Faibussowitsch {
533ccb205f8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
534d0f46423SBarry Smith   PetscInt      i, j, mbs = A->rmap->n / A->rmap->bs;
535ccb205f8SBarry Smith 
536ccb205f8SBarry Smith   PetscFunctionBegin;
53748a46eb9SPierre Jolivet   if (!a->diag) PetscCall(PetscMalloc1(mbs, &a->diag));
538ccb205f8SBarry Smith   for (i = 0; i < mbs; i++) {
539ccb205f8SBarry Smith     a->diag[i] = a->i[i + 1];
540ccb205f8SBarry Smith     for (j = a->i[i]; j < a->i[i + 1]; j++) {
541ccb205f8SBarry Smith       if (a->j[j] == i) {
542ccb205f8SBarry Smith         a->diag[i] = j;
543ccb205f8SBarry Smith         break;
544ccb205f8SBarry Smith       }
545ccb205f8SBarry Smith     }
546ccb205f8SBarry Smith   }
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
548ccb205f8SBarry Smith }
549ccb205f8SBarry Smith 
550d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A, IS isrow, IS iscol, MatReuse scall, Mat *B)
551d71ae5a4SJacob Faibussowitsch {
552ccb205f8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
553ccb205f8SBarry Smith   Mat_SeqAIJ   *c;
554ccb205f8SBarry Smith   PetscInt      i, k, first, step, lensi, nrows, ncols;
555d2a212eaSBarry Smith   PetscInt     *j_new, *i_new, *aj = a->j, *ailen = a->ilen;
5561620fd73SBarry Smith   PetscScalar  *a_new;
557ccb205f8SBarry Smith   Mat           C, *aa = a->a;
558ace3abfcSBarry Smith   PetscBool     stride, equal;
559ccb205f8SBarry Smith 
560ccb205f8SBarry Smith   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(ISEqual(isrow, iscol, &equal));
56228b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for identical column and row indices");
5639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride));
56428b400f6SJacob Faibussowitsch   PetscCheck(stride, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for stride indices");
5659566063dSJacob Faibussowitsch   PetscCall(ISStrideGetInfo(iscol, &first, &step));
56608401ef6SPierre Jolivet   PetscCheck(step == A->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only select one entry from each block");
567ccb205f8SBarry Smith 
5689566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow, &nrows));
569ccb205f8SBarry Smith   ncols = nrows;
570ccb205f8SBarry Smith 
571ccb205f8SBarry Smith   /* create submatrix */
572ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
573ccb205f8SBarry Smith     PetscInt n_cols, n_rows;
574ccb205f8SBarry Smith     C = *B;
5759566063dSJacob Faibussowitsch     PetscCall(MatGetSize(C, &n_rows, &n_cols));
576aed4548fSBarry Smith     PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size");
5779566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(C));
578ccb205f8SBarry Smith   } else {
5799566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
5809566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
581b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetType(C, MATSEQSBAIJ));
582b94d7dedSBarry Smith     else PetscCall(MatSetType(C, MATSEQAIJ));
5839566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(C, 0, ailen));
5849566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(C, 1, 0, ailen));
585ccb205f8SBarry Smith   }
586ccb205f8SBarry Smith   c = (Mat_SeqAIJ *)C->data;
587ccb205f8SBarry Smith 
588ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
589ccb205f8SBarry Smith   a_new = c->a;
590ccb205f8SBarry Smith   j_new = c->j;
591ccb205f8SBarry Smith   i_new = c->i;
592ccb205f8SBarry Smith 
593ccb205f8SBarry Smith   for (i = 0; i < nrows; i++) {
594ccb205f8SBarry Smith     lensi = ailen[i];
595ccb205f8SBarry Smith     for (k = 0; k < lensi; k++) {
596ccb205f8SBarry Smith       *j_new++ = *aj++;
5979566063dSJacob Faibussowitsch       PetscCall(MatGetValue(*aa++, first, first, a_new++));
598ccb205f8SBarry Smith     }
599ccb205f8SBarry Smith     i_new[i + 1] = i_new[i] + lensi;
600ccb205f8SBarry Smith     c->ilen[i]   = lensi;
601ccb205f8SBarry Smith   }
602ccb205f8SBarry Smith 
6039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
6049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
605ccb205f8SBarry Smith   *B = C;
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607ccb205f8SBarry Smith }
608ccb205f8SBarry Smith 
609d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A, MatAssemblyType mode)
610d71ae5a4SJacob Faibussowitsch {
611421e10b8SBarry Smith   Mat_BlockMat *a      = (Mat_BlockMat *)A->data;
612421e10b8SBarry Smith   PetscInt      fshift = 0, i, j, *ai = a->i, *aj = a->j, *imax = a->imax;
613ccb205f8SBarry Smith   PetscInt      m = a->mbs, *ip, N, *ailen = a->ilen, rmax = 0;
614421e10b8SBarry Smith   Mat          *aa = a->a, *ap;
615421e10b8SBarry Smith 
616421e10b8SBarry Smith   PetscFunctionBegin;
6173ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
618421e10b8SBarry Smith 
619421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
620421e10b8SBarry Smith   for (i = 1; i < m; i++) {
621421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
622421e10b8SBarry Smith     fshift += imax[i - 1] - ailen[i - 1];
623421e10b8SBarry Smith     rmax = PetscMax(rmax, ailen[i]);
624421e10b8SBarry Smith     if (fshift) {
625421e10b8SBarry Smith       ip = aj + ai[i];
626421e10b8SBarry Smith       ap = aa + ai[i];
627421e10b8SBarry Smith       N  = ailen[i];
628421e10b8SBarry Smith       for (j = 0; j < N; j++) {
629421e10b8SBarry Smith         ip[j - fshift] = ip[j];
630421e10b8SBarry Smith         ap[j - fshift] = ap[j];
631421e10b8SBarry Smith       }
632421e10b8SBarry Smith     }
633421e10b8SBarry Smith     ai[i] = ai[i - 1] + ailen[i - 1];
634421e10b8SBarry Smith   }
635421e10b8SBarry Smith   if (m) {
636421e10b8SBarry Smith     fshift += imax[m - 1] - ailen[m - 1];
637421e10b8SBarry Smith     ai[m] = ai[m - 1] + ailen[m - 1];
638421e10b8SBarry Smith   }
639421e10b8SBarry Smith   /* reset ilen and imax for each row */
640ad540459SPierre Jolivet   for (i = 0; i < m; i++) ailen[i] = imax[i] = ai[i + 1] - ai[i];
641421e10b8SBarry Smith   a->nz = ai[m];
642ccb205f8SBarry Smith   for (i = 0; i < a->nz; i++) {
6436bdcaf15SBarry Smith     PetscAssert(aa[i], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Null matrix at location %" PetscInt_FMT " column %" PetscInt_FMT " nz %" PetscInt_FMT, i, aj[i], a->nz);
6449566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(aa[i], MAT_FINAL_ASSEMBLY));
6459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(aa[i], MAT_FINAL_ASSEMBLY));
646ccb205f8SBarry Smith   }
6479566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded,%" PetscInt_FMT " used\n", m, A->cmap->n / A->cmap->bs, fshift, a->nz));
6489566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
6499566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax));
6502205254eSKarl Rupp 
6518e58a170SBarry Smith   A->info.mallocs += a->reallocs;
652421e10b8SBarry Smith   a->reallocs         = 0;
653421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
654421e10b8SBarry Smith   a->rmax             = rmax;
6559566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_BlockMat(A));
6563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
657421e10b8SBarry Smith }
658421e10b8SBarry Smith 
659d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetOption_BlockMat(Mat A, MatOption opt, PetscBool flg)
660d71ae5a4SJacob Faibussowitsch {
661290bbb0aSBarry Smith   PetscFunctionBegin;
6624e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
66341f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
664290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
665290bbb0aSBarry Smith   } else {
6669566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Unused matrix option %s\n", MatOptions[opt]));
667290bbb0aSBarry Smith   }
6683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
669290bbb0aSBarry Smith }
670290bbb0aSBarry Smith 
6713964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
672f4259b30SLisandro Dalcin                                        NULL,
673f4259b30SLisandro Dalcin                                        NULL,
674421e10b8SBarry Smith                                        MatMult_BlockMat,
675421e10b8SBarry Smith                                        /*  4*/ MatMultAdd_BlockMat,
676421e10b8SBarry Smith                                        MatMultTranspose_BlockMat,
677421e10b8SBarry Smith                                        MatMultTransposeAdd_BlockMat,
678f4259b30SLisandro Dalcin                                        NULL,
679f4259b30SLisandro Dalcin                                        NULL,
680f4259b30SLisandro Dalcin                                        NULL,
681f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
682f4259b30SLisandro Dalcin                                        NULL,
683f4259b30SLisandro Dalcin                                        NULL,
68441f059aeSBarry Smith                                        MatSOR_BlockMat,
685f4259b30SLisandro Dalcin                                        NULL,
686f4259b30SLisandro Dalcin                                        /* 15*/ NULL,
687f4259b30SLisandro Dalcin                                        NULL,
688f4259b30SLisandro Dalcin                                        NULL,
689f4259b30SLisandro Dalcin                                        NULL,
690f4259b30SLisandro Dalcin                                        NULL,
691f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
692421e10b8SBarry Smith                                        MatAssemblyEnd_BlockMat,
693290bbb0aSBarry Smith                                        MatSetOption_BlockMat,
694f4259b30SLisandro Dalcin                                        NULL,
695f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
696f4259b30SLisandro Dalcin                                        NULL,
697f4259b30SLisandro Dalcin                                        NULL,
698f4259b30SLisandro Dalcin                                        NULL,
699f4259b30SLisandro Dalcin                                        NULL,
700f4259b30SLisandro Dalcin                                        /* 29*/ NULL,
701f4259b30SLisandro Dalcin                                        NULL,
702f4259b30SLisandro Dalcin                                        NULL,
703f4259b30SLisandro Dalcin                                        NULL,
704f4259b30SLisandro Dalcin                                        NULL,
705f4259b30SLisandro Dalcin                                        /* 34*/ NULL,
706f4259b30SLisandro Dalcin                                        NULL,
707f4259b30SLisandro Dalcin                                        NULL,
708f4259b30SLisandro Dalcin                                        NULL,
709f4259b30SLisandro Dalcin                                        NULL,
710f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
711f4259b30SLisandro Dalcin                                        NULL,
712f4259b30SLisandro Dalcin                                        NULL,
713f4259b30SLisandro Dalcin                                        NULL,
714f4259b30SLisandro Dalcin                                        NULL,
715f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
716f4259b30SLisandro Dalcin                                        NULL,
7177d68702bSBarry Smith                                        MatShift_Basic,
718f4259b30SLisandro Dalcin                                        NULL,
719f4259b30SLisandro Dalcin                                        NULL,
720f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
721f4259b30SLisandro Dalcin                                        NULL,
722f4259b30SLisandro Dalcin                                        NULL,
723f4259b30SLisandro Dalcin                                        NULL,
724f4259b30SLisandro Dalcin                                        NULL,
725f4259b30SLisandro Dalcin                                        /* 54*/ NULL,
726f4259b30SLisandro Dalcin                                        NULL,
727f4259b30SLisandro Dalcin                                        NULL,
728f4259b30SLisandro Dalcin                                        NULL,
729f4259b30SLisandro Dalcin                                        NULL,
7307dae84e0SHong Zhang                                        /* 59*/ MatCreateSubMatrix_BlockMat,
731421e10b8SBarry Smith                                        MatDestroy_BlockMat,
732ccb205f8SBarry Smith                                        MatView_BlockMat,
733f4259b30SLisandro Dalcin                                        NULL,
734f4259b30SLisandro Dalcin                                        NULL,
735f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
736f4259b30SLisandro Dalcin                                        NULL,
737f4259b30SLisandro Dalcin                                        NULL,
738f4259b30SLisandro Dalcin                                        NULL,
739f4259b30SLisandro Dalcin                                        NULL,
740f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
741f4259b30SLisandro Dalcin                                        NULL,
742f4259b30SLisandro Dalcin                                        NULL,
743f4259b30SLisandro Dalcin                                        NULL,
744f4259b30SLisandro Dalcin                                        NULL,
745f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
746f4259b30SLisandro Dalcin                                        NULL,
747f4259b30SLisandro Dalcin                                        NULL,
748f4259b30SLisandro Dalcin                                        NULL,
749f4259b30SLisandro Dalcin                                        NULL,
750f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
751f4259b30SLisandro Dalcin                                        NULL,
752f4259b30SLisandro Dalcin                                        NULL,
753f4259b30SLisandro Dalcin                                        NULL,
7545bba2384SShri Abhyankar                                        MatLoad_BlockMat,
755f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
756f4259b30SLisandro Dalcin                                        NULL,
757f4259b30SLisandro Dalcin                                        NULL,
758f4259b30SLisandro Dalcin                                        NULL,
759f4259b30SLisandro Dalcin                                        NULL,
760f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
761f4259b30SLisandro Dalcin                                        NULL,
762f4259b30SLisandro Dalcin                                        NULL,
763f4259b30SLisandro Dalcin                                        NULL,
764f4259b30SLisandro Dalcin                                        NULL,
765f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
766f4259b30SLisandro Dalcin                                        NULL,
767f4259b30SLisandro Dalcin                                        NULL,
768f4259b30SLisandro Dalcin                                        NULL,
769f4259b30SLisandro Dalcin                                        NULL,
770f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
771f4259b30SLisandro Dalcin                                        NULL,
772f4259b30SLisandro Dalcin                                        NULL,
773f4259b30SLisandro Dalcin                                        NULL,
774f4259b30SLisandro Dalcin                                        NULL,
775f4259b30SLisandro Dalcin                                        /*104*/ NULL,
776f4259b30SLisandro Dalcin                                        NULL,
777f4259b30SLisandro Dalcin                                        NULL,
778f4259b30SLisandro Dalcin                                        NULL,
779f4259b30SLisandro Dalcin                                        NULL,
780f4259b30SLisandro Dalcin                                        /*109*/ NULL,
781f4259b30SLisandro Dalcin                                        NULL,
782f4259b30SLisandro Dalcin                                        NULL,
783f4259b30SLisandro Dalcin                                        NULL,
784f4259b30SLisandro Dalcin                                        NULL,
785f4259b30SLisandro Dalcin                                        /*114*/ NULL,
786f4259b30SLisandro Dalcin                                        NULL,
787f4259b30SLisandro Dalcin                                        NULL,
788f4259b30SLisandro Dalcin                                        NULL,
789f4259b30SLisandro Dalcin                                        NULL,
790f4259b30SLisandro Dalcin                                        /*119*/ NULL,
791f4259b30SLisandro Dalcin                                        NULL,
792f4259b30SLisandro Dalcin                                        NULL,
793f4259b30SLisandro Dalcin                                        NULL,
794f4259b30SLisandro Dalcin                                        NULL,
795f4259b30SLisandro Dalcin                                        /*124*/ NULL,
796f4259b30SLisandro Dalcin                                        NULL,
797f4259b30SLisandro Dalcin                                        NULL,
798f4259b30SLisandro Dalcin                                        NULL,
799f4259b30SLisandro Dalcin                                        NULL,
800f4259b30SLisandro Dalcin                                        /*129*/ NULL,
801f4259b30SLisandro Dalcin                                        NULL,
802f4259b30SLisandro Dalcin                                        NULL,
803f4259b30SLisandro Dalcin                                        NULL,
804f4259b30SLisandro Dalcin                                        NULL,
805f4259b30SLisandro Dalcin                                        /*134*/ NULL,
806f4259b30SLisandro Dalcin                                        NULL,
807f4259b30SLisandro Dalcin                                        NULL,
808f4259b30SLisandro Dalcin                                        NULL,
809f4259b30SLisandro Dalcin                                        NULL,
810f4259b30SLisandro Dalcin                                        /*139*/ NULL,
811f4259b30SLisandro Dalcin                                        NULL,
812d70f29a3SPierre Jolivet                                        NULL,
813d70f29a3SPierre Jolivet                                        NULL,
814d70f29a3SPierre Jolivet                                        NULL,
815d70f29a3SPierre Jolivet                                        /*144*/ NULL,
816d70f29a3SPierre Jolivet                                        NULL,
817d70f29a3SPierre Jolivet                                        NULL,
81899a7f59eSMark Adams                                        NULL,
81999a7f59eSMark Adams                                        NULL,
8207fb60732SBarry Smith                                        NULL,
821dec0b466SHong Zhang                                        /*150*/ NULL,
822eede4a3fSMark Adams                                        NULL,
823*4cc2b5b5SPierre Jolivet                                        NULL,
824dec0b466SHong Zhang                                        NULL};
825421e10b8SBarry Smith 
8268cdf2d9bSBarry Smith /*@C
8278cdf2d9bSBarry Smith   MatBlockMatSetPreallocation - For good matrix assembly performance
8288cdf2d9bSBarry Smith   the user should preallocate the matrix storage by setting the parameter nz
8298cdf2d9bSBarry Smith   (or the array nnz).  By setting these parameters accurately, performance
8308cdf2d9bSBarry Smith   during matrix assembly can be increased by more than a factor of 50.
8318cdf2d9bSBarry Smith 
832d083f849SBarry Smith   Collective
8338cdf2d9bSBarry Smith 
8348cdf2d9bSBarry Smith   Input Parameters:
8358cdf2d9bSBarry Smith + B   - The matrix
8368cdf2d9bSBarry Smith . bs  - size of each block in matrix
8378cdf2d9bSBarry Smith . nz  - number of nonzeros per block row (same for all rows)
8388cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows
8392ef1f0ffSBarry Smith          (possibly different for each row) or `NULL`
8408cdf2d9bSBarry Smith 
8418cdf2d9bSBarry Smith   Level: intermediate
8428cdf2d9bSBarry Smith 
8432ef1f0ffSBarry Smith   Notes:
8442ef1f0ffSBarry Smith   If `nnz` is given then `nz` is ignored
8452ef1f0ffSBarry Smith 
8462ef1f0ffSBarry Smith   Specify the preallocated storage with either `nz` or `nnz` (not both).
8472ef1f0ffSBarry Smith   Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
8482ef1f0ffSBarry Smith   allocation.
8492ef1f0ffSBarry Smith 
8501cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateBlockMat()`, `MatSetValues()`
8518cdf2d9bSBarry Smith @*/
852d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBlockMatSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[])
853d71ae5a4SJacob Faibussowitsch {
8548cdf2d9bSBarry Smith   PetscFunctionBegin;
855cac4c232SBarry Smith   PetscTryMethod(B, "MatBlockMatSetPreallocation_C", (Mat, PetscInt, PetscInt, const PetscInt[]), (B, bs, nz, nnz));
8563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8578cdf2d9bSBarry Smith }
8588cdf2d9bSBarry Smith 
859d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatBlockMatSetPreallocation_BlockMat(Mat A, PetscInt bs, PetscInt nz, PetscInt *nnz)
860d71ae5a4SJacob Faibussowitsch {
8618cdf2d9bSBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
8628cdf2d9bSBarry Smith   PetscInt      i;
8638cdf2d9bSBarry Smith 
8648cdf2d9bSBarry Smith   PetscFunctionBegin;
8659566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(A->rmap, bs));
8669566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(A->cmap, bs));
8679566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
8689566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
8699566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(A->rmap, &bs));
87034ef9618SShri Abhyankar 
8718cdf2d9bSBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
87208401ef6SPierre Jolivet   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
8738cdf2d9bSBarry Smith   if (nnz) {
874d0f46423SBarry Smith     for (i = 0; i < A->rmap->n / bs; i++) {
87508401ef6SPierre Jolivet       PetscCheck(nnz[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT, i, nnz[i]);
87608401ef6SPierre Jolivet       PetscCheck(nnz[i] <= A->cmap->n / bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be greater than row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT, i, nnz[i], A->cmap->n / bs);
8778cdf2d9bSBarry Smith     }
8788cdf2d9bSBarry Smith   }
879d0f46423SBarry Smith   bmat->mbs = A->rmap->n / bs;
8808cdf2d9bSBarry Smith 
8819566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs, NULL, &bmat->right));
8829566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs, NULL, &bmat->middle));
8839566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, bs, &bmat->left));
8848cdf2d9bSBarry Smith 
885aa624791SPierre Jolivet   if (!bmat->imax) PetscCall(PetscMalloc2(A->rmap->n, &bmat->imax, A->rmap->n, &bmat->ilen));
886c0aa6a63SJacob Faibussowitsch   if (PetscLikely(nnz)) {
8878cdf2d9bSBarry Smith     nz = 0;
888d0f46423SBarry Smith     for (i = 0; i < A->rmap->n / A->rmap->bs; i++) {
8898cdf2d9bSBarry Smith       bmat->imax[i] = nnz[i];
8908cdf2d9bSBarry Smith       nz += nnz[i];
8918cdf2d9bSBarry Smith     }
892f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently requires block row by row preallocation");
8938cdf2d9bSBarry Smith 
8948cdf2d9bSBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
8959566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(bmat->ilen, bmat->mbs));
8968cdf2d9bSBarry Smith 
8978cdf2d9bSBarry Smith   /* allocate the matrix space */
8989566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(A, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i));
8999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nz, &bmat->a, nz, &bmat->j, A->rmap->n + 1, &bmat->i));
9008cdf2d9bSBarry Smith   bmat->i[0] = 0;
901ad540459SPierre Jolivet   for (i = 1; i < bmat->mbs + 1; i++) bmat->i[i] = bmat->i[i - 1] + bmat->imax[i - 1];
9028cdf2d9bSBarry Smith   bmat->singlemalloc = PETSC_TRUE;
9038cdf2d9bSBarry Smith   bmat->free_a       = PETSC_TRUE;
9048cdf2d9bSBarry Smith   bmat->free_ij      = PETSC_TRUE;
9058cdf2d9bSBarry Smith 
9068cdf2d9bSBarry Smith   bmat->nz            = 0;
9078cdf2d9bSBarry Smith   bmat->maxnz         = nz;
9088cdf2d9bSBarry Smith   A->info.nz_unneeded = (double)bmat->maxnz;
9099566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
9103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9118cdf2d9bSBarry Smith }
9128cdf2d9bSBarry Smith 
913421e10b8SBarry Smith /*MC
91411a5261eSBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of `Mat`'s that represents a sparse block matrix
915421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
916421e10b8SBarry Smith 
917421e10b8SBarry Smith   Level: advanced
918421e10b8SBarry Smith 
9191cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateBlockMat()`
920421e10b8SBarry Smith M*/
921421e10b8SBarry Smith 
922d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A)
923d71ae5a4SJacob Faibussowitsch {
924421e10b8SBarry Smith   Mat_BlockMat *b;
925421e10b8SBarry Smith 
926421e10b8SBarry Smith   PetscFunctionBegin;
9274dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
928421e10b8SBarry Smith   A->data         = (void *)b;
929aea10558SJacob Faibussowitsch   A->ops[0]       = MatOps_Values;
930421e10b8SBarry Smith   A->assembled    = PETSC_TRUE;
931421e10b8SBarry Smith   A->preallocated = PETSC_FALSE;
9329566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATBLOCKMAT));
933290bbb0aSBarry Smith 
9349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatBlockMatSetPreallocation_C", MatBlockMatSetPreallocation_BlockMat));
9353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
936421e10b8SBarry Smith }
937421e10b8SBarry Smith 
938421e10b8SBarry Smith /*@C
93911a5261eSBarry Smith   MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential `Mat` object
940421e10b8SBarry Smith 
941d083f849SBarry Smith   Collective
942421e10b8SBarry Smith 
943421e10b8SBarry Smith   Input Parameters:
944421e10b8SBarry Smith + comm - MPI communicator
945421e10b8SBarry Smith . m    - number of rows
946421e10b8SBarry Smith . n    - number of columns
9478cdf2d9bSBarry Smith . bs   - size of each submatrix
94811a5261eSBarry Smith . nz   - expected maximum number of nonzero blocks in row (use `PETSC_DEFAULT` if not known)
9492ef1f0ffSBarry Smith - nnz  - expected number of nonzers per block row if known (use `NULL` otherwise)
950421e10b8SBarry Smith 
951421e10b8SBarry Smith   Output Parameter:
952421e10b8SBarry Smith . A - the matrix
953421e10b8SBarry Smith 
954421e10b8SBarry Smith   Level: intermediate
955421e10b8SBarry Smith 
95695452b02SPatrick Sanan   Notes:
95711a5261eSBarry Smith   Matrices of this type are nominally-sparse matrices in which each "entry" is a `Mat` object.  Each `Mat` must
958bbd3f336SJed Brown   have the same size and be sequential.  The local and global sizes must be compatible with this decomposition.
959bbd3f336SJed Brown 
96011a5261eSBarry Smith   For matrices containing parallel submatrices and variable block sizes, see `MATNEST`.
96111a5261eSBarry Smith 
962fe59aa6dSJacob Faibussowitsch   Developer Notes:
96311a5261eSBarry Smith   I don't like the name, it is not `MATNESTMAT`
964421e10b8SBarry Smith 
9651cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATBLOCKMAT`, `MatCreateNest()`
966421e10b8SBarry Smith @*/
967d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateBlockMat(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt bs, PetscInt nz, PetscInt *nnz, Mat *A)
968d71ae5a4SJacob Faibussowitsch {
969421e10b8SBarry Smith   PetscFunctionBegin;
9709566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
9719566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
9729566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATBLOCKMAT));
9739566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(*A, bs, nz, nnz));
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
975421e10b8SBarry Smith }
976