1421e10b8SBarry Smith 2421e10b8SBarry Smith /* 3421e10b8SBarry Smith This provides a matrix that consists of Mats 4421e10b8SBarry Smith */ 5421e10b8SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> /* use the common AIJ data-structure */ 8421e10b8SBarry Smith 9421e10b8SBarry Smith typedef struct { 10421e10b8SBarry Smith SEQAIJHEADER(Mat); 11421e10b8SBarry Smith SEQBAIJHEADER; 12ccb205f8SBarry Smith Mat *diags; 13421e10b8SBarry Smith 146e63c7a1SBarry Smith Vec left, right, middle, workb; /* dummy vectors to perform local parts of product */ 15421e10b8SBarry Smith } Mat_BlockMat; 16421e10b8SBarry Smith 179371c9d4SSatish Balay static PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) { 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)); 112290bbb0aSBarry Smith PetscFunctionReturn(0); 113290bbb0aSBarry Smith } 114290bbb0aSBarry Smith 1159371c9d4SSatish Balay static PetscErrorCode MatSOR_BlockMat(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) { 116ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat *)A->data; 117ccb205f8SBarry Smith PetscScalar *x; 118a2ea699eSBarry Smith const Mat *v; 119ccb205f8SBarry Smith const PetscScalar *b; 120d0f46423SBarry Smith PetscInt n = A->cmap->n, i, mbs = n / A->rmap->bs, j, bs = A->rmap->bs; 121ccb205f8SBarry Smith const PetscInt *idx; 122ccb205f8SBarry Smith IS row, col; 123ccb205f8SBarry Smith MatFactorInfo info; 124ccb205f8SBarry Smith Vec left = a->left, right = a->right; 125ccb205f8SBarry Smith Mat *diag; 126ccb205f8SBarry Smith 127ccb205f8SBarry Smith PetscFunctionBegin; 128ccb205f8SBarry Smith its = its * lits; 12908401ef6SPierre 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); 130aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat"); 13108401ef6SPierre Jolivet PetscCheck(omega == 1.0, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for omega not equal to 1.0"); 13228b400f6SJacob Faibussowitsch PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for fshift"); 133ccb205f8SBarry Smith 134ccb205f8SBarry Smith if (!a->diags) { 1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mbs, &a->diags)); 1369566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 137ccb205f8SBarry Smith for (i = 0; i < mbs; i++) { 1389566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(a->a[a->diag[i]], MATORDERINGND, &row, &col)); 1399566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(a->diags[i], a->a[a->diag[i]], row, col, &info)); 1409566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(a->diags[i], a->a[a->diag[i]], &info)); 1419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 143ccb205f8SBarry Smith } 144ccb205f8SBarry Smith } 145ccb205f8SBarry Smith diag = a->diags; 146ccb205f8SBarry Smith 1479566063dSJacob Faibussowitsch PetscCall(VecSet(xx, 0.0)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 150ccb205f8SBarry Smith 15141f059aeSBarry Smith /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */ 152ccb205f8SBarry Smith while (its--) { 153ccb205f8SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 154ccb205f8SBarry Smith for (i = 0; i < mbs; i++) { 155ccb205f8SBarry Smith n = a->i[i + 1] - a->i[i]; 156ccb205f8SBarry Smith idx = a->j + a->i[i]; 157ccb205f8SBarry Smith v = a->a + a->i[i]; 158ccb205f8SBarry Smith 1599566063dSJacob Faibussowitsch PetscCall(VecSet(left, 0.0)); 160ccb205f8SBarry Smith for (j = 0; j < n; j++) { 161ccb205f8SBarry Smith if (idx[j] != i) { 1629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right, x + idx[j] * bs)); 1639566063dSJacob Faibussowitsch PetscCall(MatMultAdd(v[j], right, left, left)); 1649566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 165ccb205f8SBarry Smith } 166ccb205f8SBarry Smith } 1679566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right, b + i * bs)); 1689566063dSJacob Faibussowitsch PetscCall(VecAYPX(left, -1.0, right)); 1699566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 170ccb205f8SBarry Smith 1719566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right, x + i * bs)); 1729566063dSJacob Faibussowitsch PetscCall(MatSolve(diag[i], left, right)); 1739566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 174ccb205f8SBarry Smith } 175ccb205f8SBarry Smith } 176ccb205f8SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 177ccb205f8SBarry Smith for (i = mbs - 1; i >= 0; i--) { 178ccb205f8SBarry Smith n = a->i[i + 1] - a->i[i]; 179ccb205f8SBarry Smith idx = a->j + a->i[i]; 180ccb205f8SBarry Smith v = a->a + a->i[i]; 181ccb205f8SBarry Smith 1829566063dSJacob Faibussowitsch PetscCall(VecSet(left, 0.0)); 183ccb205f8SBarry Smith for (j = 0; j < n; j++) { 184ccb205f8SBarry Smith if (idx[j] != i) { 1859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right, x + idx[j] * bs)); 1869566063dSJacob Faibussowitsch PetscCall(MatMultAdd(v[j], right, left, left)); 1879566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 188ccb205f8SBarry Smith } 189ccb205f8SBarry Smith } 1909566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right, b + i * bs)); 1919566063dSJacob Faibussowitsch PetscCall(VecAYPX(left, -1.0, right)); 1929566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 193ccb205f8SBarry Smith 1949566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(right, x + i * bs)); 1959566063dSJacob Faibussowitsch PetscCall(MatSolve(diag[i], left, right)); 1969566063dSJacob Faibussowitsch PetscCall(VecResetArray(right)); 197ccb205f8SBarry Smith } 198ccb205f8SBarry Smith } 199ccb205f8SBarry Smith } 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 202ccb205f8SBarry Smith PetscFunctionReturn(0); 203ccb205f8SBarry Smith } 204ccb205f8SBarry Smith 2059371c9d4SSatish Balay static PetscErrorCode MatSetValues_BlockMat(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is) { 206421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat *)A->data; 207421e10b8SBarry Smith PetscInt *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N, lastcol = -1; 208421e10b8SBarry Smith PetscInt *imax = a->imax, *ai = a->i, *ailen = a->ilen; 209d0f46423SBarry Smith PetscInt *aj = a->j, nonew = a->nonew, bs = A->rmap->bs, brow, bcol; 210421e10b8SBarry Smith PetscInt ridx, cidx; 211ace3abfcSBarry Smith PetscBool roworiented = a->roworiented; 212421e10b8SBarry Smith MatScalar value; 213421e10b8SBarry Smith Mat *ap, *aa = a->a; 214421e10b8SBarry Smith 215421e10b8SBarry Smith PetscFunctionBegin; 216421e10b8SBarry Smith for (k = 0; k < m; k++) { /* loop over added rows */ 217421e10b8SBarry Smith row = im[k]; 218421e10b8SBarry Smith brow = row / bs; 219421e10b8SBarry Smith if (row < 0) continue; 2206bdcaf15SBarry 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); 221421e10b8SBarry Smith rp = aj + ai[brow]; 222421e10b8SBarry Smith ap = aa + ai[brow]; 223421e10b8SBarry Smith rmax = imax[brow]; 224421e10b8SBarry Smith nrow = ailen[brow]; 225421e10b8SBarry Smith low = 0; 226421e10b8SBarry Smith high = nrow; 227421e10b8SBarry Smith for (l = 0; l < n; l++) { /* loop over added columns */ 228421e10b8SBarry Smith if (in[l] < 0) continue; 2296bdcaf15SBarry 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); 2309371c9d4SSatish Balay col = in[l]; 2319371c9d4SSatish Balay bcol = col / bs; 232b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE && brow > bcol) continue; 2339371c9d4SSatish Balay ridx = row % bs; 2349371c9d4SSatish Balay cidx = col % bs; 2352205254eSKarl Rupp if (roworiented) value = v[l + k * n]; 2362205254eSKarl Rupp else value = v[k + l * m]; 2372205254eSKarl Rupp 2382205254eSKarl Rupp if (col <= lastcol) low = 0; 2392205254eSKarl Rupp else high = nrow; 240421e10b8SBarry Smith lastcol = col; 241421e10b8SBarry Smith while (high - low > 7) { 242421e10b8SBarry Smith t = (low + high) / 2; 243421e10b8SBarry Smith if (rp[t] > bcol) high = t; 244421e10b8SBarry Smith else low = t; 245421e10b8SBarry Smith } 246421e10b8SBarry Smith for (i = low; i < high; i++) { 247421e10b8SBarry Smith if (rp[i] > bcol) break; 2482205254eSKarl Rupp if (rp[i] == bcol) goto noinsert1; 249421e10b8SBarry Smith } 250421e10b8SBarry Smith if (nonew == 1) goto noinsert1; 25108401ef6SPierre Jolivet PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col); 252fef13f97SBarry Smith MatSeqXAIJReallocateAIJ(A, a->mbs, 1, nrow, brow, bcol, rmax, aa, ai, aj, rp, ap, imax, nonew, Mat); 2539371c9d4SSatish Balay N = nrow++ - 1; 2549371c9d4SSatish Balay high++; 255421e10b8SBarry Smith /* shift up all the later entries in this row */ 256421e10b8SBarry Smith for (ii = N; ii >= i; ii--) { 257421e10b8SBarry Smith rp[ii + 1] = rp[ii]; 258421e10b8SBarry Smith ap[ii + 1] = ap[ii]; 259421e10b8SBarry Smith } 260f4259b30SLisandro Dalcin if (N >= i) ap[i] = NULL; 261421e10b8SBarry Smith rp[i] = bcol; 262421e10b8SBarry Smith a->nz++; 263e56f5c9eSBarry Smith A->nonzerostate++; 264421e10b8SBarry Smith noinsert1:; 265*48a46eb9SPierre Jolivet if (!*(ap + i)) PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, NULL, ap + i)); 2669566063dSJacob Faibussowitsch PetscCall(MatSetValues(ap[i], 1, &ridx, 1, &cidx, &value, is)); 267421e10b8SBarry Smith low = i; 268421e10b8SBarry Smith } 269421e10b8SBarry Smith ailen[brow] = nrow; 270421e10b8SBarry Smith } 271421e10b8SBarry Smith PetscFunctionReturn(0); 272421e10b8SBarry Smith } 273421e10b8SBarry Smith 2749371c9d4SSatish Balay static PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer) { 275a5973382SShri Abhyankar Mat tmpA; 276a5973382SShri Abhyankar PetscInt i, j, m, n, bs = 1, ncols, *lens, currentcol, mbs, **ii, *ilens, nextcol, *llens, cnt = 0; 277a5973382SShri Abhyankar const PetscInt *cols; 278a5973382SShri Abhyankar const PetscScalar *values; 279ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE, notdone; 280a5973382SShri Abhyankar Mat_SeqAIJ *a; 281a5973382SShri Abhyankar Mat_BlockMat *amat; 282a5973382SShri Abhyankar 283a5973382SShri Abhyankar PetscFunctionBegin; 284c98fd787SBarry Smith /* force binary viewer to load .info file if it has not yet done so */ 2859566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(viewer)); 2869566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &tmpA)); 2879566063dSJacob Faibussowitsch PetscCall(MatSetType(tmpA, MATSEQAIJ)); 2889566063dSJacob Faibussowitsch PetscCall(MatLoad_SeqAIJ(tmpA, viewer)); 289a5973382SShri Abhyankar 2909566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(tmpA, &m, &n)); 291d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, NULL, "Options for loading BlockMat matrix 1", "Mat"); 2929566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-matload_block_size", "Set the blocksize used to store the matrix", "MatLoad", bs, &bs, NULL)); 2939566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-matload_symmetric", "Store the matrix as symmetric", "MatLoad", flg, &flg, NULL)); 294d0609cedSBarry Smith PetscOptionsEnd(); 295a5973382SShri Abhyankar 296a5973382SShri Abhyankar /* Determine number of nonzero blocks for each block row */ 297a5973382SShri Abhyankar a = (Mat_SeqAIJ *)tmpA->data; 298a5973382SShri Abhyankar mbs = m / bs; 2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(mbs, &lens, bs, &ii, bs, &ilens)); 3009566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(lens, mbs)); 301a5973382SShri Abhyankar 302a5973382SShri Abhyankar for (i = 0; i < mbs; i++) { 303a5973382SShri Abhyankar for (j = 0; j < bs; j++) { 304a5973382SShri Abhyankar ii[j] = a->j + a->i[i * bs + j]; 305a5973382SShri Abhyankar ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j]; 306a5973382SShri Abhyankar } 307a5973382SShri Abhyankar 308a5973382SShri Abhyankar currentcol = -1; 309a5973382SShri Abhyankar while (PETSC_TRUE) { 310a5973382SShri Abhyankar notdone = PETSC_FALSE; 311a5973382SShri Abhyankar nextcol = 1000000000; 312a5973382SShri Abhyankar for (j = 0; j < bs; j++) { 313a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0] / bs <= currentcol)) { 314a5973382SShri Abhyankar ii[j]++; 315a5973382SShri Abhyankar ilens[j]--; 316a5973382SShri Abhyankar } 317a5973382SShri Abhyankar if (ilens[j] > 0) { 318a5973382SShri Abhyankar notdone = PETSC_TRUE; 319a5973382SShri Abhyankar nextcol = PetscMin(nextcol, ii[j][0] / bs); 320a5973382SShri Abhyankar } 321a5973382SShri Abhyankar } 322a5973382SShri Abhyankar if (!notdone) break; 323a5973382SShri Abhyankar if (!flg || (nextcol >= i)) lens[i]++; 324a5973382SShri Abhyankar currentcol = nextcol; 325a5973382SShri Abhyankar } 326a5973382SShri Abhyankar } 327a5973382SShri Abhyankar 328*48a46eb9SPierre 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)); 3299566063dSJacob Faibussowitsch PetscCall(MatBlockMatSetPreallocation(newmat, bs, 0, lens)); 3301baa6e33SBarry Smith if (flg) PetscCall(MatSetOption(newmat, MAT_SYMMETRIC, PETSC_TRUE)); 331a5973382SShri Abhyankar amat = (Mat_BlockMat *)(newmat)->data; 332a5973382SShri Abhyankar 333a5973382SShri Abhyankar /* preallocate the submatrices */ 3349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &llens)); 335a5973382SShri Abhyankar for (i = 0; i < mbs; i++) { /* loops for block rows */ 336a5973382SShri Abhyankar for (j = 0; j < bs; j++) { 337a5973382SShri Abhyankar ii[j] = a->j + a->i[i * bs + j]; 338a5973382SShri Abhyankar ilens[j] = a->i[i * bs + j + 1] - a->i[i * bs + j]; 339a5973382SShri Abhyankar } 340a5973382SShri Abhyankar 341a5973382SShri Abhyankar currentcol = 1000000000; 342a5973382SShri Abhyankar for (j = 0; j < bs; j++) { /* loop over rows in block finding first nonzero block */ 3439371c9d4SSatish Balay if (ilens[j] > 0) { currentcol = PetscMin(currentcol, ii[j][0] / bs); } 344a5973382SShri Abhyankar } 345a5973382SShri Abhyankar 346a5973382SShri Abhyankar while (PETSC_TRUE) { /* loops over blocks in block row */ 347a5973382SShri Abhyankar notdone = PETSC_FALSE; 348a5973382SShri Abhyankar nextcol = 1000000000; 3499566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(llens, bs)); 350a5973382SShri Abhyankar for (j = 0; j < bs; j++) { /* loop over rows in block */ 351a5973382SShri Abhyankar while ((ilens[j] > 0 && ii[j][0] / bs <= currentcol)) { /* loop over columns in row */ 352a5973382SShri Abhyankar ii[j]++; 353a5973382SShri Abhyankar ilens[j]--; 354a5973382SShri Abhyankar llens[j]++; 355a5973382SShri Abhyankar } 356a5973382SShri Abhyankar if (ilens[j] > 0) { 357a5973382SShri Abhyankar notdone = PETSC_TRUE; 358a5973382SShri Abhyankar nextcol = PetscMin(nextcol, ii[j][0] / bs); 359a5973382SShri Abhyankar } 360a5973382SShri Abhyankar } 36108401ef6SPierre Jolivet PetscCheck(cnt < amat->maxnz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of blocks found greater than expected %" PetscInt_FMT, cnt); 362a5973382SShri Abhyankar if (!flg || currentcol >= i) { 363a5973382SShri Abhyankar amat->j[cnt] = currentcol; 3649566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, bs, bs, 0, llens, amat->a + cnt++)); 365a5973382SShri Abhyankar } 366a5973382SShri Abhyankar 367a5973382SShri Abhyankar if (!notdone) break; 368a5973382SShri Abhyankar currentcol = nextcol; 369a5973382SShri Abhyankar } 370a5973382SShri Abhyankar amat->ilen[i] = lens[i]; 371a5973382SShri Abhyankar } 372a5973382SShri Abhyankar 3739566063dSJacob Faibussowitsch PetscCall(PetscFree3(lens, ii, ilens)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(llens)); 375a5973382SShri Abhyankar 376a5973382SShri Abhyankar /* copy over the matrix, one row at a time */ 377a5973382SShri Abhyankar for (i = 0; i < m; i++) { 3789566063dSJacob Faibussowitsch PetscCall(MatGetRow(tmpA, i, &ncols, &cols, &values)); 3799566063dSJacob Faibussowitsch PetscCall(MatSetValues(newmat, 1, &i, ncols, cols, values, INSERT_VALUES)); 3809566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(tmpA, i, &ncols, &cols, &values)); 381a5973382SShri Abhyankar } 3829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(newmat, MAT_FINAL_ASSEMBLY)); 3839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(newmat, MAT_FINAL_ASSEMBLY)); 384a5973382SShri Abhyankar PetscFunctionReturn(0); 385a5973382SShri Abhyankar } 386a5973382SShri Abhyankar 3879371c9d4SSatish Balay static PetscErrorCode MatView_BlockMat(Mat A, PetscViewer viewer) { 38864075487SBarry Smith Mat_BlockMat *a = (Mat_BlockMat *)A->data; 389ccb205f8SBarry Smith const char *name; 390ccb205f8SBarry Smith PetscViewerFormat format; 391ccb205f8SBarry Smith 392ccb205f8SBarry Smith PetscFunctionBegin; 3939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)A, &name)); 3949566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 395ccb205f8SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { 3969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Nonzero block matrices = %" PetscInt_FMT " \n", a->nz)); 397b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(PetscViewerASCIIPrintf(viewer, "Only upper triangular part of symmetric matrix is stored\n")); 398ccb205f8SBarry Smith } 399ccb205f8SBarry Smith PetscFunctionReturn(0); 400ccb205f8SBarry Smith } 401421e10b8SBarry Smith 4029371c9d4SSatish Balay static PetscErrorCode MatDestroy_BlockMat(Mat mat) { 403421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat *)mat->data; 404ccb205f8SBarry Smith PetscInt i; 405421e10b8SBarry Smith 406421e10b8SBarry Smith PetscFunctionBegin; 4079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->right)); 4089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->left)); 4099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->middle)); 4109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bmat->workb)); 411ccb205f8SBarry Smith if (bmat->diags) { 412*48a46eb9SPierre Jolivet for (i = 0; i < mat->rmap->n / mat->rmap->bs; i++) PetscCall(MatDestroy(&bmat->diags[i])); 413ccb205f8SBarry Smith } 414ccb205f8SBarry Smith if (bmat->a) { 415*48a46eb9SPierre Jolivet for (i = 0; i < bmat->nz; i++) PetscCall(MatDestroy(&bmat->a[i])); 416ccb205f8SBarry Smith } 4179566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(mat, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i)); 4189566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 419421e10b8SBarry Smith PetscFunctionReturn(0); 420421e10b8SBarry Smith } 421421e10b8SBarry Smith 4229371c9d4SSatish Balay static PetscErrorCode MatMult_BlockMat(Mat A, Vec x, Vec y) { 423421e10b8SBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat *)A->data; 424421e10b8SBarry Smith PetscScalar *xx, *yy; 425d0f46423SBarry Smith PetscInt *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j; 426421e10b8SBarry Smith Mat *aa; 427421e10b8SBarry Smith 428421e10b8SBarry Smith PetscFunctionBegin; 429421e10b8SBarry Smith /* 430421e10b8SBarry Smith Standard CSR multiply except each entry is a Mat 431421e10b8SBarry Smith */ 4329566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 433ccb205f8SBarry Smith 4349566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4359566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 436421e10b8SBarry Smith aj = bmat->j; 437421e10b8SBarry Smith aa = bmat->a; 438421e10b8SBarry Smith ii = bmat->i; 439421e10b8SBarry Smith for (i = 0; i < m; i++) { 440421e10b8SBarry Smith jrow = ii[i]; 4419566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->left, yy + bs * i)); 442421e10b8SBarry Smith n = ii[i + 1] - jrow; 443421e10b8SBarry Smith for (j = 0; j < n; j++) { 4449566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow])); 4459566063dSJacob Faibussowitsch PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left)); 4469566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 447421e10b8SBarry Smith jrow++; 448421e10b8SBarry Smith } 4499566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->left)); 450421e10b8SBarry Smith } 4519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xx)); 4529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 453421e10b8SBarry Smith PetscFunctionReturn(0); 454421e10b8SBarry Smith } 455421e10b8SBarry Smith 4569371c9d4SSatish Balay PetscErrorCode MatMult_BlockMat_Symmetric(Mat A, Vec x, Vec y) { 457290bbb0aSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat *)A->data; 458290bbb0aSBarry Smith PetscScalar *xx, *yy; 459d0f46423SBarry Smith PetscInt *aj, i, *ii, jrow, m = A->rmap->n / A->rmap->bs, bs = A->rmap->bs, n, j; 460290bbb0aSBarry Smith Mat *aa; 461290bbb0aSBarry Smith 462290bbb0aSBarry Smith PetscFunctionBegin; 463290bbb0aSBarry Smith /* 464290bbb0aSBarry Smith Standard CSR multiply except each entry is a Mat 465290bbb0aSBarry Smith */ 4669566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 467290bbb0aSBarry Smith 4689566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4699566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 470290bbb0aSBarry Smith aj = bmat->j; 471290bbb0aSBarry Smith aa = bmat->a; 472290bbb0aSBarry Smith ii = bmat->i; 473290bbb0aSBarry Smith for (i = 0; i < m; i++) { 474290bbb0aSBarry Smith jrow = ii[i]; 475290bbb0aSBarry Smith n = ii[i + 1] - jrow; 4769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->left, yy + bs * i)); 4779566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->middle, xx + bs * i)); 478290bbb0aSBarry Smith /* if we ALWAYS required a diagonal entry then could remove this if test */ 479290bbb0aSBarry Smith if (aj[jrow] == i) { 4809566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow])); 4819566063dSJacob Faibussowitsch PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left)); 4829566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 483290bbb0aSBarry Smith jrow++; 484290bbb0aSBarry Smith n--; 485290bbb0aSBarry Smith } 486290bbb0aSBarry Smith for (j = 0; j < n; j++) { 4879566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right, xx + bs * aj[jrow])); /* upper triangular part */ 4889566063dSJacob Faibussowitsch PetscCall(MatMultAdd(aa[jrow], bmat->right, bmat->left, bmat->left)); 4899566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 490290bbb0aSBarry Smith 4919566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bmat->right, yy + bs * aj[jrow])); /* lower triangular part */ 4929566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(aa[jrow], bmat->middle, bmat->right, bmat->right)); 4939566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->right)); 494290bbb0aSBarry Smith jrow++; 495290bbb0aSBarry Smith } 4969566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->left)); 4979566063dSJacob Faibussowitsch PetscCall(VecResetArray(bmat->middle)); 498290bbb0aSBarry Smith } 4999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xx)); 5009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 501290bbb0aSBarry Smith PetscFunctionReturn(0); 502290bbb0aSBarry Smith } 503290bbb0aSBarry Smith 5049371c9d4SSatish Balay static PetscErrorCode MatMultAdd_BlockMat(Mat A, Vec x, Vec y, Vec z) { 505421e10b8SBarry Smith PetscFunctionBegin; 506421e10b8SBarry Smith PetscFunctionReturn(0); 507421e10b8SBarry Smith } 508421e10b8SBarry Smith 5099371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_BlockMat(Mat A, Vec x, Vec y) { 510421e10b8SBarry Smith PetscFunctionBegin; 511421e10b8SBarry Smith PetscFunctionReturn(0); 512421e10b8SBarry Smith } 513421e10b8SBarry Smith 5149371c9d4SSatish Balay static PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A, Vec x, Vec y, Vec z) { 515421e10b8SBarry Smith PetscFunctionBegin; 516421e10b8SBarry Smith PetscFunctionReturn(0); 517421e10b8SBarry Smith } 518421e10b8SBarry Smith 519ccb205f8SBarry Smith /* 520ccb205f8SBarry Smith Adds diagonal pointers to sparse matrix structure. 521ccb205f8SBarry Smith */ 5229371c9d4SSatish Balay static PetscErrorCode MatMarkDiagonal_BlockMat(Mat A) { 523ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat *)A->data; 524d0f46423SBarry Smith PetscInt i, j, mbs = A->rmap->n / A->rmap->bs; 525ccb205f8SBarry Smith 526ccb205f8SBarry Smith PetscFunctionBegin; 527*48a46eb9SPierre Jolivet if (!a->diag) PetscCall(PetscMalloc1(mbs, &a->diag)); 528ccb205f8SBarry Smith for (i = 0; i < mbs; i++) { 529ccb205f8SBarry Smith a->diag[i] = a->i[i + 1]; 530ccb205f8SBarry Smith for (j = a->i[i]; j < a->i[i + 1]; j++) { 531ccb205f8SBarry Smith if (a->j[j] == i) { 532ccb205f8SBarry Smith a->diag[i] = j; 533ccb205f8SBarry Smith break; 534ccb205f8SBarry Smith } 535ccb205f8SBarry Smith } 536ccb205f8SBarry Smith } 537ccb205f8SBarry Smith PetscFunctionReturn(0); 538ccb205f8SBarry Smith } 539ccb205f8SBarry Smith 5409371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrix_BlockMat(Mat A, IS isrow, IS iscol, MatReuse scall, Mat *B) { 541ccb205f8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat *)A->data; 542ccb205f8SBarry Smith Mat_SeqAIJ *c; 543ccb205f8SBarry Smith PetscInt i, k, first, step, lensi, nrows, ncols; 544d2a212eaSBarry Smith PetscInt *j_new, *i_new, *aj = a->j, *ailen = a->ilen; 5451620fd73SBarry Smith PetscScalar *a_new; 546ccb205f8SBarry Smith Mat C, *aa = a->a; 547ace3abfcSBarry Smith PetscBool stride, equal; 548ccb205f8SBarry Smith 549ccb205f8SBarry Smith PetscFunctionBegin; 5509566063dSJacob Faibussowitsch PetscCall(ISEqual(isrow, iscol, &equal)); 55128b400f6SJacob Faibussowitsch PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for identical column and row indices"); 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISSTRIDE, &stride)); 55328b400f6SJacob Faibussowitsch PetscCheck(stride, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for stride indices"); 5549566063dSJacob Faibussowitsch PetscCall(ISStrideGetInfo(iscol, &first, &step)); 55508401ef6SPierre Jolivet PetscCheck(step == A->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only select one entry from each block"); 556ccb205f8SBarry Smith 5579566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow, &nrows)); 558ccb205f8SBarry Smith ncols = nrows; 559ccb205f8SBarry Smith 560ccb205f8SBarry Smith /* create submatrix */ 561ccb205f8SBarry Smith if (scall == MAT_REUSE_MATRIX) { 562ccb205f8SBarry Smith PetscInt n_cols, n_rows; 563ccb205f8SBarry Smith C = *B; 5649566063dSJacob Faibussowitsch PetscCall(MatGetSize(C, &n_rows, &n_cols)); 565aed4548fSBarry Smith PetscCheck(n_rows == nrows && n_cols == ncols, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reused submatrix wrong size"); 5669566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(C)); 567ccb205f8SBarry Smith } else { 5689566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 5699566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, nrows, ncols, PETSC_DETERMINE, PETSC_DETERMINE)); 570b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetType(C, MATSEQSBAIJ)); 571b94d7dedSBarry Smith else PetscCall(MatSetType(C, MATSEQAIJ)); 5729566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(C, 0, ailen)); 5739566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(C, 1, 0, ailen)); 574ccb205f8SBarry Smith } 575ccb205f8SBarry Smith c = (Mat_SeqAIJ *)C->data; 576ccb205f8SBarry Smith 577ccb205f8SBarry Smith /* loop over rows inserting into submatrix */ 578ccb205f8SBarry Smith a_new = c->a; 579ccb205f8SBarry Smith j_new = c->j; 580ccb205f8SBarry Smith i_new = c->i; 581ccb205f8SBarry Smith 582ccb205f8SBarry Smith for (i = 0; i < nrows; i++) { 583ccb205f8SBarry Smith lensi = ailen[i]; 584ccb205f8SBarry Smith for (k = 0; k < lensi; k++) { 585ccb205f8SBarry Smith *j_new++ = *aj++; 5869566063dSJacob Faibussowitsch PetscCall(MatGetValue(*aa++, first, first, a_new++)); 587ccb205f8SBarry Smith } 588ccb205f8SBarry Smith i_new[i + 1] = i_new[i] + lensi; 589ccb205f8SBarry Smith c->ilen[i] = lensi; 590ccb205f8SBarry Smith } 591ccb205f8SBarry Smith 5929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 5939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 594ccb205f8SBarry Smith *B = C; 595ccb205f8SBarry Smith PetscFunctionReturn(0); 596ccb205f8SBarry Smith } 597ccb205f8SBarry Smith 5989371c9d4SSatish Balay static PetscErrorCode MatAssemblyEnd_BlockMat(Mat A, MatAssemblyType mode) { 599421e10b8SBarry Smith Mat_BlockMat *a = (Mat_BlockMat *)A->data; 600421e10b8SBarry Smith PetscInt fshift = 0, i, j, *ai = a->i, *aj = a->j, *imax = a->imax; 601ccb205f8SBarry Smith PetscInt m = a->mbs, *ip, N, *ailen = a->ilen, rmax = 0; 602421e10b8SBarry Smith Mat *aa = a->a, *ap; 603421e10b8SBarry Smith 604421e10b8SBarry Smith PetscFunctionBegin; 605421e10b8SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 606421e10b8SBarry Smith 607421e10b8SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 608421e10b8SBarry Smith for (i = 1; i < m; i++) { 609421e10b8SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 610421e10b8SBarry Smith fshift += imax[i - 1] - ailen[i - 1]; 611421e10b8SBarry Smith rmax = PetscMax(rmax, ailen[i]); 612421e10b8SBarry Smith if (fshift) { 613421e10b8SBarry Smith ip = aj + ai[i]; 614421e10b8SBarry Smith ap = aa + ai[i]; 615421e10b8SBarry Smith N = ailen[i]; 616421e10b8SBarry Smith for (j = 0; j < N; j++) { 617421e10b8SBarry Smith ip[j - fshift] = ip[j]; 618421e10b8SBarry Smith ap[j - fshift] = ap[j]; 619421e10b8SBarry Smith } 620421e10b8SBarry Smith } 621421e10b8SBarry Smith ai[i] = ai[i - 1] + ailen[i - 1]; 622421e10b8SBarry Smith } 623421e10b8SBarry Smith if (m) { 624421e10b8SBarry Smith fshift += imax[m - 1] - ailen[m - 1]; 625421e10b8SBarry Smith ai[m] = ai[m - 1] + ailen[m - 1]; 626421e10b8SBarry Smith } 627421e10b8SBarry Smith /* reset ilen and imax for each row */ 6289371c9d4SSatish Balay for (i = 0; i < m; i++) { ailen[i] = imax[i] = ai[i + 1] - ai[i]; } 629421e10b8SBarry Smith a->nz = ai[m]; 630ccb205f8SBarry Smith for (i = 0; i < a->nz; i++) { 6316bdcaf15SBarry 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); 6329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(aa[i], MAT_FINAL_ASSEMBLY)); 6339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(aa[i], MAT_FINAL_ASSEMBLY)); 634ccb205f8SBarry Smith } 6359566063dSJacob 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)); 6369566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs)); 6379566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", rmax)); 6382205254eSKarl Rupp 6398e58a170SBarry Smith A->info.mallocs += a->reallocs; 640421e10b8SBarry Smith a->reallocs = 0; 641421e10b8SBarry Smith A->info.nz_unneeded = (double)fshift; 642421e10b8SBarry Smith a->rmax = rmax; 6439566063dSJacob Faibussowitsch PetscCall(MatMarkDiagonal_BlockMat(A)); 644421e10b8SBarry Smith PetscFunctionReturn(0); 645421e10b8SBarry Smith } 646421e10b8SBarry Smith 6479371c9d4SSatish Balay static PetscErrorCode MatSetOption_BlockMat(Mat A, MatOption opt, PetscBool flg) { 648290bbb0aSBarry Smith PetscFunctionBegin; 6494e0d8c25SBarry Smith if (opt == MAT_SYMMETRIC && flg) { 65041f059aeSBarry Smith A->ops->sor = MatSOR_BlockMat_Symmetric; 651290bbb0aSBarry Smith A->ops->mult = MatMult_BlockMat_Symmetric; 652290bbb0aSBarry Smith } else { 6539566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Unused matrix option %s\n", MatOptions[opt])); 654290bbb0aSBarry Smith } 655290bbb0aSBarry Smith PetscFunctionReturn(0); 656290bbb0aSBarry Smith } 657290bbb0aSBarry Smith 6583964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_BlockMat, 659f4259b30SLisandro Dalcin NULL, 660f4259b30SLisandro Dalcin NULL, 661421e10b8SBarry Smith MatMult_BlockMat, 662421e10b8SBarry Smith /* 4*/ MatMultAdd_BlockMat, 663421e10b8SBarry Smith MatMultTranspose_BlockMat, 664421e10b8SBarry Smith MatMultTransposeAdd_BlockMat, 665f4259b30SLisandro Dalcin NULL, 666f4259b30SLisandro Dalcin NULL, 667f4259b30SLisandro Dalcin NULL, 668f4259b30SLisandro Dalcin /* 10*/ NULL, 669f4259b30SLisandro Dalcin NULL, 670f4259b30SLisandro Dalcin NULL, 67141f059aeSBarry Smith MatSOR_BlockMat, 672f4259b30SLisandro Dalcin NULL, 673f4259b30SLisandro Dalcin /* 15*/ NULL, 674f4259b30SLisandro Dalcin NULL, 675f4259b30SLisandro Dalcin NULL, 676f4259b30SLisandro Dalcin NULL, 677f4259b30SLisandro Dalcin NULL, 678f4259b30SLisandro Dalcin /* 20*/ NULL, 679421e10b8SBarry Smith MatAssemblyEnd_BlockMat, 680290bbb0aSBarry Smith MatSetOption_BlockMat, 681f4259b30SLisandro Dalcin NULL, 682f4259b30SLisandro Dalcin /* 24*/ NULL, 683f4259b30SLisandro Dalcin NULL, 684f4259b30SLisandro Dalcin NULL, 685f4259b30SLisandro Dalcin NULL, 686f4259b30SLisandro Dalcin NULL, 687f4259b30SLisandro Dalcin /* 29*/ NULL, 688f4259b30SLisandro Dalcin NULL, 689f4259b30SLisandro Dalcin NULL, 690f4259b30SLisandro Dalcin NULL, 691f4259b30SLisandro Dalcin NULL, 692f4259b30SLisandro Dalcin /* 34*/ NULL, 693f4259b30SLisandro Dalcin NULL, 694f4259b30SLisandro Dalcin NULL, 695f4259b30SLisandro Dalcin NULL, 696f4259b30SLisandro Dalcin NULL, 697f4259b30SLisandro Dalcin /* 39*/ NULL, 698f4259b30SLisandro Dalcin NULL, 699f4259b30SLisandro Dalcin NULL, 700f4259b30SLisandro Dalcin NULL, 701f4259b30SLisandro Dalcin NULL, 702f4259b30SLisandro Dalcin /* 44*/ NULL, 703f4259b30SLisandro Dalcin NULL, 7047d68702bSBarry Smith MatShift_Basic, 705f4259b30SLisandro Dalcin NULL, 706f4259b30SLisandro Dalcin NULL, 707f4259b30SLisandro Dalcin /* 49*/ NULL, 708f4259b30SLisandro Dalcin NULL, 709f4259b30SLisandro Dalcin NULL, 710f4259b30SLisandro Dalcin NULL, 711f4259b30SLisandro Dalcin NULL, 712f4259b30SLisandro Dalcin /* 54*/ NULL, 713f4259b30SLisandro Dalcin NULL, 714f4259b30SLisandro Dalcin NULL, 715f4259b30SLisandro Dalcin NULL, 716f4259b30SLisandro Dalcin NULL, 7177dae84e0SHong Zhang /* 59*/ MatCreateSubMatrix_BlockMat, 718421e10b8SBarry Smith MatDestroy_BlockMat, 719ccb205f8SBarry Smith MatView_BlockMat, 720f4259b30SLisandro Dalcin NULL, 721f4259b30SLisandro Dalcin NULL, 722f4259b30SLisandro Dalcin /* 64*/ NULL, 723f4259b30SLisandro Dalcin NULL, 724f4259b30SLisandro Dalcin NULL, 725f4259b30SLisandro Dalcin NULL, 726f4259b30SLisandro Dalcin NULL, 727f4259b30SLisandro Dalcin /* 69*/ NULL, 728f4259b30SLisandro Dalcin NULL, 729f4259b30SLisandro Dalcin NULL, 730f4259b30SLisandro Dalcin NULL, 731f4259b30SLisandro Dalcin NULL, 732f4259b30SLisandro Dalcin /* 74*/ NULL, 733f4259b30SLisandro Dalcin NULL, 734f4259b30SLisandro Dalcin NULL, 735f4259b30SLisandro Dalcin NULL, 736f4259b30SLisandro Dalcin NULL, 737f4259b30SLisandro Dalcin /* 79*/ NULL, 738f4259b30SLisandro Dalcin NULL, 739f4259b30SLisandro Dalcin NULL, 740f4259b30SLisandro Dalcin NULL, 7415bba2384SShri Abhyankar MatLoad_BlockMat, 742f4259b30SLisandro Dalcin /* 84*/ NULL, 743f4259b30SLisandro Dalcin NULL, 744f4259b30SLisandro Dalcin NULL, 745f4259b30SLisandro Dalcin NULL, 746f4259b30SLisandro Dalcin NULL, 747f4259b30SLisandro Dalcin /* 89*/ NULL, 748f4259b30SLisandro Dalcin NULL, 749f4259b30SLisandro Dalcin NULL, 750f4259b30SLisandro Dalcin NULL, 751f4259b30SLisandro Dalcin NULL, 752f4259b30SLisandro Dalcin /* 94*/ NULL, 753f4259b30SLisandro Dalcin NULL, 754f4259b30SLisandro Dalcin NULL, 755f4259b30SLisandro Dalcin NULL, 756f4259b30SLisandro Dalcin NULL, 757f4259b30SLisandro Dalcin /* 99*/ NULL, 758f4259b30SLisandro Dalcin NULL, 759f4259b30SLisandro Dalcin NULL, 760f4259b30SLisandro Dalcin NULL, 761f4259b30SLisandro Dalcin NULL, 762f4259b30SLisandro Dalcin /*104*/ NULL, 763f4259b30SLisandro Dalcin NULL, 764f4259b30SLisandro Dalcin NULL, 765f4259b30SLisandro Dalcin NULL, 766f4259b30SLisandro Dalcin NULL, 767f4259b30SLisandro Dalcin /*109*/ NULL, 768f4259b30SLisandro Dalcin NULL, 769f4259b30SLisandro Dalcin NULL, 770f4259b30SLisandro Dalcin NULL, 771f4259b30SLisandro Dalcin NULL, 772f4259b30SLisandro Dalcin /*114*/ NULL, 773f4259b30SLisandro Dalcin NULL, 774f4259b30SLisandro Dalcin NULL, 775f4259b30SLisandro Dalcin NULL, 776f4259b30SLisandro Dalcin NULL, 777f4259b30SLisandro Dalcin /*119*/ NULL, 778f4259b30SLisandro Dalcin NULL, 779f4259b30SLisandro Dalcin NULL, 780f4259b30SLisandro Dalcin NULL, 781f4259b30SLisandro Dalcin NULL, 782f4259b30SLisandro Dalcin /*124*/ NULL, 783f4259b30SLisandro Dalcin NULL, 784f4259b30SLisandro Dalcin NULL, 785f4259b30SLisandro Dalcin NULL, 786f4259b30SLisandro Dalcin NULL, 787f4259b30SLisandro Dalcin /*129*/ NULL, 788f4259b30SLisandro Dalcin NULL, 789f4259b30SLisandro Dalcin NULL, 790f4259b30SLisandro Dalcin NULL, 791f4259b30SLisandro Dalcin NULL, 792f4259b30SLisandro Dalcin /*134*/ NULL, 793f4259b30SLisandro Dalcin NULL, 794f4259b30SLisandro Dalcin NULL, 795f4259b30SLisandro Dalcin NULL, 796f4259b30SLisandro Dalcin NULL, 797f4259b30SLisandro Dalcin /*139*/ NULL, 798f4259b30SLisandro Dalcin NULL, 799d70f29a3SPierre Jolivet NULL, 800d70f29a3SPierre Jolivet NULL, 801d70f29a3SPierre Jolivet NULL, 802d70f29a3SPierre Jolivet /*144*/ NULL, 803d70f29a3SPierre Jolivet NULL, 804d70f29a3SPierre Jolivet NULL, 80599a7f59eSMark Adams NULL, 80699a7f59eSMark Adams NULL, 8077fb60732SBarry Smith NULL, 8089371c9d4SSatish Balay /*150*/ NULL}; 809421e10b8SBarry Smith 8108cdf2d9bSBarry Smith /*@C 8118cdf2d9bSBarry Smith MatBlockMatSetPreallocation - For good matrix assembly performance 8128cdf2d9bSBarry Smith the user should preallocate the matrix storage by setting the parameter nz 8138cdf2d9bSBarry Smith (or the array nnz). By setting these parameters accurately, performance 8148cdf2d9bSBarry Smith during matrix assembly can be increased by more than a factor of 50. 8158cdf2d9bSBarry Smith 816d083f849SBarry Smith Collective 8178cdf2d9bSBarry Smith 8188cdf2d9bSBarry Smith Input Parameters: 8198cdf2d9bSBarry Smith + B - The matrix 8208cdf2d9bSBarry Smith . bs - size of each block in matrix 8218cdf2d9bSBarry Smith . nz - number of nonzeros per block row (same for all rows) 8228cdf2d9bSBarry Smith - nnz - array containing the number of nonzeros in the various block rows 8230298fd71SBarry Smith (possibly different for each row) or NULL 8248cdf2d9bSBarry Smith 8258cdf2d9bSBarry Smith Notes: 8268cdf2d9bSBarry Smith If nnz is given then nz is ignored 8278cdf2d9bSBarry Smith 8288cdf2d9bSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 8290298fd71SBarry Smith Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory 8308cdf2d9bSBarry Smith allocation. For large problems you MUST preallocate memory or you 8318cdf2d9bSBarry Smith will get TERRIBLE performance, see the users' manual chapter on matrices. 8328cdf2d9bSBarry Smith 8338cdf2d9bSBarry Smith Level: intermediate 8348cdf2d9bSBarry Smith 835db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateBlockMat()`, `MatSetValues()` 8368cdf2d9bSBarry Smith 8378cdf2d9bSBarry Smith @*/ 8389371c9d4SSatish Balay PetscErrorCode MatBlockMatSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[]) { 8398cdf2d9bSBarry Smith PetscFunctionBegin; 840cac4c232SBarry Smith PetscTryMethod(B, "MatBlockMatSetPreallocation_C", (Mat, PetscInt, PetscInt, const PetscInt[]), (B, bs, nz, nnz)); 8418cdf2d9bSBarry Smith PetscFunctionReturn(0); 8428cdf2d9bSBarry Smith } 8438cdf2d9bSBarry Smith 8449371c9d4SSatish Balay static PetscErrorCode MatBlockMatSetPreallocation_BlockMat(Mat A, PetscInt bs, PetscInt nz, PetscInt *nnz) { 8458cdf2d9bSBarry Smith Mat_BlockMat *bmat = (Mat_BlockMat *)A->data; 8468cdf2d9bSBarry Smith PetscInt i; 8478cdf2d9bSBarry Smith 8488cdf2d9bSBarry Smith PetscFunctionBegin; 8499566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->rmap, bs)); 8509566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->cmap, bs)); 8519566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 8529566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 8539566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetBlockSize(A->rmap, &bs)); 85434ef9618SShri Abhyankar 8558cdf2d9bSBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 85608401ef6SPierre Jolivet PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz); 8578cdf2d9bSBarry Smith if (nnz) { 858d0f46423SBarry Smith for (i = 0; i < A->rmap->n / bs; i++) { 85908401ef6SPierre 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]); 86008401ef6SPierre 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); 8618cdf2d9bSBarry Smith } 8628cdf2d9bSBarry Smith } 863d0f46423SBarry Smith bmat->mbs = A->rmap->n / bs; 8648cdf2d9bSBarry Smith 8659566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs, NULL, &bmat->right)); 8669566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, bs, NULL, &bmat->middle)); 8679566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, bs, &bmat->left)); 8688cdf2d9bSBarry Smith 8692ee49352SLisandro Dalcin if (!bmat->imax) { 8709566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->rmap->n, &bmat->imax, A->rmap->n, &bmat->ilen)); 8719566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)A, 2 * A->rmap->n * sizeof(PetscInt))); 8722ee49352SLisandro Dalcin } 873c0aa6a63SJacob Faibussowitsch if (PetscLikely(nnz)) { 8748cdf2d9bSBarry Smith nz = 0; 875d0f46423SBarry Smith for (i = 0; i < A->rmap->n / A->rmap->bs; i++) { 8768cdf2d9bSBarry Smith bmat->imax[i] = nnz[i]; 8778cdf2d9bSBarry Smith nz += nnz[i]; 8788cdf2d9bSBarry Smith } 879f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently requires block row by row preallocation"); 8808cdf2d9bSBarry Smith 8818cdf2d9bSBarry Smith /* bmat->ilen will count nonzeros in each row so far. */ 8829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(bmat->ilen, bmat->mbs)); 8838cdf2d9bSBarry Smith 8848cdf2d9bSBarry Smith /* allocate the matrix space */ 8859566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, (PetscScalar **)&bmat->a, &bmat->j, &bmat->i)); 8869566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nz, &bmat->a, nz, &bmat->j, A->rmap->n + 1, &bmat->i)); 8879566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)A, (A->rmap->n + 1) * sizeof(PetscInt) + nz * (sizeof(PetscScalar) + sizeof(PetscInt)))); 8888cdf2d9bSBarry Smith bmat->i[0] = 0; 8899371c9d4SSatish Balay for (i = 1; i < bmat->mbs + 1; i++) { bmat->i[i] = bmat->i[i - 1] + bmat->imax[i - 1]; } 8908cdf2d9bSBarry Smith bmat->singlemalloc = PETSC_TRUE; 8918cdf2d9bSBarry Smith bmat->free_a = PETSC_TRUE; 8928cdf2d9bSBarry Smith bmat->free_ij = PETSC_TRUE; 8938cdf2d9bSBarry Smith 8948cdf2d9bSBarry Smith bmat->nz = 0; 8958cdf2d9bSBarry Smith bmat->maxnz = nz; 8968cdf2d9bSBarry Smith A->info.nz_unneeded = (double)bmat->maxnz; 8979566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 8988cdf2d9bSBarry Smith PetscFunctionReturn(0); 8998cdf2d9bSBarry Smith } 9008cdf2d9bSBarry Smith 901421e10b8SBarry Smith /*MC 902421e10b8SBarry Smith MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix 903421e10b8SBarry Smith consisting of (usually) sparse blocks. 904421e10b8SBarry Smith 905421e10b8SBarry Smith Level: advanced 906421e10b8SBarry Smith 907db781477SPatrick Sanan .seealso: `MatCreateBlockMat()` 908421e10b8SBarry Smith 909421e10b8SBarry Smith M*/ 910421e10b8SBarry Smith 9119371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A) { 912421e10b8SBarry Smith Mat_BlockMat *b; 913421e10b8SBarry Smith 914421e10b8SBarry Smith PetscFunctionBegin; 9159566063dSJacob Faibussowitsch PetscCall(PetscNewLog(A, &b)); 916421e10b8SBarry Smith A->data = (void *)b; 9179566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A->ops, &MatOps_Values, sizeof(struct _MatOps))); 918421e10b8SBarry Smith 919421e10b8SBarry Smith A->assembled = PETSC_TRUE; 920421e10b8SBarry Smith A->preallocated = PETSC_FALSE; 9219566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATBLOCKMAT)); 922290bbb0aSBarry Smith 9239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatBlockMatSetPreallocation_C", MatBlockMatSetPreallocation_BlockMat)); 924421e10b8SBarry Smith PetscFunctionReturn(0); 925421e10b8SBarry Smith } 926421e10b8SBarry Smith 927421e10b8SBarry Smith /*@C 928bbd3f336SJed Brown MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object 929421e10b8SBarry Smith 930d083f849SBarry Smith Collective 931421e10b8SBarry Smith 932421e10b8SBarry Smith Input Parameters: 933421e10b8SBarry Smith + comm - MPI communicator 934421e10b8SBarry Smith . m - number of rows 935421e10b8SBarry Smith . n - number of columns 9368cdf2d9bSBarry Smith . bs - size of each submatrix 9378cdf2d9bSBarry Smith . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known) 9380298fd71SBarry Smith - nnz - expected number of nonzers per block row if known (use NULL otherwise) 939421e10b8SBarry Smith 940421e10b8SBarry Smith Output Parameter: 941421e10b8SBarry Smith . A - the matrix 942421e10b8SBarry Smith 943421e10b8SBarry Smith Level: intermediate 944421e10b8SBarry Smith 94595452b02SPatrick Sanan Notes: 94695452b02SPatrick Sanan Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object. Each Mat must 947bbd3f336SJed Brown have the same size and be sequential. The local and global sizes must be compatible with this decomposition. 948bbd3f336SJed Brown 949bbd3f336SJed Brown For matrices containing parallel submatrices and variable block sizes, see MATNEST. 950421e10b8SBarry Smith 951db781477SPatrick Sanan .seealso: `MATBLOCKMAT`, `MatCreateNest()` 952421e10b8SBarry Smith @*/ 9539371c9d4SSatish Balay PetscErrorCode MatCreateBlockMat(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt bs, PetscInt nz, PetscInt *nnz, Mat *A) { 954421e10b8SBarry Smith PetscFunctionBegin; 9559566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, A)); 9569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A, m, n, PETSC_DETERMINE, PETSC_DETERMINE)); 9579566063dSJacob Faibussowitsch PetscCall(MatSetType(*A, MATBLOCKMAT)); 9589566063dSJacob Faibussowitsch PetscCall(MatBlockMatSetPreallocation(*A, bs, nz, nnz)); 959421e10b8SBarry Smith PetscFunctionReturn(0); 960421e10b8SBarry Smith } 961