xref: /petsc/src/mat/impls/blockmat/seq/blockmat.c (revision 5b12497f4ab2fae16b3dbcff7d3c608aa116b75c)
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 
16*421480d9SBarry Smith MatGetDiagonalMarkers(BlockMat, A->rmap->bs)
17*421480d9SBarry Smith 
MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)18d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
19d71ae5a4SJacob Faibussowitsch {
20290bbb0aSBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat *)A->data;
21290bbb0aSBarry Smith   PetscScalar       *x;
22a2ea699eSBarry Smith   const Mat         *v;
23290bbb0aSBarry Smith   const PetscScalar *b;
24d0f46423SBarry Smith   PetscInt           n = A->cmap->n, i, mbs = n / A->rmap->bs, j, bs = A->rmap->bs;
25290bbb0aSBarry Smith   const PetscInt    *idx;
26290bbb0aSBarry Smith   IS                 row, col;
27290bbb0aSBarry Smith   MatFactorInfo      info;
286e63c7a1SBarry Smith   Vec                left = a->left, right = a->right, middle = a->middle;
29290bbb0aSBarry Smith   Mat               *diag;
30290bbb0aSBarry Smith 
31290bbb0aSBarry Smith   PetscFunctionBegin;
32290bbb0aSBarry Smith   its = its * lits;
3308401ef6SPierre 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);
34aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
3508401ef6SPierre Jolivet   PetscCheck(omega == 1.0, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for omega not equal to 1.0");
3628b400f6SJacob Faibussowitsch   PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for fshift");
377a46b595SBarry 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");
38290bbb0aSBarry Smith 
39290bbb0aSBarry Smith   if (!a->diags) {
40*421480d9SBarry Smith     const PetscInt *adiag;
41*421480d9SBarry Smith 
42*421480d9SBarry Smith     PetscCall(MatGetDiagonalMarkers_BlockMat(A, &adiag, NULL));
439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs, &a->diags));
449566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
45290bbb0aSBarry Smith     for (i = 0; i < mbs; i++) {
469566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND, &row, &col));
47*421480d9SBarry Smith       PetscCall(MatCholeskyFactorSymbolic(a->diags[i], a->a[adiag[i]], row, &info));
48*421480d9SBarry Smith       PetscCall(MatCholeskyFactorNumeric(a->diags[i], a->a[adiag[i]], &info));
499566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&row));
509566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&col));
51290bbb0aSBarry Smith     }
529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(bb, &a->workb));
53290bbb0aSBarry Smith   }
54290bbb0aSBarry Smith   diag = a->diags;
55290bbb0aSBarry Smith 
569566063dSJacob Faibussowitsch   PetscCall(VecSet(xx, 0.0));
579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
58dd8e379bSPierre Jolivet   /* copy right-hand side because it must be modified during iteration */
599566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, a->workb));
609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(a->workb, &b));
61290bbb0aSBarry Smith 
6241f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
63290bbb0aSBarry Smith   while (its--) {
64290bbb0aSBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
65290bbb0aSBarry Smith       for (i = 0; i < mbs; i++) {
666e63c7a1SBarry Smith         n   = a->i[i + 1] - a->i[i] - 1;
676e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
686e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
69290bbb0aSBarry Smith 
709566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
71290bbb0aSBarry Smith         for (j = 0; j < n; j++) {
729566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(right, x + idx[j] * bs));
739566063dSJacob Faibussowitsch           PetscCall(MatMultAdd(v[j], right, left, left));
749566063dSJacob Faibussowitsch           PetscCall(VecResetArray(right));
75290bbb0aSBarry Smith         }
769566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
779566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
789566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
79290bbb0aSBarry Smith 
809566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
819566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
826e63c7a1SBarry Smith 
83dd8e379bSPierre Jolivet         /* now adjust right-hand side, see MatSOR_SeqSBAIJ */
846e63c7a1SBarry Smith         for (j = 0; j < n; j++) {
859566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(v[j], right, left));
869566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(middle, b + idx[j] * bs));
879566063dSJacob Faibussowitsch           PetscCall(VecAXPY(middle, -1.0, left));
889566063dSJacob Faibussowitsch           PetscCall(VecResetArray(middle));
896e63c7a1SBarry Smith         }
909566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
91290bbb0aSBarry Smith       }
92290bbb0aSBarry Smith     }
93290bbb0aSBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
94290bbb0aSBarry Smith       for (i = mbs - 1; i >= 0; i--) {
956e63c7a1SBarry Smith         n   = a->i[i + 1] - a->i[i] - 1;
966e63c7a1SBarry Smith         idx = a->j + a->i[i] + 1;
976e63c7a1SBarry Smith         v   = a->a + a->i[i] + 1;
98290bbb0aSBarry Smith 
999566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
100290bbb0aSBarry Smith         for (j = 0; j < n; j++) {
1019566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(right, x + idx[j] * bs));
1029566063dSJacob Faibussowitsch           PetscCall(MatMultAdd(v[j], right, left, left));
1039566063dSJacob Faibussowitsch           PetscCall(VecResetArray(right));
104290bbb0aSBarry Smith         }
1059566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
1069566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
1079566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
108290bbb0aSBarry Smith 
1099566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
1109566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
1119566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
112290bbb0aSBarry Smith       }
113290bbb0aSBarry Smith     }
114290bbb0aSBarry Smith   }
1159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(a->workb, &b));
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118290bbb0aSBarry Smith }
119290bbb0aSBarry Smith 
MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)120d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSOR_BlockMat(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
121d71ae5a4SJacob Faibussowitsch {
122ccb205f8SBarry Smith   Mat_BlockMat      *a = (Mat_BlockMat *)A->data;
123ccb205f8SBarry Smith   PetscScalar       *x;
124a2ea699eSBarry Smith   const Mat         *v;
125ccb205f8SBarry Smith   const PetscScalar *b;
126d0f46423SBarry Smith   PetscInt           n = A->cmap->n, i, mbs = n / A->rmap->bs, j, bs = A->rmap->bs;
127ccb205f8SBarry Smith   const PetscInt    *idx;
128ccb205f8SBarry Smith   IS                 row, col;
129ccb205f8SBarry Smith   MatFactorInfo      info;
130ccb205f8SBarry Smith   Vec                left = a->left, right = a->right;
131ccb205f8SBarry Smith   Mat               *diag;
132ccb205f8SBarry Smith 
133ccb205f8SBarry Smith   PetscFunctionBegin;
134ccb205f8SBarry Smith   its = its * lits;
13508401ef6SPierre 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);
136aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
13708401ef6SPierre Jolivet   PetscCheck(omega == 1.0, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for omega not equal to 1.0");
13828b400f6SJacob Faibussowitsch   PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for fshift");
139ccb205f8SBarry Smith 
140ccb205f8SBarry Smith   if (!a->diags) {
141*421480d9SBarry Smith     const PetscInt *adiag;
142*421480d9SBarry Smith 
143*421480d9SBarry Smith     PetscCall(MatGetDiagonalMarkers_BlockMat(A, &adiag, NULL));
1449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mbs, &a->diags));
1459566063dSJacob Faibussowitsch     PetscCall(MatFactorInfoInitialize(&info));
146ccb205f8SBarry Smith     for (i = 0; i < mbs; i++) {
147*421480d9SBarry Smith       PetscCall(MatGetOrdering(a->a[adiag[i]], MATORDERINGND, &row, &col));
148*421480d9SBarry Smith       PetscCall(MatLUFactorSymbolic(a->diags[i], a->a[adiag[i]], row, col, &info));
149*421480d9SBarry Smith       PetscCall(MatLUFactorNumeric(a->diags[i], a->a[adiag[i]], &info));
1509566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&row));
1519566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&col));
152ccb205f8SBarry Smith     }
153ccb205f8SBarry Smith   }
154ccb205f8SBarry Smith   diag = a->diags;
155ccb205f8SBarry Smith 
1569566063dSJacob Faibussowitsch   PetscCall(VecSet(xx, 0.0));
1579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
159ccb205f8SBarry Smith 
16041f059aeSBarry Smith   /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
161ccb205f8SBarry Smith   while (its--) {
162ccb205f8SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
163ccb205f8SBarry Smith       for (i = 0; i < mbs; i++) {
164ccb205f8SBarry Smith         n   = a->i[i + 1] - a->i[i];
165ccb205f8SBarry Smith         idx = a->j + a->i[i];
166ccb205f8SBarry Smith         v   = a->a + a->i[i];
167ccb205f8SBarry Smith 
1689566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
169ccb205f8SBarry Smith         for (j = 0; j < n; j++) {
170ccb205f8SBarry Smith           if (idx[j] != i) {
1719566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(right, x + idx[j] * bs));
1729566063dSJacob Faibussowitsch             PetscCall(MatMultAdd(v[j], right, left, left));
1739566063dSJacob Faibussowitsch             PetscCall(VecResetArray(right));
174ccb205f8SBarry Smith           }
175ccb205f8SBarry Smith         }
1769566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
1779566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
1789566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
179ccb205f8SBarry Smith 
1809566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
1819566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
1829566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
183ccb205f8SBarry Smith       }
184ccb205f8SBarry Smith     }
185ccb205f8SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
186ccb205f8SBarry Smith       for (i = mbs - 1; i >= 0; i--) {
187ccb205f8SBarry Smith         n   = a->i[i + 1] - a->i[i];
188ccb205f8SBarry Smith         idx = a->j + a->i[i];
189ccb205f8SBarry Smith         v   = a->a + a->i[i];
190ccb205f8SBarry Smith 
1919566063dSJacob Faibussowitsch         PetscCall(VecSet(left, 0.0));
192ccb205f8SBarry Smith         for (j = 0; j < n; j++) {
193ccb205f8SBarry Smith           if (idx[j] != i) {
1949566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(right, x + idx[j] * bs));
1959566063dSJacob Faibussowitsch             PetscCall(MatMultAdd(v[j], right, left, left));
1969566063dSJacob Faibussowitsch             PetscCall(VecResetArray(right));
197ccb205f8SBarry Smith           }
198ccb205f8SBarry Smith         }
1999566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, b + i * bs));
2009566063dSJacob Faibussowitsch         PetscCall(VecAYPX(left, -1.0, right));
2019566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
202ccb205f8SBarry Smith 
2039566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(right, x + i * bs));
2049566063dSJacob Faibussowitsch         PetscCall(MatSolve(diag[i], left, right));
2059566063dSJacob Faibussowitsch         PetscCall(VecResetArray(right));
206ccb205f8SBarry Smith       }
207ccb205f8SBarry Smith     }
208ccb205f8SBarry Smith   }
2099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212ccb205f8SBarry Smith }
213ccb205f8SBarry Smith 
MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)214d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValues_BlockMat(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
215d71ae5a4SJacob Faibussowitsch {
216421e10b8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
217ea0465efSStefano Zampini   PetscInt     *rp, k, low, high, t, row, nrow, i, col, l, lastcol = -1;
218ea0465efSStefano Zampini   PetscInt     *ai = a->i, *ailen = a->ilen;
219d0f46423SBarry Smith   PetscInt     *aj = a->j, nonew = a->nonew, bs = A->rmap->bs, brow, bcol;
220421e10b8SBarry Smith   PetscInt      ridx, cidx;
221ace3abfcSBarry Smith   PetscBool     roworiented = a->roworiented;
222421e10b8SBarry Smith   MatScalar     value;
223421e10b8SBarry Smith   Mat          *ap, *aa = a->a;
224421e10b8SBarry Smith 
225421e10b8SBarry Smith   PetscFunctionBegin;
226421e10b8SBarry Smith   for (k = 0; k < m; k++) { /* loop over added rows */
227421e10b8SBarry Smith     row  = im[k];
228421e10b8SBarry Smith     brow = row / bs;
229421e10b8SBarry Smith     if (row < 0) continue;
2306bdcaf15SBarry 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);
231421e10b8SBarry Smith     rp   = aj + ai[brow];
232421e10b8SBarry Smith     ap   = aa + ai[brow];
233421e10b8SBarry Smith     nrow = ailen[brow];
234421e10b8SBarry Smith     low  = 0;
235421e10b8SBarry Smith     high = nrow;
236421e10b8SBarry Smith     for (l = 0; l < n; l++) { /* loop over added columns */
237421e10b8SBarry Smith       if (in[l] < 0) continue;
2386bdcaf15SBarry 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);
2399371c9d4SSatish Balay       col  = in[l];
2409371c9d4SSatish Balay       bcol = col / bs;
241b94d7dedSBarry Smith       if (A->symmetric == PETSC_BOOL3_TRUE && brow > bcol) continue;
2429371c9d4SSatish Balay       ridx = row % bs;
2439371c9d4SSatish Balay       cidx = col % bs;
2442205254eSKarl Rupp       if (roworiented) value = v[l + k * n];
2452205254eSKarl Rupp       else value = v[k + l * m];
2462205254eSKarl Rupp 
2472205254eSKarl Rupp       if (col <= lastcol) low = 0;
2482205254eSKarl Rupp       else high = nrow;
249421e10b8SBarry Smith       lastcol = col;
250421e10b8SBarry Smith       while (high - low > 7) {
251421e10b8SBarry Smith         t = (low + high) / 2;
252421e10b8SBarry Smith         if (rp[t] > bcol) high = t;
253421e10b8SBarry Smith         else low = t;
254421e10b8SBarry Smith       }
255421e10b8SBarry Smith       for (i = low; i < high; i++) {
256421e10b8SBarry Smith         if (rp[i] > bcol) break;
2572205254eSKarl Rupp         if (rp[i] == bcol) goto noinsert1;
258421e10b8SBarry Smith       }
259421e10b8SBarry Smith       if (nonew == 1) goto noinsert1;
260ea0465efSStefano Zampini       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the block matrix", row, col);
261ea0465efSStefano Zampini #if 0
262fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A, a->mbs, 1, nrow, brow, bcol, rmax, aa, ai, aj, rp, ap, imax, nonew, Mat);
2639371c9d4SSatish Balay       N = nrow++ - 1;
2649371c9d4SSatish Balay       high++;
265421e10b8SBarry Smith       /* shift up all the later entries in this row */
266421e10b8SBarry Smith       for (ii = N; ii >= i; ii--) {
267421e10b8SBarry Smith         rp[ii + 1] = rp[ii];
268421e10b8SBarry Smith         ap[ii + 1] = ap[ii];
269421e10b8SBarry Smith       }
270f4259b30SLisandro Dalcin       if (N >= i) ap[i] = NULL;
271421e10b8SBarry Smith       rp[i] = bcol;
272421e10b8SBarry Smith       a->nz++;
273ea0465efSStefano Zampini #endif
274421e10b8SBarry Smith     noinsert1:;
27548a46eb9SPierre Jolivet       if (!*(ap + i)) PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, NULL, ap + i));
2769566063dSJacob Faibussowitsch       PetscCall(MatSetValues(ap[i], 1, &ridx, 1, &cidx, &value, is));
277421e10b8SBarry Smith       low = i;
278421e10b8SBarry Smith     }
279421e10b8SBarry Smith     ailen[brow] = nrow;
280421e10b8SBarry Smith   }
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
282421e10b8SBarry Smith }
283421e10b8SBarry Smith 
MatLoad_BlockMat(Mat newmat,PetscViewer viewer)284d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
285d71ae5a4SJacob Faibussowitsch {
286a5973382SShri Abhyankar   Mat                tmpA;
287a5973382SShri Abhyankar   PetscInt           i, j, m, n, bs = 1, ncols, *lens, currentcol, mbs, **ii, *ilens, nextcol, *llens, cnt = 0;
288a5973382SShri Abhyankar   const PetscInt    *cols;
289a5973382SShri Abhyankar   const PetscScalar *values;
290ace3abfcSBarry Smith   PetscBool          flg = PETSC_FALSE, notdone;
291a5973382SShri Abhyankar   Mat_SeqAIJ        *a;
292a5973382SShri Abhyankar   Mat_BlockMat      *amat;
293a5973382SShri Abhyankar 
294a5973382SShri Abhyankar   PetscFunctionBegin;
295c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
2969566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(viewer));
2979566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &tmpA));
2989566063dSJacob Faibussowitsch   PetscCall(MatSetType(tmpA, MATSEQAIJ));
2999566063dSJacob Faibussowitsch   PetscCall(MatLoad_SeqAIJ(tmpA, viewer));
300a5973382SShri Abhyankar 
3019566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(tmpA, &m, &n));
302d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, NULL, "Options for loading BlockMat matrix 1", "Mat");
3039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-matload_block_size", "Set the blocksize used to store the matrix", "MatLoad", bs, &bs, NULL));
3049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-matload_symmetric", "Store the matrix as symmetric", "MatLoad", flg, &flg, NULL));
305d0609cedSBarry Smith   PetscOptionsEnd();
306a5973382SShri Abhyankar 
307a5973382SShri Abhyankar   /* Determine number of nonzero blocks for each block row */
308a5973382SShri Abhyankar   a   = (Mat_SeqAIJ *)tmpA->data;
309a5973382SShri Abhyankar   mbs = m / bs;
3109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(mbs, &lens, bs, &ii, bs, &ilens));
3119566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lens, mbs));
312a5973382SShri Abhyankar 
313a5973382SShri Abhyankar   for (i = 0; i < mbs; i++) {
314a5973382SShri Abhyankar     for (j = 0; j < bs; j++) {
315a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i * bs + j];
316a5973382SShri Abhyankar       ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j];
317a5973382SShri Abhyankar     }
318a5973382SShri Abhyankar 
319a5973382SShri Abhyankar     currentcol = -1;
320a5973382SShri Abhyankar     while (PETSC_TRUE) {
321a5973382SShri Abhyankar       notdone = PETSC_FALSE;
322a5973382SShri Abhyankar       nextcol = 1000000000;
323a5973382SShri Abhyankar       for (j = 0; j < bs; j++) {
324f4f49eeaSPierre Jolivet         while (ilens[j] > 0 && ii[j][0] / bs <= currentcol) {
325a5973382SShri Abhyankar           ii[j]++;
326a5973382SShri Abhyankar           ilens[j]--;
327a5973382SShri Abhyankar         }
328a5973382SShri Abhyankar         if (ilens[j] > 0) {
329a5973382SShri Abhyankar           notdone = PETSC_TRUE;
330a5973382SShri Abhyankar           nextcol = PetscMin(nextcol, ii[j][0] / bs);
331a5973382SShri Abhyankar         }
332a5973382SShri Abhyankar       }
333a5973382SShri Abhyankar       if (!notdone) break;
334a5973382SShri Abhyankar       if (!flg || (nextcol >= i)) lens[i]++;
335a5973382SShri Abhyankar       currentcol = nextcol;
336a5973382SShri Abhyankar     }
337a5973382SShri Abhyankar   }
338a5973382SShri Abhyankar 
33948a46eb9SPierre 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));
3409566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(newmat, bs, 0, lens));
3411baa6e33SBarry Smith   if (flg) PetscCall(MatSetOption(newmat, MAT_SYMMETRIC, PETSC_TRUE));
34257508eceSPierre Jolivet   amat = (Mat_BlockMat *)newmat->data;
343a5973382SShri Abhyankar 
344a5973382SShri Abhyankar   /* preallocate the submatrices */
3459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &llens));
346a5973382SShri Abhyankar   for (i = 0; i < mbs; i++) { /* loops for block rows */
347a5973382SShri Abhyankar     for (j = 0; j < bs; j++) {
348a5973382SShri Abhyankar       ii[j]    = a->j + a->i[i * bs + j];
349a5973382SShri Abhyankar       ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j];
350a5973382SShri Abhyankar     }
351a5973382SShri Abhyankar 
352a5973382SShri Abhyankar     currentcol = 1000000000;
353a5973382SShri Abhyankar     for (j = 0; j < bs; j++) { /* loop over rows in block finding first nonzero block */
354ad540459SPierre Jolivet       if (ilens[j] > 0) currentcol = PetscMin(currentcol, ii[j][0] / bs);
355a5973382SShri Abhyankar     }
356a5973382SShri Abhyankar 
357a5973382SShri Abhyankar     while (PETSC_TRUE) { /* loops over blocks in block row */
358a5973382SShri Abhyankar       notdone = PETSC_FALSE;
359a5973382SShri Abhyankar       nextcol = 1000000000;
3609566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(llens, bs));
361a5973382SShri Abhyankar       for (j = 0; j < bs; j++) {                              /* loop over rows in block */
362f4f49eeaSPierre Jolivet         while (ilens[j] > 0 && ii[j][0] / bs <= currentcol) { /* loop over columns in row */
363a5973382SShri Abhyankar           ii[j]++;
364a5973382SShri Abhyankar           ilens[j]--;
365a5973382SShri Abhyankar           llens[j]++;
366a5973382SShri Abhyankar         }
367a5973382SShri Abhyankar         if (ilens[j] > 0) {
368a5973382SShri Abhyankar           notdone = PETSC_TRUE;
369a5973382SShri Abhyankar           nextcol = PetscMin(nextcol, ii[j][0] / bs);
370a5973382SShri Abhyankar         }
371a5973382SShri Abhyankar       }
37208401ef6SPierre Jolivet       PetscCheck(cnt < amat->maxnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of blocks found greater than expected %" PetscInt_FMT, cnt);
373a5973382SShri Abhyankar       if (!flg || currentcol >= i) {
374a5973382SShri Abhyankar         amat->j[cnt] = currentcol;
3759566063dSJacob Faibussowitsch         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, llens, amat->a + cnt++));
376a5973382SShri Abhyankar       }
377a5973382SShri Abhyankar 
378a5973382SShri Abhyankar       if (!notdone) break;
379a5973382SShri Abhyankar       currentcol = nextcol;
380a5973382SShri Abhyankar     }
381a5973382SShri Abhyankar     amat->ilen[i] = lens[i];
382a5973382SShri Abhyankar   }
383a5973382SShri Abhyankar 
3849566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lens, ii, ilens));
3859566063dSJacob Faibussowitsch   PetscCall(PetscFree(llens));
386a5973382SShri Abhyankar 
387a5973382SShri Abhyankar   /* copy over the matrix, one row at a time */
388a5973382SShri Abhyankar   for (i = 0; i < m; i++) {
3899566063dSJacob Faibussowitsch     PetscCall(MatGetRow(tmpA, i, &ncols, &cols, &values));
3909566063dSJacob Faibussowitsch     PetscCall(MatSetValues(newmat, 1, &i, ncols, cols, values, INSERT_VALUES));
3919566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(tmpA, i, &ncols, &cols, &values));
392a5973382SShri Abhyankar   }
3939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(newmat, MAT_FINAL_ASSEMBLY));
3949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(newmat, MAT_FINAL_ASSEMBLY));
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
396a5973382SShri Abhyankar }
397a5973382SShri Abhyankar 
MatView_BlockMat(Mat A,PetscViewer viewer)398d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_BlockMat(Mat A, PetscViewer viewer)
399d71ae5a4SJacob Faibussowitsch {
40064075487SBarry Smith   Mat_BlockMat     *a = (Mat_BlockMat *)A->data;
401ccb205f8SBarry Smith   const char       *name;
402ccb205f8SBarry Smith   PetscViewerFormat format;
403ccb205f8SBarry Smith 
404ccb205f8SBarry Smith   PetscFunctionBegin;
4059566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)A, &name));
4069566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
407ccb205f8SBarry Smith   if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
4089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Nonzero block matrices = %" PetscInt_FMT " \n", a->nz));
409b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(PetscViewerASCIIPrintf(viewer, "Only upper triangular part of symmetric matrix is stored\n"));
410ccb205f8SBarry Smith   }
4113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
412ccb205f8SBarry Smith }
413421e10b8SBarry Smith 
MatDestroy_BlockMat(Mat mat)414d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_BlockMat(Mat mat)
415d71ae5a4SJacob Faibussowitsch {
416421e10b8SBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)mat->data;
417ccb205f8SBarry Smith   PetscInt      i;
418421e10b8SBarry Smith 
419421e10b8SBarry Smith   PetscFunctionBegin;
4209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->right));
4219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->left));
4229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->middle));
4239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bmat->workb));
424ccb205f8SBarry Smith   if (bmat->diags) {
42548a46eb9SPierre Jolivet     for (i = 0; i < mat->rmap->n / mat->rmap->bs; i++) PetscCall(MatDestroy(&bmat->diags[i]));
426ccb205f8SBarry Smith   }
427ccb205f8SBarry Smith   if (bmat->a) {
42848a46eb9SPierre Jolivet     for (i = 0; i < bmat->nz; i++) PetscCall(MatDestroy(&bmat->a[i]));
429ccb205f8SBarry Smith   }
4309566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(mat, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i));
4319566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
4323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
433421e10b8SBarry Smith }
434421e10b8SBarry Smith 
MatMult_BlockMat(Mat A,Vec x,Vec y)435d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_BlockMat(Mat A, Vec x, Vec y)
436d71ae5a4SJacob Faibussowitsch {
437421e10b8SBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
438421e10b8SBarry Smith   PetscScalar  *xx, *yy;
439d0f46423SBarry Smith   PetscInt     *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j;
440421e10b8SBarry Smith   Mat          *aa;
441421e10b8SBarry Smith 
442421e10b8SBarry Smith   PetscFunctionBegin;
443421e10b8SBarry Smith   /*
444421e10b8SBarry Smith      Standard CSR multiply except each entry is a Mat
445421e10b8SBarry Smith   */
4469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
447ccb205f8SBarry Smith 
4489566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
450421e10b8SBarry Smith   aj = bmat->j;
451421e10b8SBarry Smith   aa = bmat->a;
452421e10b8SBarry Smith   ii = bmat->i;
453421e10b8SBarry Smith   for (i = 0; i < m; i++) {
454421e10b8SBarry Smith     jrow = ii[i];
4559566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left, yy + bs * i));
456421e10b8SBarry Smith     n = ii[i + 1] - jrow;
457421e10b8SBarry Smith     for (j = 0; j < n; j++) {
4589566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow]));
4599566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4609566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
461421e10b8SBarry Smith       jrow++;
462421e10b8SBarry Smith     }
4639566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
464421e10b8SBarry Smith   }
4659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
4669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
4673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
468421e10b8SBarry Smith }
469421e10b8SBarry Smith 
MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)47066976f2fSJacob Faibussowitsch static PetscErrorCode MatMult_BlockMat_Symmetric(Mat A, Vec x, Vec y)
471d71ae5a4SJacob Faibussowitsch {
472290bbb0aSBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
473290bbb0aSBarry Smith   PetscScalar  *xx, *yy;
474d0f46423SBarry Smith   PetscInt     *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j;
475290bbb0aSBarry Smith   Mat          *aa;
476290bbb0aSBarry Smith 
477290bbb0aSBarry Smith   PetscFunctionBegin;
478290bbb0aSBarry Smith   /*
479290bbb0aSBarry Smith      Standard CSR multiply except each entry is a Mat
480290bbb0aSBarry Smith   */
4819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &xx));
482290bbb0aSBarry Smith 
4839566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
485290bbb0aSBarry Smith   aj = bmat->j;
486290bbb0aSBarry Smith   aa = bmat->a;
487290bbb0aSBarry Smith   ii = bmat->i;
488290bbb0aSBarry Smith   for (i = 0; i < m; i++) {
489290bbb0aSBarry Smith     jrow = ii[i];
490290bbb0aSBarry Smith     n    = ii[i + 1] - jrow;
4919566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->left, yy + bs * i));
4929566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bmat->middle, xx + bs * i));
493290bbb0aSBarry Smith     /* if we ALWAYS required a diagonal entry then could remove this if test */
494290bbb0aSBarry Smith     if (aj[jrow] == i) {
4959566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow]));
4969566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
4979566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
498290bbb0aSBarry Smith       jrow++;
499290bbb0aSBarry Smith       n--;
500290bbb0aSBarry Smith     }
501290bbb0aSBarry Smith     for (j = 0; j < n; j++) {
5029566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow])); /* upper triangular part */
5039566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left));
5049566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
505290bbb0aSBarry Smith 
5069566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(bmat->right, yy + bs * aj[jrow])); /* lower triangular part */
5079566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(aa[jrow], bmat->middle, bmat->right, bmat->right));
5089566063dSJacob Faibussowitsch       PetscCall(VecResetArray(bmat->right));
509290bbb0aSBarry Smith       jrow++;
510290bbb0aSBarry Smith     }
5119566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->left));
5129566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bmat->middle));
513290bbb0aSBarry Smith   }
5149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &xx));
5159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
517290bbb0aSBarry Smith }
518290bbb0aSBarry Smith 
MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)519d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_BlockMat(Mat A, Vec x, Vec y, Vec z)
520d71ae5a4SJacob Faibussowitsch {
521421e10b8SBarry Smith   PetscFunctionBegin;
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
523421e10b8SBarry Smith }
524421e10b8SBarry Smith 
MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)525d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_BlockMat(Mat A, Vec x, Vec y)
526d71ae5a4SJacob Faibussowitsch {
527421e10b8SBarry Smith   PetscFunctionBegin;
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
529421e10b8SBarry Smith }
530421e10b8SBarry Smith 
MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)531d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A, Vec x, Vec y, Vec z)
532d71ae5a4SJacob Faibussowitsch {
533421e10b8SBarry Smith   PetscFunctionBegin;
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535421e10b8SBarry Smith }
536421e10b8SBarry Smith 
MatCreateSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat * B)537d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A, IS isrow, IS iscol, MatReuse scall, Mat *B)
538d71ae5a4SJacob Faibussowitsch {
539ccb205f8SBarry Smith   Mat_BlockMat *a = (Mat_BlockMat *)A->data;
540ccb205f8SBarry Smith   Mat_SeqAIJ   *c;
541ccb205f8SBarry Smith   PetscInt      i, k, first, step, lensi, nrows, ncols;
542d2a212eaSBarry Smith   PetscInt     *j_new, *i_new, *aj = a->j, *ailen = a->ilen;
5431620fd73SBarry Smith   PetscScalar  *a_new;
544ccb205f8SBarry Smith   Mat           C, *aa = a->a;
545ace3abfcSBarry Smith   PetscBool     stride, equal;
546ccb205f8SBarry Smith 
547ccb205f8SBarry Smith   PetscFunctionBegin;
5489566063dSJacob Faibussowitsch   PetscCall(ISEqual(isrow, iscol, &equal));
54928b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for identical column and row indices");
5509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride));
55128b400f6SJacob Faibussowitsch   PetscCheck(stride, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for stride indices");
5529566063dSJacob Faibussowitsch   PetscCall(ISStrideGetInfo(iscol, &first, &step));
55308401ef6SPierre Jolivet   PetscCheck(step == A->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only select one entry from each block");
554ccb205f8SBarry Smith 
5559566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow, &nrows));
556ccb205f8SBarry Smith   ncols = nrows;
557ccb205f8SBarry Smith 
558ccb205f8SBarry Smith   /* create submatrix */
559ccb205f8SBarry Smith   if (scall == MAT_REUSE_MATRIX) {
560ccb205f8SBarry Smith     PetscInt n_cols, n_rows;
561ccb205f8SBarry Smith     C = *B;
5629566063dSJacob Faibussowitsch     PetscCall(MatGetSize(C, &n_rows, &n_cols));
563aed4548fSBarry Smith     PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size");
5649566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(C));
565ccb205f8SBarry Smith   } else {
5669566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
5679566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE));
568b94d7dedSBarry Smith     if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetType(C, MATSEQSBAIJ));
569b94d7dedSBarry Smith     else PetscCall(MatSetType(C, MATSEQAIJ));
5709566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(C, 0, ailen));
5719566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(C, 1, 0, ailen));
572ccb205f8SBarry Smith   }
573ccb205f8SBarry Smith   c = (Mat_SeqAIJ *)C->data;
574ccb205f8SBarry Smith 
575ccb205f8SBarry Smith   /* loop over rows inserting into submatrix */
576ccb205f8SBarry Smith   a_new = c->a;
577ccb205f8SBarry Smith   j_new = c->j;
578ccb205f8SBarry Smith   i_new = c->i;
579ccb205f8SBarry Smith 
580ccb205f8SBarry Smith   for (i = 0; i < nrows; i++) {
581ccb205f8SBarry Smith     lensi = ailen[i];
582ccb205f8SBarry Smith     for (k = 0; k < lensi; k++) {
583ccb205f8SBarry Smith       *j_new++ = *aj++;
5849566063dSJacob Faibussowitsch       PetscCall(MatGetValue(*aa++, first, first, a_new++));
585ccb205f8SBarry Smith     }
586ccb205f8SBarry Smith     i_new[i + 1] = i_new[i] + lensi;
587ccb205f8SBarry Smith     c->ilen[i]   = lensi;
588ccb205f8SBarry Smith   }
589ccb205f8SBarry Smith 
5909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
5919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
592ccb205f8SBarry Smith   *B = C;
5933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
594ccb205f8SBarry Smith }
595ccb205f8SBarry Smith 
MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)596d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A, MatAssemblyType mode)
597d71ae5a4SJacob Faibussowitsch {
598421e10b8SBarry Smith   Mat_BlockMat *a      = (Mat_BlockMat *)A->data;
599421e10b8SBarry Smith   PetscInt      fshift = 0, i, j, *ai = a->i, *aj = a->j, *imax = a->imax;
600ccb205f8SBarry Smith   PetscInt      m = a->mbs, *ip, N, *ailen = a->ilen, rmax = 0;
601421e10b8SBarry Smith   Mat          *aa = a->a, *ap;
602421e10b8SBarry Smith 
603421e10b8SBarry Smith   PetscFunctionBegin;
6043ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
605421e10b8SBarry Smith 
606421e10b8SBarry Smith   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
607421e10b8SBarry Smith   for (i = 1; i < m; i++) {
608421e10b8SBarry Smith     /* move each row back by the amount of empty slots (fshift) before it*/
609421e10b8SBarry Smith     fshift += imax[i - 1] - ailen[i - 1];
610421e10b8SBarry Smith     rmax = PetscMax(rmax, ailen[i]);
611421e10b8SBarry Smith     if (fshift) {
612421e10b8SBarry Smith       ip = aj + ai[i];
613421e10b8SBarry Smith       ap = aa + ai[i];
614421e10b8SBarry Smith       N  = ailen[i];
615421e10b8SBarry Smith       for (j = 0; j < N; j++) {
616421e10b8SBarry Smith         ip[j - fshift] = ip[j];
617421e10b8SBarry Smith         ap[j - fshift] = ap[j];
618421e10b8SBarry Smith       }
619421e10b8SBarry Smith     }
620421e10b8SBarry Smith     ai[i] = ai[i - 1] + ailen[i - 1];
621421e10b8SBarry Smith   }
622421e10b8SBarry Smith   if (m) {
623421e10b8SBarry Smith     fshift += imax[m - 1] - ailen[m - 1];
624421e10b8SBarry Smith     ai[m] = ai[m - 1] + ailen[m - 1];
625421e10b8SBarry Smith   }
626421e10b8SBarry Smith   /* reset ilen and imax for each row */
627ad540459SPierre Jolivet   for (i = 0; i < m; i++) ailen[i] = imax[i] = ai[i + 1] - ai[i];
628421e10b8SBarry Smith   a->nz = ai[m];
629ccb205f8SBarry Smith   for (i = 0; i < a->nz; i++) {
6306bdcaf15SBarry 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);
6319566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(aa[i], MAT_FINAL_ASSEMBLY));
6329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(aa[i], MAT_FINAL_ASSEMBLY));
633ccb205f8SBarry Smith   }
6349566063dSJacob 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));
6359566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
6369566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax));
6372205254eSKarl Rupp 
6388e58a170SBarry Smith   A->info.mallocs += a->reallocs;
639421e10b8SBarry Smith   a->reallocs         = 0;
640421e10b8SBarry Smith   A->info.nz_unneeded = (double)fshift;
641421e10b8SBarry Smith   a->rmax             = rmax;
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
643421e10b8SBarry Smith }
644421e10b8SBarry Smith 
MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg)645d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetOption_BlockMat(Mat A, MatOption opt, PetscBool flg)
646d71ae5a4SJacob Faibussowitsch {
647290bbb0aSBarry Smith   PetscFunctionBegin;
6484e0d8c25SBarry Smith   if (opt == MAT_SYMMETRIC && flg) {
64941f059aeSBarry Smith     A->ops->sor  = MatSOR_BlockMat_Symmetric;
650290bbb0aSBarry Smith     A->ops->mult = MatMult_BlockMat_Symmetric;
651290bbb0aSBarry Smith   } else {
6529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Unused matrix option %s\n", MatOptions[opt]));
653290bbb0aSBarry Smith   }
6543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
655290bbb0aSBarry Smith }
656290bbb0aSBarry Smith 
6573964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
658f4259b30SLisandro Dalcin                                        NULL,
659f4259b30SLisandro Dalcin                                        NULL,
660421e10b8SBarry Smith                                        MatMult_BlockMat,
661421e10b8SBarry Smith                                        /*  4*/ MatMultAdd_BlockMat,
662421e10b8SBarry Smith                                        MatMultTranspose_BlockMat,
663421e10b8SBarry Smith                                        MatMultTransposeAdd_BlockMat,
664f4259b30SLisandro Dalcin                                        NULL,
665f4259b30SLisandro Dalcin                                        NULL,
666f4259b30SLisandro Dalcin                                        NULL,
667f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
668f4259b30SLisandro Dalcin                                        NULL,
669f4259b30SLisandro Dalcin                                        NULL,
67041f059aeSBarry Smith                                        MatSOR_BlockMat,
671f4259b30SLisandro Dalcin                                        NULL,
672f4259b30SLisandro Dalcin                                        /* 15*/ NULL,
673f4259b30SLisandro Dalcin                                        NULL,
674f4259b30SLisandro Dalcin                                        NULL,
675f4259b30SLisandro Dalcin                                        NULL,
676f4259b30SLisandro Dalcin                                        NULL,
677f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
678421e10b8SBarry Smith                                        MatAssemblyEnd_BlockMat,
679290bbb0aSBarry Smith                                        MatSetOption_BlockMat,
680f4259b30SLisandro Dalcin                                        NULL,
681f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
682f4259b30SLisandro Dalcin                                        NULL,
683f4259b30SLisandro Dalcin                                        NULL,
684f4259b30SLisandro Dalcin                                        NULL,
685f4259b30SLisandro Dalcin                                        NULL,
686f4259b30SLisandro Dalcin                                        /* 29*/ NULL,
687f4259b30SLisandro Dalcin                                        NULL,
688f4259b30SLisandro Dalcin                                        NULL,
689f4259b30SLisandro Dalcin                                        NULL,
690f4259b30SLisandro Dalcin                                        NULL,
691f4259b30SLisandro Dalcin                                        /* 34*/ NULL,
692f4259b30SLisandro Dalcin                                        NULL,
693f4259b30SLisandro Dalcin                                        NULL,
694f4259b30SLisandro Dalcin                                        NULL,
695f4259b30SLisandro Dalcin                                        NULL,
696f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
697f4259b30SLisandro Dalcin                                        NULL,
698f4259b30SLisandro Dalcin                                        NULL,
699f4259b30SLisandro Dalcin                                        NULL,
700f4259b30SLisandro Dalcin                                        NULL,
701f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
702f4259b30SLisandro Dalcin                                        NULL,
7037d68702bSBarry Smith                                        MatShift_Basic,
704f4259b30SLisandro Dalcin                                        NULL,
705f4259b30SLisandro Dalcin                                        NULL,
706f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
707f4259b30SLisandro Dalcin                                        NULL,
708f4259b30SLisandro Dalcin                                        NULL,
709f4259b30SLisandro Dalcin                                        NULL,
710f4259b30SLisandro Dalcin                                        NULL,
711f4259b30SLisandro Dalcin                                        /* 54*/ NULL,
712f4259b30SLisandro Dalcin                                        NULL,
713f4259b30SLisandro Dalcin                                        NULL,
714f4259b30SLisandro Dalcin                                        NULL,
715f4259b30SLisandro Dalcin                                        NULL,
7167dae84e0SHong Zhang                                        /* 59*/ MatCreateSubMatrix_BlockMat,
717421e10b8SBarry Smith                                        MatDestroy_BlockMat,
718ccb205f8SBarry Smith                                        MatView_BlockMat,
719f4259b30SLisandro Dalcin                                        NULL,
720f4259b30SLisandro Dalcin                                        NULL,
721f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
722f4259b30SLisandro Dalcin                                        NULL,
723f4259b30SLisandro Dalcin                                        NULL,
724f4259b30SLisandro Dalcin                                        NULL,
725f4259b30SLisandro Dalcin                                        NULL,
726f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
727f4259b30SLisandro Dalcin                                        NULL,
728f4259b30SLisandro Dalcin                                        NULL,
729f4259b30SLisandro Dalcin                                        NULL,
730f4259b30SLisandro Dalcin                                        NULL,
731f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
732f4259b30SLisandro Dalcin                                        NULL,
733f4259b30SLisandro Dalcin                                        NULL,
734f4259b30SLisandro Dalcin                                        NULL,
7358bb0f5c6SPierre Jolivet                                        MatLoad_BlockMat,
736f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
737f4259b30SLisandro Dalcin                                        NULL,
738f4259b30SLisandro Dalcin                                        NULL,
739f4259b30SLisandro Dalcin                                        NULL,
7408bb0f5c6SPierre Jolivet                                        NULL,
741f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
742f4259b30SLisandro Dalcin                                        NULL,
743f4259b30SLisandro Dalcin                                        NULL,
744f4259b30SLisandro Dalcin                                        NULL,
745f4259b30SLisandro Dalcin                                        NULL,
746f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
747f4259b30SLisandro Dalcin                                        NULL,
748f4259b30SLisandro Dalcin                                        NULL,
749f4259b30SLisandro Dalcin                                        NULL,
750f4259b30SLisandro Dalcin                                        NULL,
751f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
752f4259b30SLisandro Dalcin                                        NULL,
753f4259b30SLisandro Dalcin                                        NULL,
754f4259b30SLisandro Dalcin                                        NULL,
755f4259b30SLisandro Dalcin                                        NULL,
756f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
757f4259b30SLisandro Dalcin                                        NULL,
758f4259b30SLisandro Dalcin                                        NULL,
759f4259b30SLisandro Dalcin                                        NULL,
760f4259b30SLisandro Dalcin                                        NULL,
761f4259b30SLisandro Dalcin                                        /*104*/ NULL,
762f4259b30SLisandro Dalcin                                        NULL,
763f4259b30SLisandro Dalcin                                        NULL,
764f4259b30SLisandro Dalcin                                        NULL,
765f4259b30SLisandro Dalcin                                        NULL,
766f4259b30SLisandro Dalcin                                        /*109*/ NULL,
767f4259b30SLisandro Dalcin                                        NULL,
768f4259b30SLisandro Dalcin                                        NULL,
769f4259b30SLisandro Dalcin                                        NULL,
770f4259b30SLisandro Dalcin                                        NULL,
771f4259b30SLisandro Dalcin                                        /*114*/ NULL,
772f4259b30SLisandro Dalcin                                        NULL,
773f4259b30SLisandro Dalcin                                        NULL,
774f4259b30SLisandro Dalcin                                        NULL,
775f4259b30SLisandro Dalcin                                        NULL,
776f4259b30SLisandro Dalcin                                        /*119*/ NULL,
777f4259b30SLisandro Dalcin                                        NULL,
778f4259b30SLisandro Dalcin                                        NULL,
779f4259b30SLisandro Dalcin                                        NULL,
780f4259b30SLisandro Dalcin                                        NULL,
781f4259b30SLisandro Dalcin                                        /*124*/ NULL,
782f4259b30SLisandro Dalcin                                        NULL,
783f4259b30SLisandro Dalcin                                        NULL,
784f4259b30SLisandro Dalcin                                        NULL,
785f4259b30SLisandro Dalcin                                        NULL,
786f4259b30SLisandro Dalcin                                        /*129*/ NULL,
787f4259b30SLisandro Dalcin                                        NULL,
788f4259b30SLisandro Dalcin                                        NULL,
789f4259b30SLisandro Dalcin                                        NULL,
790f4259b30SLisandro Dalcin                                        NULL,
791f4259b30SLisandro Dalcin                                        /*134*/ NULL,
792f4259b30SLisandro Dalcin                                        NULL,
793f4259b30SLisandro Dalcin                                        NULL,
794f4259b30SLisandro Dalcin                                        NULL,
795f4259b30SLisandro Dalcin                                        NULL,
796f4259b30SLisandro Dalcin                                        /*139*/ NULL,
797f4259b30SLisandro Dalcin                                        NULL,
798d70f29a3SPierre Jolivet                                        NULL,
79903db1824SAlex Lindsay                                        NULL,
800dec0b466SHong Zhang                                        NULL};
801421e10b8SBarry Smith 
8028cdf2d9bSBarry Smith /*@C
8038cdf2d9bSBarry Smith   MatBlockMatSetPreallocation - For good matrix assembly performance
8048cdf2d9bSBarry Smith   the user should preallocate the matrix storage by setting the parameter nz
8058cdf2d9bSBarry Smith   (or the array nnz).  By setting these parameters accurately, performance
8068cdf2d9bSBarry Smith   during matrix assembly can be increased by more than a factor of 50.
8078cdf2d9bSBarry Smith 
808d083f849SBarry Smith   Collective
8098cdf2d9bSBarry Smith 
8108cdf2d9bSBarry Smith   Input Parameters:
8118cdf2d9bSBarry Smith + B   - The matrix
8128cdf2d9bSBarry Smith . bs  - size of each block in matrix
8138cdf2d9bSBarry Smith . nz  - number of nonzeros per block row (same for all rows)
8148cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows
8152ef1f0ffSBarry Smith          (possibly different for each row) or `NULL`
8168cdf2d9bSBarry Smith 
8178cdf2d9bSBarry Smith   Level: intermediate
8188cdf2d9bSBarry Smith 
8192ef1f0ffSBarry Smith   Notes:
8202ef1f0ffSBarry Smith   If `nnz` is given then `nz` is ignored
8212ef1f0ffSBarry Smith 
8222ef1f0ffSBarry Smith   Specify the preallocated storage with either `nz` or `nnz` (not both).
8232ef1f0ffSBarry Smith   Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
8242ef1f0ffSBarry Smith   allocation.
8252ef1f0ffSBarry Smith 
8261cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateBlockMat()`, `MatSetValues()`
8278cdf2d9bSBarry Smith @*/
MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])828d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBlockMatSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[])
829d71ae5a4SJacob Faibussowitsch {
8308cdf2d9bSBarry Smith   PetscFunctionBegin;
831cac4c232SBarry Smith   PetscTryMethod(B, "MatBlockMatSetPreallocation_C", (Mat, PetscInt, PetscInt, const PetscInt[]), (B, bs, nz, nnz));
8323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8338cdf2d9bSBarry Smith }
8348cdf2d9bSBarry Smith 
MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt * nnz)835d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatBlockMatSetPreallocation_BlockMat(Mat A, PetscInt bs, PetscInt nz, PetscInt *nnz)
836d71ae5a4SJacob Faibussowitsch {
8378cdf2d9bSBarry Smith   Mat_BlockMat *bmat = (Mat_BlockMat *)A->data;
8388cdf2d9bSBarry Smith   PetscInt      i;
8398cdf2d9bSBarry Smith 
8408cdf2d9bSBarry Smith   PetscFunctionBegin;
8419566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(A->rmap, bs));
8429566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(A->cmap, bs));
8439566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
8449566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
8459566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(A->rmap, &bs));
84634ef9618SShri Abhyankar 
8478cdf2d9bSBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
84808401ef6SPierre Jolivet   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
8498cdf2d9bSBarry Smith   if (nnz) {
850d0f46423SBarry Smith     for (i = 0; i < A->rmap->n / bs; i++) {
85108401ef6SPierre 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]);
85208401ef6SPierre 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);
8538cdf2d9bSBarry Smith     }
8548cdf2d9bSBarry Smith   }
855d0f46423SBarry Smith   bmat->mbs = A->rmap->n / bs;
8568cdf2d9bSBarry Smith 
8579566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs, NULL, &bmat->right));
8589566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs, NULL, &bmat->middle));
8599566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, bs, &bmat->left));
8608cdf2d9bSBarry Smith 
861aa624791SPierre Jolivet   if (!bmat->imax) PetscCall(PetscMalloc2(A->rmap->n, &bmat->imax, A->rmap->n, &bmat->ilen));
862966bd95aSPierre Jolivet   PetscCheck(PetscLikely(nnz), PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently requires block row by row preallocation");
8638cdf2d9bSBarry Smith   nz = 0;
864d0f46423SBarry Smith   for (i = 0; i < A->rmap->n / A->rmap->bs; i++) {
8658cdf2d9bSBarry Smith     bmat->imax[i] = nnz[i];
8668cdf2d9bSBarry Smith     nz += nnz[i];
8678cdf2d9bSBarry Smith   }
8688cdf2d9bSBarry Smith 
8698cdf2d9bSBarry Smith   /* bmat->ilen will count nonzeros in each row so far. */
8709566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(bmat->ilen, bmat->mbs));
8718cdf2d9bSBarry Smith 
8728cdf2d9bSBarry Smith   /* allocate the matrix space */
8739566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(A, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i));
8749f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscScalar), (void **)&bmat->a));
8759f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscInt), (void **)&bmat->j));
8769f0612e4SBarry Smith   PetscCall(PetscShmgetAllocateArray(A->rmap->n + 1, sizeof(PetscInt), (void **)&bmat->i));
8778cdf2d9bSBarry Smith   bmat->free_a  = PETSC_TRUE;
8788cdf2d9bSBarry Smith   bmat->free_ij = PETSC_TRUE;
8798cdf2d9bSBarry Smith 
8809f0612e4SBarry Smith   bmat->i[0] = 0;
8819f0612e4SBarry Smith   for (i = 1; i < bmat->mbs + 1; i++) bmat->i[i] = bmat->i[i - 1] + bmat->imax[i - 1];
8829f0612e4SBarry Smith 
8838cdf2d9bSBarry Smith   bmat->nz            = 0;
8848cdf2d9bSBarry Smith   bmat->maxnz         = nz;
8858cdf2d9bSBarry Smith   A->info.nz_unneeded = (double)bmat->maxnz;
8869566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
8873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8888cdf2d9bSBarry Smith }
8898cdf2d9bSBarry Smith 
890421e10b8SBarry Smith /*MC
89111a5261eSBarry Smith    MATBLOCKMAT - A matrix that is defined by a set of `Mat`'s that represents a sparse block matrix
892421e10b8SBarry Smith                  consisting of (usually) sparse blocks.
893421e10b8SBarry Smith 
894421e10b8SBarry Smith   Level: advanced
895421e10b8SBarry Smith 
8961cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateBlockMat()`
897421e10b8SBarry Smith M*/
898421e10b8SBarry Smith 
MatCreate_BlockMat(Mat A)899d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A)
900d71ae5a4SJacob Faibussowitsch {
901421e10b8SBarry Smith   Mat_BlockMat *b;
902421e10b8SBarry Smith 
903421e10b8SBarry Smith   PetscFunctionBegin;
9044dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
905421e10b8SBarry Smith   A->data         = (void *)b;
906aea10558SJacob Faibussowitsch   A->ops[0]       = MatOps_Values;
907421e10b8SBarry Smith   A->assembled    = PETSC_TRUE;
908421e10b8SBarry Smith   A->preallocated = PETSC_FALSE;
9099566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATBLOCKMAT));
910290bbb0aSBarry Smith 
9119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatBlockMatSetPreallocation_C", MatBlockMatSetPreallocation_BlockMat));
9123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
913421e10b8SBarry Smith }
914421e10b8SBarry Smith 
915421e10b8SBarry Smith /*@C
91611a5261eSBarry Smith   MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential `Mat` object
917421e10b8SBarry Smith 
918d083f849SBarry Smith   Collective
919421e10b8SBarry Smith 
920421e10b8SBarry Smith   Input Parameters:
921421e10b8SBarry Smith + comm - MPI communicator
922421e10b8SBarry Smith . m    - number of rows
923421e10b8SBarry Smith . n    - number of columns
9248cdf2d9bSBarry Smith . bs   - size of each submatrix
92511a5261eSBarry Smith . nz   - expected maximum number of nonzero blocks in row (use `PETSC_DEFAULT` if not known)
9262ef1f0ffSBarry Smith - nnz  - expected number of nonzers per block row if known (use `NULL` otherwise)
927421e10b8SBarry Smith 
928421e10b8SBarry Smith   Output Parameter:
929421e10b8SBarry Smith . A - the matrix
930421e10b8SBarry Smith 
931421e10b8SBarry Smith   Level: intermediate
932421e10b8SBarry Smith 
93395452b02SPatrick Sanan   Notes:
93411a5261eSBarry Smith   Matrices of this type are nominally-sparse matrices in which each "entry" is a `Mat` object.  Each `Mat` must
935bbd3f336SJed Brown   have the same size and be sequential.  The local and global sizes must be compatible with this decomposition.
936bbd3f336SJed Brown 
93711a5261eSBarry Smith   For matrices containing parallel submatrices and variable block sizes, see `MATNEST`.
93811a5261eSBarry Smith 
939fe59aa6dSJacob Faibussowitsch   Developer Notes:
94011a5261eSBarry Smith   I don't like the name, it is not `MATNESTMAT`
941421e10b8SBarry Smith 
9421cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATBLOCKMAT`, `MatCreateNest()`
943421e10b8SBarry Smith @*/
MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt * nnz,Mat * A)944d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateBlockMat(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt bs, PetscInt nz, PetscInt *nnz, Mat *A)
945d71ae5a4SJacob Faibussowitsch {
946421e10b8SBarry Smith   PetscFunctionBegin;
9479566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
9489566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
9499566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATBLOCKMAT));
9509566063dSJacob Faibussowitsch   PetscCall(MatBlockMatSetPreallocation(*A, bs, nz, nnz));
9513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
952421e10b8SBarry Smith }
953