xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
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));
53290bbb0aSBarry Smith   /* 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 
7841f059aeSBarry Smith         /* 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++;
265e56f5c9eSBarry Smith       A->nonzerostate++;
266421e10b8SBarry Smith     noinsert1:;
26748a46eb9SPierre Jolivet       if (!*(ap + i)) PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, NULL, ap + i));
2689566063dSJacob Faibussowitsch       PetscCall(MatSetValues(ap[i], 1, &ridx, 1, &cidx, &value, is));
269421e10b8SBarry Smith       low = i;
270421e10b8SBarry Smith     }
271421e10b8SBarry Smith     ailen[brow] = nrow;
272421e10b8SBarry Smith   }
2733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
274421e10b8SBarry Smith }
275421e10b8SBarry Smith 
276d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
277d71ae5a4SJacob Faibussowitsch {
278a5973382SShri Abhyankar   Mat                tmpA;
279a5973382SShri Abhyankar   PetscInt           i, j, m, n, bs = 1, ncols, *lens, currentcol, mbs, **ii, *ilens, nextcol, *llens, cnt = 0;
280a5973382SShri Abhyankar   const PetscInt    *cols;
281a5973382SShri Abhyankar   const PetscScalar *values;
282ace3abfcSBarry Smith   PetscBool          flg = PETSC_FALSE, notdone;
283a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
284a5973382SShri Abhyankar   Mat_BlockMat      *amat;
285a5973382SShri Abhyankar 
286a5973382SShri Abhyankar   PetscFunctionBegin;
287c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
2889566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
2899566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &tmpA));
2909566063dSJacob Faibussowitsch   PetscCall(MatSetType(tmpA, MATSEQAIJ));
2919566063dSJacob Faibussowitsch   PetscCall(MatLoad_SeqAIJ(tmpA, viewer));
292a5973382SShri Abhyankar 
2939566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(tmpA, &m, &n));
294d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, NULL, "Options for loading BlockMat matrix 1", "Mat");
2959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-matload_block_size", "Set the blocksize used to store the matrix", "MatLoad", bs, &bs, NULL));
2969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-matload_symmetric", "Store the matrix as symmetric", "MatLoad", flg, &flg, NULL));
297d0609cedSBarry Smith   PetscOptionsEnd();
298a5973382SShri Abhyankar 
299a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
300a5973382SShri Abhyankar   a   = (Mat_SeqAIJ *)tmpA->data;
301a5973382SShri Abhyankar   mbs = m / bs;
3029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &lens, bs, &ii, bs, &ilens));
3039566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lens, mbs));
304a5973382SShri Abhyankar 
305a5973382SShri Abhyankar   for (i = 0; i < mbs; i++) {
306a5973382SShri Abhyankar     for (j = 0; j < bs; j++) {
307a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i * bs + j];
308a5973382SShri Abhyankar       ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j];
309a5973382SShri Abhyankar     }
310a5973382SShri Abhyankar 
311a5973382SShri Abhyankar     currentcol = -1;
312a5973382SShri Abhyankar     while (PETSC_TRUE) {
313a5973382SShri Abhyankar       notdone = PETSC_FALSE;
314a5973382SShri Abhyankar       nextcol = 1000000000;
315a5973382SShri Abhyankar       for (j = 0; j < bs; j++) {
316*f4f49eeaSPierre Jolivet         while (ilens[j] > 0 && ii[j][0] / bs <= currentcol) {
317a5973382SShri Abhyankar           ii[j]++;
318a5973382SShri Abhyankar           ilens[j]--;
319a5973382SShri Abhyankar         }
320a5973382SShri Abhyankar         if (ilens[j] > 0) {
321a5973382SShri Abhyankar           notdone = PETSC_TRUE;
322a5973382SShri Abhyankar           nextcol = PetscMin(nextcol, ii[j][0] / bs);
323a5973382SShri Abhyankar         }
324a5973382SShri Abhyankar       }
325a5973382SShri Abhyankar       if (!notdone) break;
326a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
327a5973382SShri Abhyankar       currentcol = nextcol;
328a5973382SShri Abhyankar     }
329a5973382SShri Abhyankar   }
330a5973382SShri Abhyankar 
33148a46eb9SPierre 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));
3329566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(newmat, bs, 0, lens));
3331baa6e33SBarry Smith   if (flg) PetscCall(MatSetOption(newmat, MAT_SYMMETRIC, PETSC_TRUE));
334a5973382SShri Abhyankar   amat = (Mat_BlockMat *)(newmat)->data;
335a5973382SShri Abhyankar 
336a5973382SShri Abhyankar   /* preallocate the submatrices */
3379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &llens));
338a5973382SShri Abhyankar   for (i = 0; i < mbs; i++) { /* loops for block rows */
339a5973382SShri Abhyankar     for (j = 0; j < bs; j++) {
340a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i * bs + j];
341a5973382SShri Abhyankar       ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j];
342a5973382SShri Abhyankar     }
343a5973382SShri Abhyankar 
344a5973382SShri Abhyankar     currentcol = 1000000000;
345a5973382SShri Abhyankar     for (j = 0; j < bs; j++) { /* loop over rows in block finding first nonzero block */
346ad540459SPierre Jolivet       if (ilens[j] > 0) currentcol = PetscMin(currentcol, ii[j][0] / bs);
347a5973382SShri Abhyankar     }
348a5973382SShri Abhyankar 
349a5973382SShri Abhyankar     while (PETSC_TRUE) { /* loops over blocks in block row */
350a5973382SShri Abhyankar       notdone = PETSC_FALSE;
351a5973382SShri Abhyankar       nextcol = 1000000000;
3529566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(llens, bs));
353a5973382SShri Abhyankar       for (j = 0; j < bs; j++) {                              /* loop over rows in block */
354*f4f49eeaSPierre Jolivet         while (ilens[j] > 0 && ii[j][0] / bs <= currentcol) { /* loop over columns in row */
355a5973382SShri Abhyankar           ii[j]++;
356a5973382SShri Abhyankar           ilens[j]--;
357a5973382SShri Abhyankar           llens[j]++;
358a5973382SShri Abhyankar         }
359a5973382SShri Abhyankar         if (ilens[j] > 0) {
360a5973382SShri Abhyankar           notdone = PETSC_TRUE;
361a5973382SShri Abhyankar           nextcol = PetscMin(nextcol, ii[j][0] / bs);
362a5973382SShri Abhyankar         }
363a5973382SShri Abhyankar       }
36408401ef6SPierre Jolivet       PetscCheck(cnt < amat->maxnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of blocks found greater than expected %" PetscInt_FMT, cnt);
365a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
366a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
3679566063dSJacob Faibussowitsch         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, llens, amat->a + cnt++));
368a5973382SShri Abhyankar       }
369a5973382SShri Abhyankar 
370a5973382SShri Abhyankar       if (!notdone) break;
371a5973382SShri Abhyankar       currentcol = nextcol;
372a5973382SShri Abhyankar     }
373a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
374a5973382SShri Abhyankar   }
375a5973382SShri Abhyankar 
3769566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lens, ii, ilens));
3779566063dSJacob Faibussowitsch   PetscCall(PetscFree(llens));
378a5973382SShri Abhyankar 
379a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
380a5973382SShri Abhyankar   for (i = 0; i < m; i++) {
3819566063dSJacob Faibussowitsch     PetscCall(MatGetRow(tmpA, i, &ncols, &cols, &values));
3829566063dSJacob Faibussowitsch     PetscCall(MatSetValues(newmat, 1, &i, ncols, cols, values, INSERT_VALUES));
3839566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(tmpA, i, &ncols, &cols, &values));
384a5973382SShri Abhyankar   }
3859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(newmat, MAT_FINAL_ASSEMBLY));
3869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(newmat, MAT_FINAL_ASSEMBLY));
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388a5973382SShri Abhyankar }
389a5973382SShri Abhyankar 
390d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_BlockMat(Mat A, PetscViewer viewer)
391d71ae5a4SJacob Faibussowitsch {
39264075487SBarry Smith   Mat_BlockMat     *a = (Mat_BlockMat *)A->data;
393ccb205f8SBarry Smith   const char       *name;
394ccb205f8SBarry Smith   PetscViewerFormat format;
395ccb205f8SBarry Smith 
396ccb205f8SBarry Smith   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)A, &name));
3989566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
399ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
4009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Nonzero block matrices = %" PetscInt_FMT " \n", a->nz));
401b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(PetscViewerASCIIPrintf(viewer, "Only upper triangular part of symmetric matrix is stored\n"));
402ccb205f8SBarry Smith   }
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
404ccb205f8SBarry Smith }
405421e10b8SBarry Smith 
406d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_BlockMat(Mat mat)
407d71ae5a4SJacob Faibussowitsch {
408421e10b8SBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)mat->data;
409ccb205f8SBarry Smith   PetscInt      i;
410421e10b8SBarry Smith 
411421e10b8SBarry Smith   PetscFunctionBegin;
4129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->right));
4139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->left));
4149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->middle));
4159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->workb));
416ccb205f8SBarry Smith   if (bmat->diags) {
41748a46eb9SPierre Jolivet     for (i = 0; i < mat->rmap->n / mat->rmap->bs; i++) PetscCall(MatDestroy(&bmat->diags[i]));
418ccb205f8SBarry Smith   }
419ccb205f8SBarry Smith   if (bmat->a) {
42048a46eb9SPierre Jolivet     for (i = 0; i < bmat->nz; i++) PetscCall(MatDestroy(&bmat->a[i]));
421ccb205f8SBarry Smith   }
4229566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(mat, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
425421e10b8SBarry Smith }
426421e10b8SBarry Smith 
427d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_BlockMat(Mat A, Vec x, Vec y)
428d71ae5a4SJacob Faibussowitsch {
429421e10b8SBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
430421e10b8SBarry Smith   PetscScalar  *xx, *yy;
431d0f46423SBarry Smith   PetscInt     *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j;
432421e10b8SBarry Smith   Mat          *aa;
433421e10b8SBarry Smith 
434421e10b8SBarry Smith   PetscFunctionBegin;
435421e10b8SBarry Smith   /*
436421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
437421e10b8SBarry Smith   */
4389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
439ccb205f8SBarry Smith 
4409566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
442421e10b8SBarry Smith   aj = bmat->j;
443421e10b8SBarry Smith   aa = bmat->a;
444421e10b8SBarry Smith   ii = bmat->i;
445421e10b8SBarry Smith   for (i = 0; i < m; i++) {
446421e10b8SBarry Smith     jrow = ii[i];
4479566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left, yy + bs * i));
448421e10b8SBarry Smith     n = ii[i + 1] - jrow;
449421e10b8SBarry Smith     for (j = 0; j < n; j++) {
4509566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow]));
4519566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4529566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
453421e10b8SBarry Smith       jrow++;
454421e10b8SBarry Smith     }
4559566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
456421e10b8SBarry Smith   }
4579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
4589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460421e10b8SBarry Smith }
461421e10b8SBarry Smith 
46266976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_BlockMat_Symmetric(Mat A, Vec x, Vec y)
463d71ae5a4SJacob Faibussowitsch {
464290bbb0aSBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
465290bbb0aSBarry Smith   PetscScalar  *xx, *yy;
466d0f46423SBarry Smith   PetscInt     *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j;
467290bbb0aSBarry Smith   Mat          *aa;
468290bbb0aSBarry Smith 
469290bbb0aSBarry Smith   PetscFunctionBegin;
470290bbb0aSBarry Smith   /*
471290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
472290bbb0aSBarry Smith   */
4739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
474290bbb0aSBarry Smith 
4759566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
477290bbb0aSBarry Smith   aj = bmat->j;
478290bbb0aSBarry Smith   aa = bmat->a;
479290bbb0aSBarry Smith   ii = bmat->i;
480290bbb0aSBarry Smith   for (i = 0; i < m; i++) {
481290bbb0aSBarry Smith     jrow = ii[i];
482290bbb0aSBarry Smith     n    = ii[i + 1] - jrow;
4839566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left, yy + bs * i));
4849566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->middle, xx + bs * i));
485290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
486290bbb0aSBarry Smith     if (aj[jrow] == i) {
4879566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow]));
4889566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4899566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
490290bbb0aSBarry Smith       jrow++;
491290bbb0aSBarry Smith       n--;
492290bbb0aSBarry Smith     }
493290bbb0aSBarry Smith     for (j = 0; j < n; j++) {
4949566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow])); /* upper triangular part */
4959566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4969566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
497290bbb0aSBarry Smith 
4989566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, yy + bs * aj[jrow])); /* lower triangular part */
4999566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(aa[jrow], bmat->middle, bmat->right, bmat->right));
5009566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
501290bbb0aSBarry Smith       jrow++;
502290bbb0aSBarry Smith     }
5039566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
5049566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->middle));
505290bbb0aSBarry Smith   }
5069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
5079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
5083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
509290bbb0aSBarry Smith }
510290bbb0aSBarry Smith 
511d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_BlockMat(Mat A, Vec x, Vec y, Vec z)
512d71ae5a4SJacob Faibussowitsch {
513421e10b8SBarry Smith   PetscFunctionBegin;
5143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515421e10b8SBarry Smith }
516421e10b8SBarry Smith 
517d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_BlockMat(Mat A, Vec x, Vec y)
518d71ae5a4SJacob Faibussowitsch {
519421e10b8SBarry Smith   PetscFunctionBegin;
5203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
521421e10b8SBarry Smith }
522421e10b8SBarry Smith 
523d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A, Vec x, Vec y, Vec z)
524d71ae5a4SJacob Faibussowitsch {
525421e10b8SBarry Smith   PetscFunctionBegin;
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527421e10b8SBarry Smith }
528421e10b8SBarry Smith 
529ccb205f8SBarry Smith /*
530ccb205f8SBarry Smith      Adds diagonal pointers to sparse matrix structure.
531ccb205f8SBarry Smith */
532d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
533d71ae5a4SJacob Faibussowitsch {
534ccb205f8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
535d0f46423SBarry Smith   PetscInt      i, j, mbs = A->rmap->n / A->rmap->bs;
536ccb205f8SBarry Smith 
537ccb205f8SBarry Smith   PetscFunctionBegin;
53848a46eb9SPierre Jolivet   if (!a->diag) PetscCall(PetscMalloc1(mbs, &a->diag));
539ccb205f8SBarry Smith   for (i = 0; i < mbs; i++) {
540ccb205f8SBarry Smith     a->diag[i] = a->i[i + 1];
541ccb205f8SBarry Smith     for (j = a->i[i]; j < a->i[i + 1]; j++) {
542ccb205f8SBarry Smith       if (a->j[j] == i) {
543ccb205f8SBarry Smith         a->diag[i] = j;
544ccb205f8SBarry Smith         break;
545ccb205f8SBarry Smith       }
546ccb205f8SBarry Smith     }
547ccb205f8SBarry Smith   }
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549ccb205f8SBarry Smith }
550ccb205f8SBarry Smith 
551d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A, IS isrow, IS iscol, MatReuse scall, Mat *B)
552d71ae5a4SJacob Faibussowitsch {
553ccb205f8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
554ccb205f8SBarry Smith   Mat_SeqAIJ   *c;
555ccb205f8SBarry Smith   PetscInt      i, k, first, step, lensi, nrows, ncols;
556d2a212eaSBarry Smith   PetscInt     *j_new, *i_new, *aj = a->j, *ailen = a->ilen;
5571620fd73SBarry Smith   PetscScalar  *a_new;
558ccb205f8SBarry Smith   Mat           C, *aa = a->a;
559ace3abfcSBarry Smith   PetscBool     stride, equal;
560ccb205f8SBarry Smith 
561ccb205f8SBarry Smith   PetscFunctionBegin;
5629566063dSJacob Faibussowitsch   PetscCall(ISEqual(isrow, iscol, &equal));
56328b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for identical column and row indices");
5649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride));
56528b400f6SJacob Faibussowitsch   PetscCheck(stride, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for stride indices");
5669566063dSJacob Faibussowitsch   PetscCall(ISStrideGetInfo(iscol, &first, &step));
56708401ef6SPierre Jolivet   PetscCheck(step == A->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only select one entry from each block");
568ccb205f8SBarry Smith 
5699566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow, &nrows));
570ccb205f8SBarry Smith   ncols = nrows;
571ccb205f8SBarry Smith 
572ccb205f8SBarry Smith   /* create submatrix */
573ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
574ccb205f8SBarry Smith     PetscInt n_cols, n_rows;
575ccb205f8SBarry Smith     C = *B;
5769566063dSJacob Faibussowitsch     PetscCall(MatGetSize(C, &n_rows, &n_cols));
577aed4548fSBarry Smith     PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size");
5789566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(C));
579ccb205f8SBarry Smith   } else {
5809566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
5819566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
582b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetType(C, MATSEQSBAIJ));
583b94d7dedSBarry Smith     else PetscCall(MatSetType(C, MATSEQAIJ));
5849566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(C, 0, ailen));
5859566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(C, 1, 0, ailen));
586ccb205f8SBarry Smith   }
587ccb205f8SBarry Smith   c = (Mat_SeqAIJ *)C->data;
588ccb205f8SBarry Smith 
589ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
590ccb205f8SBarry Smith   a_new = c->a;
591ccb205f8SBarry Smith   j_new = c->j;
592ccb205f8SBarry Smith   i_new = c->i;
593ccb205f8SBarry Smith 
594ccb205f8SBarry Smith   for (i = 0; i < nrows; i++) {
595ccb205f8SBarry Smith     lensi = ailen[i];
596ccb205f8SBarry Smith     for (k = 0; k < lensi; k++) {
597ccb205f8SBarry Smith       *j_new++ = *aj++;
5989566063dSJacob Faibussowitsch       PetscCall(MatGetValue(*aa++, first, first, a_new++));
599ccb205f8SBarry Smith     }
600ccb205f8SBarry Smith     i_new[i + 1] = i_new[i] + lensi;
601ccb205f8SBarry Smith     c->ilen[i]   = lensi;
602ccb205f8SBarry Smith   }
603ccb205f8SBarry Smith 
6049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
6059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
606ccb205f8SBarry Smith   *B = C;
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
608ccb205f8SBarry Smith }
609ccb205f8SBarry Smith 
610d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A, MatAssemblyType mode)
611d71ae5a4SJacob Faibussowitsch {
612421e10b8SBarry Smith   Mat_BlockMat *a      = (Mat_BlockMat *)A->data;
613421e10b8SBarry Smith   PetscInt      fshift = 0, i, j, *ai = a->i, *aj = a->j, *imax = a->imax;
614ccb205f8SBarry Smith   PetscInt      m = a->mbs, *ip, N, *ailen = a->ilen, rmax = 0;
615421e10b8SBarry Smith   Mat          *aa = a->a, *ap;
616421e10b8SBarry Smith 
617421e10b8SBarry Smith   PetscFunctionBegin;
6183ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
619421e10b8SBarry Smith 
620421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
621421e10b8SBarry Smith   for (i = 1; i < m; i++) {
622421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
623421e10b8SBarry Smith     fshift += imax[i - 1] - ailen[i - 1];
624421e10b8SBarry Smith     rmax = PetscMax(rmax, ailen[i]);
625421e10b8SBarry Smith     if (fshift) {
626421e10b8SBarry Smith       ip = aj + ai[i];
627421e10b8SBarry Smith       ap = aa + ai[i];
628421e10b8SBarry Smith       N  = ailen[i];
629421e10b8SBarry Smith       for (j = 0; j < N; j++) {
630421e10b8SBarry Smith         ip[j - fshift] = ip[j];
631421e10b8SBarry Smith         ap[j - fshift] = ap[j];
632421e10b8SBarry Smith       }
633421e10b8SBarry Smith     }
634421e10b8SBarry Smith     ai[i] = ai[i - 1] + ailen[i - 1];
635421e10b8SBarry Smith   }
636421e10b8SBarry Smith   if (m) {
637421e10b8SBarry Smith     fshift += imax[m - 1] - ailen[m - 1];
638421e10b8SBarry Smith     ai[m] = ai[m - 1] + ailen[m - 1];
639421e10b8SBarry Smith   }
640421e10b8SBarry Smith   /* reset ilen and imax for each row */
641ad540459SPierre Jolivet   for (i = 0; i < m; i++) ailen[i] = imax[i] = ai[i + 1] - ai[i];
642421e10b8SBarry Smith   a->nz = ai[m];
643ccb205f8SBarry Smith   for (i = 0; i < a->nz; i++) {
6446bdcaf15SBarry 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);
6459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(aa[i], MAT_FINAL_ASSEMBLY));
6469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(aa[i], MAT_FINAL_ASSEMBLY));
647ccb205f8SBarry Smith   }
6489566063dSJacob 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));
6499566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
6509566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax));
6512205254eSKarl Rupp 
6528e58a170SBarry Smith   A->info.mallocs += a->reallocs;
653421e10b8SBarry Smith   a->reallocs         = 0;
654421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
655421e10b8SBarry Smith   a->rmax             = rmax;
6569566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_BlockMat(A));
6573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
658421e10b8SBarry Smith }
659421e10b8SBarry Smith 
660d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetOption_BlockMat(Mat A, MatOption opt, PetscBool flg)
661d71ae5a4SJacob Faibussowitsch {
662290bbb0aSBarry Smith   PetscFunctionBegin;
6634e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
66441f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
665290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
666290bbb0aSBarry Smith   } else {
6679566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Unused matrix option %s\n", MatOptions[opt]));
668290bbb0aSBarry Smith   }
6693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
670290bbb0aSBarry Smith }
671290bbb0aSBarry Smith 
6723964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
673f4259b30SLisandro Dalcin                                        NULL,
674f4259b30SLisandro Dalcin                                        NULL,
675421e10b8SBarry Smith                                        MatMult_BlockMat,
676421e10b8SBarry Smith                                        /*  4*/ MatMultAdd_BlockMat,
677421e10b8SBarry Smith                                        MatMultTranspose_BlockMat,
678421e10b8SBarry Smith                                        MatMultTransposeAdd_BlockMat,
679f4259b30SLisandro Dalcin                                        NULL,
680f4259b30SLisandro Dalcin                                        NULL,
681f4259b30SLisandro Dalcin                                        NULL,
682f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
683f4259b30SLisandro Dalcin                                        NULL,
684f4259b30SLisandro Dalcin                                        NULL,
68541f059aeSBarry Smith                                        MatSOR_BlockMat,
686f4259b30SLisandro Dalcin                                        NULL,
687f4259b30SLisandro Dalcin                                        /* 15*/ NULL,
688f4259b30SLisandro Dalcin                                        NULL,
689f4259b30SLisandro Dalcin                                        NULL,
690f4259b30SLisandro Dalcin                                        NULL,
691f4259b30SLisandro Dalcin                                        NULL,
692f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
693421e10b8SBarry Smith                                        MatAssemblyEnd_BlockMat,
694290bbb0aSBarry Smith                                        MatSetOption_BlockMat,
695f4259b30SLisandro Dalcin                                        NULL,
696f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
697f4259b30SLisandro Dalcin                                        NULL,
698f4259b30SLisandro Dalcin                                        NULL,
699f4259b30SLisandro Dalcin                                        NULL,
700f4259b30SLisandro Dalcin                                        NULL,
701f4259b30SLisandro Dalcin                                        /* 29*/ NULL,
702f4259b30SLisandro Dalcin                                        NULL,
703f4259b30SLisandro Dalcin                                        NULL,
704f4259b30SLisandro Dalcin                                        NULL,
705f4259b30SLisandro Dalcin                                        NULL,
706f4259b30SLisandro Dalcin                                        /* 34*/ NULL,
707f4259b30SLisandro Dalcin                                        NULL,
708f4259b30SLisandro Dalcin                                        NULL,
709f4259b30SLisandro Dalcin                                        NULL,
710f4259b30SLisandro Dalcin                                        NULL,
711f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
712f4259b30SLisandro Dalcin                                        NULL,
713f4259b30SLisandro Dalcin                                        NULL,
714f4259b30SLisandro Dalcin                                        NULL,
715f4259b30SLisandro Dalcin                                        NULL,
716f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
717f4259b30SLisandro Dalcin                                        NULL,
7187d68702bSBarry Smith                                        MatShift_Basic,
719f4259b30SLisandro Dalcin                                        NULL,
720f4259b30SLisandro Dalcin                                        NULL,
721f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
722f4259b30SLisandro Dalcin                                        NULL,
723f4259b30SLisandro Dalcin                                        NULL,
724f4259b30SLisandro Dalcin                                        NULL,
725f4259b30SLisandro Dalcin                                        NULL,
726f4259b30SLisandro Dalcin                                        /* 54*/ NULL,
727f4259b30SLisandro Dalcin                                        NULL,
728f4259b30SLisandro Dalcin                                        NULL,
729f4259b30SLisandro Dalcin                                        NULL,
730f4259b30SLisandro Dalcin                                        NULL,
7317dae84e0SHong Zhang                                        /* 59*/ MatCreateSubMatrix_BlockMat,
732421e10b8SBarry Smith                                        MatDestroy_BlockMat,
733ccb205f8SBarry Smith                                        MatView_BlockMat,
734f4259b30SLisandro Dalcin                                        NULL,
735f4259b30SLisandro Dalcin                                        NULL,
736f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
737f4259b30SLisandro Dalcin                                        NULL,
738f4259b30SLisandro Dalcin                                        NULL,
739f4259b30SLisandro Dalcin                                        NULL,
740f4259b30SLisandro Dalcin                                        NULL,
741f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
742f4259b30SLisandro Dalcin                                        NULL,
743f4259b30SLisandro Dalcin                                        NULL,
744f4259b30SLisandro Dalcin                                        NULL,
745f4259b30SLisandro Dalcin                                        NULL,
746f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
747f4259b30SLisandro Dalcin                                        NULL,
748f4259b30SLisandro Dalcin                                        NULL,
749f4259b30SLisandro Dalcin                                        NULL,
750f4259b30SLisandro Dalcin                                        NULL,
751f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
752f4259b30SLisandro Dalcin                                        NULL,
753f4259b30SLisandro Dalcin                                        NULL,
754f4259b30SLisandro Dalcin                                        NULL,
7555bba2384SShri Abhyankar                                        MatLoad_BlockMat,
756f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
757f4259b30SLisandro Dalcin                                        NULL,
758f4259b30SLisandro Dalcin                                        NULL,
759f4259b30SLisandro Dalcin                                        NULL,
760f4259b30SLisandro Dalcin                                        NULL,
761f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
762f4259b30SLisandro Dalcin                                        NULL,
763f4259b30SLisandro Dalcin                                        NULL,
764f4259b30SLisandro Dalcin                                        NULL,
765f4259b30SLisandro Dalcin                                        NULL,
766f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
767f4259b30SLisandro Dalcin                                        NULL,
768f4259b30SLisandro Dalcin                                        NULL,
769f4259b30SLisandro Dalcin                                        NULL,
770f4259b30SLisandro Dalcin                                        NULL,
771f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
772f4259b30SLisandro Dalcin                                        NULL,
773f4259b30SLisandro Dalcin                                        NULL,
774f4259b30SLisandro Dalcin                                        NULL,
775f4259b30SLisandro Dalcin                                        NULL,
776f4259b30SLisandro Dalcin                                        /*104*/ NULL,
777f4259b30SLisandro Dalcin                                        NULL,
778f4259b30SLisandro Dalcin                                        NULL,
779f4259b30SLisandro Dalcin                                        NULL,
780f4259b30SLisandro Dalcin                                        NULL,
781f4259b30SLisandro Dalcin                                        /*109*/ NULL,
782f4259b30SLisandro Dalcin                                        NULL,
783f4259b30SLisandro Dalcin                                        NULL,
784f4259b30SLisandro Dalcin                                        NULL,
785f4259b30SLisandro Dalcin                                        NULL,
786f4259b30SLisandro Dalcin                                        /*114*/ NULL,
787f4259b30SLisandro Dalcin                                        NULL,
788f4259b30SLisandro Dalcin                                        NULL,
789f4259b30SLisandro Dalcin                                        NULL,
790f4259b30SLisandro Dalcin                                        NULL,
791f4259b30SLisandro Dalcin                                        /*119*/ NULL,
792f4259b30SLisandro Dalcin                                        NULL,
793f4259b30SLisandro Dalcin                                        NULL,
794f4259b30SLisandro Dalcin                                        NULL,
795f4259b30SLisandro Dalcin                                        NULL,
796f4259b30SLisandro Dalcin                                        /*124*/ NULL,
797f4259b30SLisandro Dalcin                                        NULL,
798f4259b30SLisandro Dalcin                                        NULL,
799f4259b30SLisandro Dalcin                                        NULL,
800f4259b30SLisandro Dalcin                                        NULL,
801f4259b30SLisandro Dalcin                                        /*129*/ NULL,
802f4259b30SLisandro Dalcin                                        NULL,
803f4259b30SLisandro Dalcin                                        NULL,
804f4259b30SLisandro Dalcin                                        NULL,
805f4259b30SLisandro Dalcin                                        NULL,
806f4259b30SLisandro Dalcin                                        /*134*/ NULL,
807f4259b30SLisandro Dalcin                                        NULL,
808f4259b30SLisandro Dalcin                                        NULL,
809f4259b30SLisandro Dalcin                                        NULL,
810f4259b30SLisandro Dalcin                                        NULL,
811f4259b30SLisandro Dalcin                                        /*139*/ NULL,
812f4259b30SLisandro Dalcin                                        NULL,
813d70f29a3SPierre Jolivet                                        NULL,
814d70f29a3SPierre Jolivet                                        NULL,
815d70f29a3SPierre Jolivet                                        NULL,
816d70f29a3SPierre Jolivet                                        /*144*/ NULL,
817d70f29a3SPierre Jolivet                                        NULL,
818d70f29a3SPierre Jolivet                                        NULL,
81999a7f59eSMark Adams                                        NULL,
82099a7f59eSMark Adams                                        NULL,
8217fb60732SBarry Smith                                        NULL,
822dec0b466SHong Zhang                                        /*150*/ NULL,
823eede4a3fSMark Adams                                        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