xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision e3c720941f06ea1a5bfc20f190bf798d81831d82)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay /*
3a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
449b5e25fSSatish Balay   matrix storage format.
549b5e25fSSatish Balay */
6c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
8c6db04a5SJed Brown #include <petscblaslapack.h>
949b5e25fSSatish Balay 
10c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/relax.h>
1170dcbbb9SBarry Smith #define USESHORT
12c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/relax.h>
1370dcbbb9SBarry Smith 
1426cec326SBarry Smith /* defines MatSetValues_Seq_Hash(), MatAssemblyEnd_Seq_Hash(), MatSetUp_Seq_Hash() */
1526cec326SBarry Smith #define TYPE SBAIJ
1626cec326SBarry Smith #define TYPE_SBAIJ
1726cec326SBarry Smith #define TYPE_BS
1826cec326SBarry Smith #include "../src/mat/impls/aij/seq/seqhashmatsetvalues.h"
1926cec326SBarry Smith #undef TYPE_BS
2026cec326SBarry Smith #define TYPE_BS _BS
2126cec326SBarry Smith #define TYPE_BS_ON
2226cec326SBarry Smith #include "../src/mat/impls/aij/seq/seqhashmatsetvalues.h"
2326cec326SBarry Smith #undef TYPE_BS
2426cec326SBarry Smith #undef TYPE_SBAIJ
2526cec326SBarry Smith #include "../src/mat/impls/aij/seq/seqhashmat.h"
2626cec326SBarry Smith #undef TYPE
2726cec326SBarry Smith #undef TYPE_BS_ON
2826cec326SBarry Smith 
296214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
30cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
316214f412SHong Zhang #endif
32d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
33d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_SBAIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
34d24d4204SJose E. Roman #endif
3528d58a37SPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat, MatType, MatReuse, Mat *);
36b5b17502SBarry Smith 
3749b5e25fSSatish Balay /*
3849b5e25fSSatish Balay      Checks for missing diagonals
3949b5e25fSSatish Balay */
40d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A, PetscBool *missing, PetscInt *dd)
41d71ae5a4SJacob Faibussowitsch {
42045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
437734d3b5SMatthew G. Knepley   PetscInt     *diag, *ii = a->i, i;
4449b5e25fSSatish Balay 
4549b5e25fSSatish Balay   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSBAIJ(A));
472af78befSBarry Smith   *missing = PETSC_FALSE;
487734d3b5SMatthew G. Knepley   if (A->rmap->n > 0 && !ii) {
49358d2f5dSShri Abhyankar     *missing = PETSC_TRUE;
50358d2f5dSShri Abhyankar     if (dd) *dd = 0;
519566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
52358d2f5dSShri Abhyankar   } else {
53358d2f5dSShri Abhyankar     diag = a->diag;
5449b5e25fSSatish Balay     for (i = 0; i < a->mbs; i++) {
557734d3b5SMatthew G. Knepley       if (diag[i] >= ii[i + 1]) {
562af78befSBarry Smith         *missing = PETSC_TRUE;
572af78befSBarry Smith         if (dd) *dd = i;
582af78befSBarry Smith         break;
592af78befSBarry Smith       }
6049b5e25fSSatish Balay     }
61358d2f5dSShri Abhyankar   }
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6349b5e25fSSatish Balay }
6449b5e25fSSatish Balay 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
66d71ae5a4SJacob Faibussowitsch {
67045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
6848dd3d27SHong Zhang   PetscInt      i, j;
6949b5e25fSSatish Balay 
7049b5e25fSSatish Balay   PetscFunctionBegin;
7109f38230SBarry Smith   if (!a->diag) {
729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(a->mbs, &a->diag));
73c760cd28SBarry Smith     a->free_diag = PETSC_TRUE;
7409f38230SBarry Smith   }
7548dd3d27SHong Zhang   for (i = 0; i < a->mbs; i++) {
7648dd3d27SHong Zhang     a->diag[i] = a->i[i + 1];
7748dd3d27SHong Zhang     for (j = a->i[i]; j < a->i[i + 1]; j++) {
7848dd3d27SHong Zhang       if (a->j[j] == i) {
7948dd3d27SHong Zhang         a->diag[i] = j;
8048dd3d27SHong Zhang         break;
8148dd3d27SHong Zhang       }
8248dd3d27SHong Zhang     }
8348dd3d27SHong Zhang   }
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8549b5e25fSSatish Balay }
8649b5e25fSSatish Balay 
87d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *nn, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done)
88d71ae5a4SJacob Faibussowitsch {
89a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
902462f5fdSStefano Zampini   PetscInt      i, j, n = a->mbs, nz = a->i[n], *tia, *tja, bs = A->rmap->bs, k, l, cnt;
912462f5fdSStefano Zampini   PetscInt    **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
9249b5e25fSSatish Balay 
9349b5e25fSSatish Balay   PetscFunctionBegin;
94d3e5a4abSHong Zhang   *nn = n;
953ba16761SJacob Faibussowitsch   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
962462f5fdSStefano Zampini   if (symmetric) {
979566063dSJacob Faibussowitsch     PetscCall(MatToSymmetricIJ_SeqAIJ(n, a->i, a->j, PETSC_FALSE, 0, 0, &tia, &tja));
982462f5fdSStefano Zampini     nz = tia[n];
992462f5fdSStefano Zampini   } else {
1009371c9d4SSatish Balay     tia = a->i;
1019371c9d4SSatish Balay     tja = a->j;
1022462f5fdSStefano Zampini   }
1032462f5fdSStefano Zampini 
1042462f5fdSStefano Zampini   if (!blockcompressed && bs > 1) {
1052462f5fdSStefano Zampini     (*nn) *= bs;
1068f7157efSSatish Balay     /* malloc & create the natural set of indices */
1079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((n + 1) * bs, ia));
1082462f5fdSStefano Zampini     if (n) {
1092462f5fdSStefano Zampini       (*ia)[0] = oshift;
110ad540459SPierre Jolivet       for (j = 1; j < bs; j++) (*ia)[j] = (tia[1] - tia[0]) * bs + (*ia)[j - 1];
1112462f5fdSStefano Zampini     }
1122462f5fdSStefano Zampini 
1132462f5fdSStefano Zampini     for (i = 1; i < n; i++) {
1142462f5fdSStefano Zampini       (*ia)[i * bs] = (tia[i] - tia[i - 1]) * bs + (*ia)[i * bs - 1];
115ad540459SPierre Jolivet       for (j = 1; j < bs; j++) (*ia)[i * bs + j] = (tia[i + 1] - tia[i]) * bs + (*ia)[i * bs + j - 1];
1162462f5fdSStefano Zampini     }
117ad540459SPierre Jolivet     if (n) (*ia)[n * bs] = (tia[n] - tia[n - 1]) * bs + (*ia)[n * bs - 1];
1182462f5fdSStefano Zampini 
1192462f5fdSStefano Zampini     if (inja) {
1209566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nz * bs * bs, ja));
1212462f5fdSStefano Zampini       cnt = 0;
1222462f5fdSStefano Zampini       for (i = 0; i < n; i++) {
1238f7157efSSatish Balay         for (j = 0; j < bs; j++) {
1242462f5fdSStefano Zampini           for (k = tia[i]; k < tia[i + 1]; k++) {
125ad540459SPierre Jolivet             for (l = 0; l < bs; l++) (*ja)[cnt++] = bs * tja[k] + l;
1268f7157efSSatish Balay           }
1278f7157efSSatish Balay         }
1288f7157efSSatish Balay       }
1298f7157efSSatish Balay     }
1302462f5fdSStefano Zampini 
1312462f5fdSStefano Zampini     if (symmetric) { /* deallocate memory allocated in MatToSymmetricIJ_SeqAIJ() */
1329566063dSJacob Faibussowitsch       PetscCall(PetscFree(tia));
1339566063dSJacob Faibussowitsch       PetscCall(PetscFree(tja));
1342462f5fdSStefano Zampini     }
1352462f5fdSStefano Zampini   } else if (oshift == 1) {
1362462f5fdSStefano Zampini     if (symmetric) {
1372462f5fdSStefano Zampini       nz = tia[A->rmap->n / bs];
1382462f5fdSStefano Zampini       /*  add 1 to i and j indices */
1392462f5fdSStefano Zampini       for (i = 0; i < A->rmap->n / bs + 1; i++) tia[i] = tia[i] + 1;
1402462f5fdSStefano Zampini       *ia = tia;
1412462f5fdSStefano Zampini       if (ja) {
1422462f5fdSStefano Zampini         for (i = 0; i < nz; i++) tja[i] = tja[i] + 1;
1432462f5fdSStefano Zampini         *ja = tja;
1442462f5fdSStefano Zampini       }
1452462f5fdSStefano Zampini     } else {
1462462f5fdSStefano Zampini       nz = a->i[A->rmap->n / bs];
1472462f5fdSStefano Zampini       /* malloc space and  add 1 to i and j indices */
1489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(A->rmap->n / bs + 1, ia));
1492462f5fdSStefano Zampini       for (i = 0; i < A->rmap->n / bs + 1; i++) (*ia)[i] = a->i[i] + 1;
1502462f5fdSStefano Zampini       if (ja) {
1519566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nz, ja));
1522462f5fdSStefano Zampini         for (i = 0; i < nz; i++) (*ja)[i] = a->j[i] + 1;
1532462f5fdSStefano Zampini       }
1542462f5fdSStefano Zampini     }
1552462f5fdSStefano Zampini   } else {
1562462f5fdSStefano Zampini     *ia = tia;
1572462f5fdSStefano Zampini     if (ja) *ja = tja;
158a6ece127SHong Zhang   }
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16049b5e25fSSatish Balay }
16149b5e25fSSatish Balay 
162d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
163d71ae5a4SJacob Faibussowitsch {
16449b5e25fSSatish Balay   PetscFunctionBegin;
1653ba16761SJacob Faibussowitsch   if (!ia) PetscFunctionReturn(PETSC_SUCCESS);
1662462f5fdSStefano Zampini   if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) {
1679566063dSJacob Faibussowitsch     PetscCall(PetscFree(*ia));
1689566063dSJacob Faibussowitsch     if (ja) PetscCall(PetscFree(*ja));
169a6ece127SHong Zhang   }
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17149b5e25fSSatish Balay }
17249b5e25fSSatish Balay 
173d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
174d71ae5a4SJacob Faibussowitsch {
17549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
17649b5e25fSSatish Balay 
17749b5e25fSSatish Balay   PetscFunctionBegin;
178b4e2f619SBarry Smith   if (A->hash_active) {
179b4e2f619SBarry Smith     PetscInt bs;
180*e3c72094SPierre Jolivet     A->ops[0] = a->cops;
181b4e2f619SBarry Smith     PetscCall(PetscHMapIJVDestroy(&a->ht));
182b4e2f619SBarry Smith     PetscCall(MatGetBlockSize(A, &bs));
183b4e2f619SBarry Smith     if (bs > 1) PetscCall(PetscHSetIJDestroy(&a->bht));
184b4e2f619SBarry Smith     PetscCall(PetscFree(a->dnz));
185b4e2f619SBarry Smith     PetscCall(PetscFree(a->bdnz));
186b4e2f619SBarry Smith     A->hash_active = PETSC_FALSE;
187b4e2f619SBarry Smith   }
188a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
1893ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->N, a->nz));
190a9f03627SSatish Balay #endif
1919566063dSJacob Faibussowitsch   PetscCall(MatSeqXAIJFreeAIJ(A, &a->a, &a->j, &a->i));
1929566063dSJacob Faibussowitsch   if (a->free_diag) PetscCall(PetscFree(a->diag));
1939566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->row));
1949566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->col));
1959566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->icol));
1969566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->idiag));
1979566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->inode.size));
1989566063dSJacob Faibussowitsch   if (a->free_imax_ilen) PetscCall(PetscFree2(a->imax, a->ilen));
1999566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->solve_work));
2009566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->sor_work));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->solves_work));
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->mult_work));
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->saved_values));
2049566063dSJacob Faibussowitsch   if (a->free_jshort) PetscCall(PetscFree(a->jshort));
2059566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->inew));
2069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&a->parent));
2079566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
208901853e0SKris Buschelman 
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
2102e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJGetArray_C", NULL));
2112e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJRestoreArray_C", NULL));
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
2139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
2149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJSetColumnIndices_C", NULL));
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_seqaij_C", NULL));
2169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_seqbaij_C", NULL));
2179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJSetPreallocation_C", NULL));
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSBAIJSetPreallocationCSR_C", NULL));
2196214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_elemental_C", NULL));
2216214f412SHong Zhang #endif
222d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsbaij_scalapack_C", NULL));
224d24d4204SJose E. Roman #endif
2252e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22749b5e25fSSatish Balay }
22849b5e25fSSatish Balay 
229d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetOption_SeqSBAIJ(Mat A, MatOption op, PetscBool flg)
230d71ae5a4SJacob Faibussowitsch {
231045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
232eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
233eb1ec7c1SStefano Zampini   PetscInt bs;
234eb1ec7c1SStefano Zampini #endif
23549b5e25fSSatish Balay 
23649b5e25fSSatish Balay   PetscFunctionBegin;
237eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
2389566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
239eb1ec7c1SStefano Zampini #endif
2404d9d31abSKris Buschelman   switch (op) {
241d71ae5a4SJacob Faibussowitsch   case MAT_ROW_ORIENTED:
242d71ae5a4SJacob Faibussowitsch     a->roworiented = flg;
243d71ae5a4SJacob Faibussowitsch     break;
244d71ae5a4SJacob Faibussowitsch   case MAT_KEEP_NONZERO_PATTERN:
245d71ae5a4SJacob Faibussowitsch     a->keepnonzeropattern = flg;
246d71ae5a4SJacob Faibussowitsch     break;
247d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATIONS:
248d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? 0 : 1);
249d71ae5a4SJacob Faibussowitsch     break;
250d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATION_ERR:
251d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -1 : 0);
252d71ae5a4SJacob Faibussowitsch     break;
253d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_ALLOCATION_ERR:
254d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -2 : 0);
255d71ae5a4SJacob Faibussowitsch     break;
256d71ae5a4SJacob Faibussowitsch   case MAT_UNUSED_NONZERO_LOCATION_ERR:
257d71ae5a4SJacob Faibussowitsch     a->nounused = (flg ? -1 : 0);
258d71ae5a4SJacob Faibussowitsch     break;
2598c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
2604d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
2614d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
262d71ae5a4SJacob Faibussowitsch   case MAT_SORTED_FULL:
263d71ae5a4SJacob Faibussowitsch     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
264d71ae5a4SJacob Faibussowitsch     break;
2659a4540c5SBarry Smith   case MAT_HERMITIAN:
266eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
267eb1ec7c1SStefano Zampini     if (flg) { /* disable transpose ops */
26808401ef6SPierre Jolivet       PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for Hermitian with block size greater than 1");
269eb1ec7c1SStefano Zampini       A->ops->multtranspose    = NULL;
270eb1ec7c1SStefano Zampini       A->ops->multtransposeadd = NULL;
271b94d7dedSBarry Smith       A->symmetric             = PETSC_BOOL3_FALSE;
272eb1ec7c1SStefano Zampini     }
2730f2140c7SStefano Zampini #endif
274eeffb40dSHong Zhang     break;
27577e54ba9SKris Buschelman   case MAT_SYMMETRIC:
276eb1ec7c1SStefano Zampini   case MAT_SPD:
277eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
278eb1ec7c1SStefano Zampini     if (flg) { /* An hermitian and symmetric matrix has zero imaginary part (restore back transpose ops) */
279eb1ec7c1SStefano Zampini       A->ops->multtranspose    = A->ops->mult;
280eb1ec7c1SStefano Zampini       A->ops->multtransposeadd = A->ops->multadd;
281eb1ec7c1SStefano Zampini     }
282eb1ec7c1SStefano Zampini #endif
283eb1ec7c1SStefano Zampini     break;
284eb1ec7c1SStefano Zampini     /* These options are handled directly by MatSetOption() */
28577e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
2869a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
287b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
288672ba085SHong Zhang   case MAT_STRUCTURE_ONLY:
289b94d7dedSBarry Smith   case MAT_SPD_ETERNAL:
2904dcd73b1SHong Zhang     /* These options are handled directly by MatSetOption() */
291290bbb0aSBarry Smith     break;
292d71ae5a4SJacob Faibussowitsch   case MAT_IGNORE_LOWER_TRIANGULAR:
293d71ae5a4SJacob Faibussowitsch     a->ignore_ltriangular = flg;
294d71ae5a4SJacob Faibussowitsch     break;
295d71ae5a4SJacob Faibussowitsch   case MAT_ERROR_LOWER_TRIANGULAR:
296d71ae5a4SJacob Faibussowitsch     a->ignore_ltriangular = flg;
297d71ae5a4SJacob Faibussowitsch     break;
298d71ae5a4SJacob Faibussowitsch   case MAT_GETROW_UPPERTRIANGULAR:
299d71ae5a4SJacob Faibussowitsch     a->getrow_utriangular = flg;
300d71ae5a4SJacob Faibussowitsch     break;
301d71ae5a4SJacob Faibussowitsch   case MAT_SUBMAT_SINGLEIS:
302d71ae5a4SJacob Faibussowitsch     break;
303d71ae5a4SJacob Faibussowitsch   default:
304d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
30549b5e25fSSatish Balay   }
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30749b5e25fSSatish Balay }
30849b5e25fSSatish Balay 
309d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRow_SeqSBAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
310d71ae5a4SJacob Faibussowitsch {
31149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
31249b5e25fSSatish Balay 
31349b5e25fSSatish Balay   PetscFunctionBegin;
31408401ef6SPierre Jolivet   PetscCheck(!A || a->getrow_utriangular, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");
31552768537SHong Zhang 
316f5edf698SHong Zhang   /* Get the upper triangular part of the row */
3179566063dSJacob Faibussowitsch   PetscCall(MatGetRow_SeqBAIJ_private(A, row, nz, idx, v, a->i, a->j, a->a));
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31949b5e25fSSatish Balay }
32049b5e25fSSatish Balay 
321d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
322d71ae5a4SJacob Faibussowitsch {
32349b5e25fSSatish Balay   PetscFunctionBegin;
324cb4a9cd9SHong Zhang   if (nz) *nz = 0;
3259566063dSJacob Faibussowitsch   if (idx) PetscCall(PetscFree(*idx));
3269566063dSJacob Faibussowitsch   if (v) PetscCall(PetscFree(*v));
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32849b5e25fSSatish Balay }
32949b5e25fSSatish Balay 
330d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
331d71ae5a4SJacob Faibussowitsch {
332f5edf698SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
333f5edf698SHong Zhang 
334f5edf698SHong Zhang   PetscFunctionBegin;
335f5edf698SHong Zhang   a->getrow_utriangular = PETSC_TRUE;
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337f5edf698SHong Zhang }
338a323099bSStefano Zampini 
339d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
340d71ae5a4SJacob Faibussowitsch {
341f5edf698SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
342f5edf698SHong Zhang 
343f5edf698SHong Zhang   PetscFunctionBegin;
344f5edf698SHong Zhang   a->getrow_utriangular = PETSC_FALSE;
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
346f5edf698SHong Zhang }
347f5edf698SHong Zhang 
348d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTranspose_SeqSBAIJ(Mat A, MatReuse reuse, Mat *B)
349d71ae5a4SJacob Faibussowitsch {
35049b5e25fSSatish Balay   PetscFunctionBegin;
3517fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
352cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
3539566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, B));
354cf37664fSBarry Smith   } else if (reuse == MAT_REUSE_MATRIX) {
3559566063dSJacob Faibussowitsch     PetscCall(MatCopy(A, *B, SAME_NONZERO_PATTERN));
356fc4dec0aSBarry Smith   }
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35849b5e25fSSatish Balay }
35949b5e25fSSatish Balay 
360d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A, PetscViewer viewer)
361d71ae5a4SJacob Faibussowitsch {
36249b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ *)A->data;
363d0f46423SBarry Smith   PetscInt          i, j, bs = A->rmap->bs, k, l, bs2 = a->bs2;
364f3ef73ceSBarry Smith   PetscViewerFormat format;
365121deb67SSatish Balay   PetscInt         *diag;
366b3a0534dSBarry Smith   const char       *matname;
36749b5e25fSSatish Balay 
36849b5e25fSSatish Balay   PetscFunctionBegin;
3699566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
370456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
3719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  block size is %" PetscInt_FMT "\n", bs));
372fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
373d2507d54SMatthew Knepley     Mat aij;
374ade3a672SBarry Smith 
375d5f3da31SBarry Smith     if (A->factortype && bs > 1) {
3769566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n"));
3773ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
37870d5e725SHong Zhang     }
3799566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &aij));
38023a3927dSBarry Smith     if (((PetscObject)A)->name) PetscCall(PetscObjectGetName((PetscObject)A, &matname));
38123a3927dSBarry Smith     if (((PetscObject)A)->name) PetscCall(PetscObjectSetName((PetscObject)aij, matname));
38223a3927dSBarry Smith     PetscCall(MatView_SeqAIJ(aij, viewer));
3839566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&aij));
384fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
385b3a0534dSBarry Smith     Mat B;
386b3a0534dSBarry Smith 
387b3a0534dSBarry Smith     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &B));
388b3a0534dSBarry Smith     if (((PetscObject)A)->name) PetscCall(PetscObjectGetName((PetscObject)A, &matname));
389b3a0534dSBarry Smith     if (((PetscObject)A)->name) PetscCall(PetscObjectSetName((PetscObject)B, matname));
390b3a0534dSBarry Smith     PetscCall(MatView_SeqAIJ(B, viewer));
391b3a0534dSBarry Smith     PetscCall(MatDestroy(&B));
392c1490034SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
3933ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
39449b5e25fSSatish Balay   } else {
3959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
3962c990fa1SHong Zhang     if (A->factortype) { /* for factored matrix */
39708401ef6SPierre Jolivet       PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "matrix is factored with bs>1. Not implemented yet");
3982c990fa1SHong Zhang 
399121deb67SSatish Balay       diag = a->diag;
400121deb67SSatish Balay       for (i = 0; i < a->mbs; i++) { /* for row block i */
4019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
4022c990fa1SHong Zhang         /* diagonal entry */
4032c990fa1SHong Zhang #if defined(PETSC_USE_COMPLEX)
4042c990fa1SHong Zhang         if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
4059566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i) ", a->j[diag[i]], (double)PetscRealPart(1.0 / a->a[diag[i]]), (double)PetscImaginaryPart(1.0 / a->a[diag[i]])));
4062c990fa1SHong Zhang         } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
4079566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i) ", a->j[diag[i]], (double)PetscRealPart(1.0 / a->a[diag[i]]), -(double)PetscImaginaryPart(1.0 / a->a[diag[i]])));
4082c990fa1SHong Zhang         } else {
4099566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[diag[i]], (double)PetscRealPart(1.0 / a->a[diag[i]])));
4102c990fa1SHong Zhang         }
4112c990fa1SHong Zhang #else
4129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[diag[i]], (double)(1.0 / a->a[diag[i]])));
4132c990fa1SHong Zhang #endif
4142c990fa1SHong Zhang         /* off-diagonal entries */
4152c990fa1SHong Zhang         for (k = a->i[i]; k < a->i[i + 1] - 1; k++) {
4162c990fa1SHong Zhang #if defined(PETSC_USE_COMPLEX)
417ca0704adSBarry Smith           if (PetscImaginaryPart(a->a[k]) > 0.0) {
4189566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i) ", bs * a->j[k], (double)PetscRealPart(a->a[k]), (double)PetscImaginaryPart(a->a[k])));
419ca0704adSBarry Smith           } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
4209566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i) ", bs * a->j[k], (double)PetscRealPart(a->a[k]), -(double)PetscImaginaryPart(a->a[k])));
4212c990fa1SHong Zhang           } else {
4229566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", bs * a->j[k], (double)PetscRealPart(a->a[k])));
4232c990fa1SHong Zhang           }
4242c990fa1SHong Zhang #else
4259566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->j[k], (double)a->a[k]));
4262c990fa1SHong Zhang #endif
4272c990fa1SHong Zhang         }
4289566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
4292c990fa1SHong Zhang       }
4302c990fa1SHong Zhang 
4312c990fa1SHong Zhang     } else {                         /* for non-factored matrix */
4320c74a584SJed Brown       for (i = 0; i < a->mbs; i++) { /* for row block i */
4330c74a584SJed Brown         for (j = 0; j < bs; j++) {   /* for row bs*i + j */
4349566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i * bs + j));
4350c74a584SJed Brown           for (k = a->i[i]; k < a->i[i + 1]; k++) { /* for column block */
4360c74a584SJed Brown             for (l = 0; l < bs; l++) {              /* for column */
43749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
43849b5e25fSSatish Balay               if (PetscImaginaryPart(a->a[bs2 * k + l * bs + j]) > 0.0) {
4399371c9d4SSatish Balay                 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i) ", bs * a->j[k] + l, (double)PetscRealPart(a->a[bs2 * k + l * bs + j]), (double)PetscImaginaryPart(a->a[bs2 * k + l * bs + j])));
44049b5e25fSSatish Balay               } else if (PetscImaginaryPart(a->a[bs2 * k + l * bs + j]) < 0.0) {
4419371c9d4SSatish Balay                 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i) ", bs * a->j[k] + l, (double)PetscRealPart(a->a[bs2 * k + l * bs + j]), -(double)PetscImaginaryPart(a->a[bs2 * k + l * bs + j])));
44249b5e25fSSatish Balay               } else {
4439566063dSJacob Faibussowitsch                 PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", bs * a->j[k] + l, (double)PetscRealPart(a->a[bs2 * k + l * bs + j])));
44449b5e25fSSatish Balay               }
44549b5e25fSSatish Balay #else
4469566063dSJacob Faibussowitsch               PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", bs * a->j[k] + l, (double)a->a[bs2 * k + l * bs + j]));
44749b5e25fSSatish Balay #endif
44849b5e25fSSatish Balay             }
44949b5e25fSSatish Balay           }
4509566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
45149b5e25fSSatish Balay         }
45249b5e25fSSatish Balay       }
4532c990fa1SHong Zhang     }
4549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
45549b5e25fSSatish Balay   }
4569566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45849b5e25fSSatish Balay }
45949b5e25fSSatish Balay 
4609804daf3SBarry Smith #include <petscdraw.h>
461d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw, void *Aa)
462d71ae5a4SJacob Faibussowitsch {
46349b5e25fSSatish Balay   Mat           A = (Mat)Aa;
46449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
465d0f46423SBarry Smith   PetscInt      row, i, j, k, l, mbs = a->mbs, color, bs = A->rmap->bs, bs2 = a->bs2;
46649b5e25fSSatish Balay   PetscReal     xl, yl, xr, yr, x_l, x_r, y_l, y_r;
46749b5e25fSSatish Balay   MatScalar    *aa;
468b0a32e0cSBarry Smith   PetscViewer   viewer;
46949b5e25fSSatish Balay 
47049b5e25fSSatish Balay   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
4729566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
47349b5e25fSSatish Balay 
47449b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
475383922c3SLisandro Dalcin 
476d0609cedSBarry Smith   PetscDrawCollectiveBegin(draw);
4779566063dSJacob Faibussowitsch   PetscCall(PetscDrawString(draw, .3 * (xl + xr), .3 * (yl + yr), PETSC_DRAW_BLACK, "symmetric"));
478383922c3SLisandro Dalcin   /* Blue for negative, Cyan for zero and  Red for positive */
479b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
48049b5e25fSSatish Balay   for (i = 0, row = 0; i < mbs; i++, row += bs) {
48149b5e25fSSatish Balay     for (j = a->i[i]; j < a->i[i + 1]; j++) {
4829371c9d4SSatish Balay       y_l = A->rmap->N - row - 1.0;
4839371c9d4SSatish Balay       y_r = y_l + 1.0;
4849371c9d4SSatish Balay       x_l = a->j[j] * bs;
4859371c9d4SSatish Balay       x_r = x_l + 1.0;
48649b5e25fSSatish Balay       aa  = a->a + j * bs2;
48749b5e25fSSatish Balay       for (k = 0; k < bs; k++) {
48849b5e25fSSatish Balay         for (l = 0; l < bs; l++) {
48949b5e25fSSatish Balay           if (PetscRealPart(*aa++) >= 0.) continue;
4909566063dSJacob Faibussowitsch           PetscCall(PetscDrawRectangle(draw, x_l + k, y_l - l, x_r + k, y_r - l, color, color, color, color));
49149b5e25fSSatish Balay         }
49249b5e25fSSatish Balay       }
49349b5e25fSSatish Balay     }
49449b5e25fSSatish Balay   }
495b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
49649b5e25fSSatish Balay   for (i = 0, row = 0; i < mbs; i++, row += bs) {
49749b5e25fSSatish Balay     for (j = a->i[i]; j < a->i[i + 1]; j++) {
4989371c9d4SSatish Balay       y_l = A->rmap->N - row - 1.0;
4999371c9d4SSatish Balay       y_r = y_l + 1.0;
5009371c9d4SSatish Balay       x_l = a->j[j] * bs;
5019371c9d4SSatish Balay       x_r = x_l + 1.0;
50249b5e25fSSatish Balay       aa  = a->a + j * bs2;
50349b5e25fSSatish Balay       for (k = 0; k < bs; k++) {
50449b5e25fSSatish Balay         for (l = 0; l < bs; l++) {
50549b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
5069566063dSJacob Faibussowitsch           PetscCall(PetscDrawRectangle(draw, x_l + k, y_l - l, x_r + k, y_r - l, color, color, color, color));
50749b5e25fSSatish Balay         }
50849b5e25fSSatish Balay       }
50949b5e25fSSatish Balay     }
51049b5e25fSSatish Balay   }
511b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
51249b5e25fSSatish Balay   for (i = 0, row = 0; i < mbs; i++, row += bs) {
51349b5e25fSSatish Balay     for (j = a->i[i]; j < a->i[i + 1]; j++) {
5149371c9d4SSatish Balay       y_l = A->rmap->N - row - 1.0;
5159371c9d4SSatish Balay       y_r = y_l + 1.0;
5169371c9d4SSatish Balay       x_l = a->j[j] * bs;
5179371c9d4SSatish Balay       x_r = x_l + 1.0;
51849b5e25fSSatish Balay       aa  = a->a + j * bs2;
51949b5e25fSSatish Balay       for (k = 0; k < bs; k++) {
52049b5e25fSSatish Balay         for (l = 0; l < bs; l++) {
52149b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
5229566063dSJacob Faibussowitsch           PetscCall(PetscDrawRectangle(draw, x_l + k, y_l - l, x_r + k, y_r - l, color, color, color, color));
52349b5e25fSSatish Balay         }
52449b5e25fSSatish Balay       }
52549b5e25fSSatish Balay     }
52649b5e25fSSatish Balay   }
527d0609cedSBarry Smith   PetscDrawCollectiveEnd(draw);
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52949b5e25fSSatish Balay }
53049b5e25fSSatish Balay 
531d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A, PetscViewer viewer)
532d71ae5a4SJacob Faibussowitsch {
53349b5e25fSSatish Balay   PetscReal xl, yl, xr, yr, w, h;
534b0a32e0cSBarry Smith   PetscDraw draw;
535ace3abfcSBarry Smith   PetscBool isnull;
53649b5e25fSSatish Balay 
53749b5e25fSSatish Balay   PetscFunctionBegin;
5389566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
5399566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
5403ba16761SJacob Faibussowitsch   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
54149b5e25fSSatish Balay 
5429371c9d4SSatish Balay   xr = A->rmap->N;
5439371c9d4SSatish Balay   yr = A->rmap->N;
5449371c9d4SSatish Balay   h  = yr / 10.0;
5459371c9d4SSatish Balay   w  = xr / 10.0;
5469371c9d4SSatish Balay   xr += w;
5479371c9d4SSatish Balay   yr += h;
5489371c9d4SSatish Balay   xl = -w;
5499371c9d4SSatish Balay   yl = -h;
5509566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
5519566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
5529566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatView_SeqSBAIJ_Draw_Zoom, A));
5539566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
5549566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
5553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55649b5e25fSSatish Balay }
55749b5e25fSSatish Balay 
558618cc2edSLisandro Dalcin /* Used for both MPIBAIJ and MPISBAIJ matrices */
559618cc2edSLisandro Dalcin #define MatView_SeqSBAIJ_Binary MatView_SeqBAIJ_Binary
560618cc2edSLisandro Dalcin 
561d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSBAIJ(Mat A, PetscViewer viewer)
562d71ae5a4SJacob Faibussowitsch {
563618cc2edSLisandro Dalcin   PetscBool iascii, isbinary, isdraw;
56449b5e25fSSatish Balay 
56549b5e25fSSatish Balay   PetscFunctionBegin;
5669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
5679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
5689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
56932077d6dSBarry Smith   if (iascii) {
5709566063dSJacob Faibussowitsch     PetscCall(MatView_SeqSBAIJ_ASCII(A, viewer));
571618cc2edSLisandro Dalcin   } else if (isbinary) {
5729566063dSJacob Faibussowitsch     PetscCall(MatView_SeqSBAIJ_Binary(A, viewer));
57349b5e25fSSatish Balay   } else if (isdraw) {
5749566063dSJacob Faibussowitsch     PetscCall(MatView_SeqSBAIJ_Draw(A, viewer));
57549b5e25fSSatish Balay   } else {
576a5e6ed63SBarry Smith     Mat         B;
577ade3a672SBarry Smith     const char *matname;
5789566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &B));
57923a3927dSBarry Smith     if (((PetscObject)A)->name) PetscCall(PetscObjectGetName((PetscObject)A, &matname));
58023a3927dSBarry Smith     if (((PetscObject)A)->name) PetscCall(PetscObjectSetName((PetscObject)B, matname));
5819566063dSJacob Faibussowitsch     PetscCall(MatView(B, viewer));
5829566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
58349b5e25fSSatish Balay   }
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58549b5e25fSSatish Balay }
58649b5e25fSSatish Balay 
587d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetValues_SeqSBAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
588d71ae5a4SJacob Faibussowitsch {
589045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
59013f74950SBarry Smith   PetscInt     *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
59113f74950SBarry Smith   PetscInt     *ai = a->i, *ailen = a->ilen;
592d0f46423SBarry Smith   PetscInt      brow, bcol, ridx, cidx, bs = A->rmap->bs, bs2 = a->bs2;
59397e567efSBarry Smith   MatScalar    *ap, *aa = a->a;
59449b5e25fSSatish Balay 
59549b5e25fSSatish Balay   PetscFunctionBegin;
59649b5e25fSSatish Balay   for (k = 0; k < m; k++) { /* loop over rows */
5979371c9d4SSatish Balay     row  = im[k];
5989371c9d4SSatish Balay     brow = row / bs;
5999371c9d4SSatish Balay     if (row < 0) {
6009371c9d4SSatish Balay       v += n;
6019371c9d4SSatish Balay       continue;
6029371c9d4SSatish Balay     } /* negative row */
60354c59aa7SJacob Faibussowitsch     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);
6049371c9d4SSatish Balay     rp   = aj + ai[brow];
6059371c9d4SSatish Balay     ap   = aa + bs2 * ai[brow];
60649b5e25fSSatish Balay     nrow = ailen[brow];
60749b5e25fSSatish Balay     for (l = 0; l < n; l++) { /* loop over columns */
6089371c9d4SSatish Balay       if (in[l] < 0) {
6099371c9d4SSatish Balay         v++;
6109371c9d4SSatish Balay         continue;
6119371c9d4SSatish Balay       } /* negative column */
61254c59aa7SJacob Faibussowitsch       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);
61349b5e25fSSatish Balay       col  = in[l];
61449b5e25fSSatish Balay       bcol = col / bs;
61549b5e25fSSatish Balay       cidx = col % bs;
61649b5e25fSSatish Balay       ridx = row % bs;
61749b5e25fSSatish Balay       high = nrow;
61849b5e25fSSatish Balay       low  = 0; /* assume unsorted */
61949b5e25fSSatish Balay       while (high - low > 5) {
62049b5e25fSSatish Balay         t = (low + high) / 2;
62149b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
62249b5e25fSSatish Balay         else low = t;
62349b5e25fSSatish Balay       }
62449b5e25fSSatish Balay       for (i = low; i < high; i++) {
62549b5e25fSSatish Balay         if (rp[i] > bcol) break;
62649b5e25fSSatish Balay         if (rp[i] == bcol) {
62749b5e25fSSatish Balay           *v++ = ap[bs2 * i + bs * cidx + ridx];
62849b5e25fSSatish Balay           goto finished;
62949b5e25fSSatish Balay         }
63049b5e25fSSatish Balay       }
63197e567efSBarry Smith       *v++ = 0.0;
63249b5e25fSSatish Balay     finished:;
63349b5e25fSSatish Balay     }
63449b5e25fSSatish Balay   }
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63649b5e25fSSatish Balay }
63749b5e25fSSatish Balay 
638d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPermute_SeqSBAIJ(Mat A, IS rowp, IS colp, Mat *B)
639d71ae5a4SJacob Faibussowitsch {
640dc29a518SPierre Jolivet   Mat C;
641dc29a518SPierre Jolivet 
642dc29a518SPierre Jolivet   PetscFunctionBegin;
6439566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSEQBAIJ, MAT_INITIAL_MATRIX, &C));
6449566063dSJacob Faibussowitsch   PetscCall(MatPermute(C, rowp, colp, B));
6459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
64648a46eb9SPierre Jolivet   if (rowp == colp) PetscCall(MatConvert(*B, MATSEQSBAIJ, MAT_INPLACE_MATRIX, B));
6473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
648dc29a518SPierre Jolivet }
64949b5e25fSSatish Balay 
650d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
651d71ae5a4SJacob Faibussowitsch {
6520880e062SHong Zhang   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ *)A->data;
653e2ee6c50SBarry Smith   PetscInt          *rp, k, low, high, t, ii, jj, row, nrow, i, col, l, rmax, N, lastcol = -1;
65413f74950SBarry Smith   PetscInt          *imax = a->imax, *ai = a->i, *ailen = a->ilen;
655d0f46423SBarry Smith   PetscInt          *aj = a->j, nonew = a->nonew, bs2 = a->bs2, bs = A->rmap->bs, stepval;
656ace3abfcSBarry Smith   PetscBool          roworiented = a->roworiented;
657dd6ea824SBarry Smith   const PetscScalar *value       = v;
658f15d580aSBarry Smith   MatScalar         *ap, *aa = a->a, *bap;
6590880e062SHong Zhang 
66049b5e25fSSatish Balay   PetscFunctionBegin;
66126fbe8dcSKarl Rupp   if (roworiented) stepval = (n - 1) * bs;
66226fbe8dcSKarl Rupp   else stepval = (m - 1) * bs;
66326fbe8dcSKarl Rupp 
6640880e062SHong Zhang   for (k = 0; k < m; k++) { /* loop over added rows */
6650880e062SHong Zhang     row = im[k];
6660880e062SHong Zhang     if (row < 0) continue;
6676bdcaf15SBarry Smith     PetscCheck(row < a->mbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Block index row too large %" PetscInt_FMT " max %" PetscInt_FMT, row, a->mbs - 1);
6680880e062SHong Zhang     rp   = aj + ai[row];
6690880e062SHong Zhang     ap   = aa + bs2 * ai[row];
6700880e062SHong Zhang     rmax = imax[row];
6710880e062SHong Zhang     nrow = ailen[row];
6720880e062SHong Zhang     low  = 0;
673818f2c47SBarry Smith     high = nrow;
6740880e062SHong Zhang     for (l = 0; l < n; l++) { /* loop over added columns */
6750880e062SHong Zhang       if (in[l] < 0) continue;
6760880e062SHong Zhang       col = in[l];
6776bdcaf15SBarry Smith       PetscCheck(col < a->nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Block index column too large %" PetscInt_FMT " max %" PetscInt_FMT, col, a->nbs - 1);
678b98bf0e1SJed Brown       if (col < row) {
67926fbe8dcSKarl Rupp         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
68026fbe8dcSKarl Rupp         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
681b98bf0e1SJed Brown       }
68226fbe8dcSKarl Rupp       if (roworiented) value = v + k * (stepval + bs) * bs + l * bs;
68326fbe8dcSKarl Rupp       else value = v + l * (stepval + bs) * bs + k * bs;
68426fbe8dcSKarl Rupp 
68526fbe8dcSKarl Rupp       if (col <= lastcol) low = 0;
68626fbe8dcSKarl Rupp       else high = nrow;
68726fbe8dcSKarl Rupp 
688e2ee6c50SBarry Smith       lastcol = col;
6890880e062SHong Zhang       while (high - low > 7) {
6900880e062SHong Zhang         t = (low + high) / 2;
6910880e062SHong Zhang         if (rp[t] > col) high = t;
6920880e062SHong Zhang         else low = t;
6930880e062SHong Zhang       }
6940880e062SHong Zhang       for (i = low; i < high; i++) {
6950880e062SHong Zhang         if (rp[i] > col) break;
6960880e062SHong Zhang         if (rp[i] == col) {
6970880e062SHong Zhang           bap = ap + bs2 * i;
6980880e062SHong Zhang           if (roworiented) {
6990880e062SHong Zhang             if (is == ADD_VALUES) {
7000880e062SHong Zhang               for (ii = 0; ii < bs; ii++, value += stepval) {
701ad540459SPierre Jolivet                 for (jj = ii; jj < bs2; jj += bs) bap[jj] += *value++;
7020880e062SHong Zhang               }
7030880e062SHong Zhang             } else {
7040880e062SHong Zhang               for (ii = 0; ii < bs; ii++, value += stepval) {
705ad540459SPierre Jolivet                 for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
7060880e062SHong Zhang               }
7070880e062SHong Zhang             }
7080880e062SHong Zhang           } else {
7090880e062SHong Zhang             if (is == ADD_VALUES) {
7100880e062SHong Zhang               for (ii = 0; ii < bs; ii++, value += stepval) {
711ad540459SPierre Jolivet                 for (jj = 0; jj < bs; jj++) *bap++ += *value++;
7120880e062SHong Zhang               }
7130880e062SHong Zhang             } else {
7140880e062SHong Zhang               for (ii = 0; ii < bs; ii++, value += stepval) {
715ad540459SPierre Jolivet                 for (jj = 0; jj < bs; jj++) *bap++ = *value++;
7160880e062SHong Zhang               }
7170880e062SHong Zhang             }
7180880e062SHong Zhang           }
7190880e062SHong Zhang           goto noinsert2;
7200880e062SHong Zhang         }
7210880e062SHong Zhang       }
7220880e062SHong Zhang       if (nonew == 1) goto noinsert2;
72308401ef6SPierre Jolivet       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new block index nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
724fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, row, col, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
7259371c9d4SSatish Balay       N = nrow++ - 1;
7269371c9d4SSatish Balay       high++;
7270880e062SHong Zhang       /* shift up all the later entries in this row */
7289566063dSJacob Faibussowitsch       PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
7299566063dSJacob Faibussowitsch       PetscCall(PetscArraymove(ap + bs2 * (i + 1), ap + bs2 * i, bs2 * (N - i + 1)));
7309566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ap + bs2 * i, bs2));
7310880e062SHong Zhang       rp[i] = col;
7320880e062SHong Zhang       bap   = ap + bs2 * i;
7330880e062SHong Zhang       if (roworiented) {
7340880e062SHong Zhang         for (ii = 0; ii < bs; ii++, value += stepval) {
735ad540459SPierre Jolivet           for (jj = ii; jj < bs2; jj += bs) bap[jj] = *value++;
7360880e062SHong Zhang         }
7370880e062SHong Zhang       } else {
7380880e062SHong Zhang         for (ii = 0; ii < bs; ii++, value += stepval) {
739ad540459SPierre Jolivet           for (jj = 0; jj < bs; jj++) *bap++ = *value++;
7400880e062SHong Zhang         }
7410880e062SHong Zhang       }
7420880e062SHong Zhang     noinsert2:;
7430880e062SHong Zhang       low = i;
7440880e062SHong Zhang     }
7450880e062SHong Zhang     ailen[row] = nrow;
7460880e062SHong Zhang   }
7473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74849b5e25fSSatish Balay }
74949b5e25fSSatish Balay 
75064831d72SBarry Smith /*
75164831d72SBarry Smith     This is not yet used
75264831d72SBarry Smith */
753d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
754d71ae5a4SJacob Faibussowitsch {
7550def2e27SBarry Smith   Mat_SeqSBAIJ   *a  = (Mat_SeqSBAIJ *)A->data;
7560def2e27SBarry Smith   const PetscInt *ai = a->i, *aj = a->j, *cols;
7570def2e27SBarry Smith   PetscInt        i = 0, j, blk_size, m = A->rmap->n, node_count = 0, nzx, nzy, *ns, row, nz, cnt, cnt2, *counts;
758ace3abfcSBarry Smith   PetscBool       flag;
7590def2e27SBarry Smith 
7600def2e27SBarry Smith   PetscFunctionBegin;
7619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &ns));
7620def2e27SBarry Smith   while (i < m) {
7630def2e27SBarry Smith     nzx = ai[i + 1] - ai[i]; /* Number of nonzeros */
7640def2e27SBarry Smith     /* Limits the number of elements in a node to 'a->inode.limit' */
7650def2e27SBarry Smith     for (j = i + 1, blk_size = 1; j < m && blk_size < a->inode.limit; ++j, ++blk_size) {
7660def2e27SBarry Smith       nzy = ai[j + 1] - ai[j];
7670def2e27SBarry Smith       if (nzy != (nzx - j + i)) break;
7689566063dSJacob Faibussowitsch       PetscCall(PetscArraycmp(aj + ai[i] + j - i, aj + ai[j], nzy, &flag));
7690def2e27SBarry Smith       if (!flag) break;
7700def2e27SBarry Smith     }
7710def2e27SBarry Smith     ns[node_count++] = blk_size;
77226fbe8dcSKarl Rupp 
7730def2e27SBarry Smith     i = j;
7740def2e27SBarry Smith   }
7750def2e27SBarry Smith   if (!a->inode.size && m && node_count > .9 * m) {
7769566063dSJacob Faibussowitsch     PetscCall(PetscFree(ns));
7779566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Found %" PetscInt_FMT " nodes out of %" PetscInt_FMT " rows. Not using Inode routines\n", node_count, m));
7780def2e27SBarry Smith   } else {
7790def2e27SBarry Smith     a->inode.node_count = node_count;
78026fbe8dcSKarl Rupp 
7819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(node_count, &a->inode.size));
7829566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->inode.size, ns, node_count));
7839566063dSJacob Faibussowitsch     PetscCall(PetscFree(ns));
7849566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Found %" PetscInt_FMT " nodes of %" PetscInt_FMT ". Limit used: %" PetscInt_FMT ". Using Inode routines\n", node_count, m, a->inode.limit));
7850def2e27SBarry Smith 
7860def2e27SBarry Smith     /* count collections of adjacent columns in each inode */
7870def2e27SBarry Smith     row = 0;
7880def2e27SBarry Smith     cnt = 0;
7890def2e27SBarry Smith     for (i = 0; i < node_count; i++) {
7900def2e27SBarry Smith       cols = aj + ai[row] + a->inode.size[i];
7910def2e27SBarry Smith       nz   = ai[row + 1] - ai[row] - a->inode.size[i];
7920def2e27SBarry Smith       for (j = 1; j < nz; j++) {
79326fbe8dcSKarl Rupp         if (cols[j] != cols[j - 1] + 1) cnt++;
7940def2e27SBarry Smith       }
7950def2e27SBarry Smith       cnt++;
7960def2e27SBarry Smith       row += a->inode.size[i];
7970def2e27SBarry Smith     }
7989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2 * cnt, &counts));
7990def2e27SBarry Smith     cnt = 0;
8000def2e27SBarry Smith     row = 0;
8010def2e27SBarry Smith     for (i = 0; i < node_count; i++) {
8020def2e27SBarry Smith       cols            = aj + ai[row] + a->inode.size[i];
8030def2e27SBarry Smith       counts[2 * cnt] = cols[0];
8040def2e27SBarry Smith       nz              = ai[row + 1] - ai[row] - a->inode.size[i];
8050def2e27SBarry Smith       cnt2            = 1;
8060def2e27SBarry Smith       for (j = 1; j < nz; j++) {
8070def2e27SBarry Smith         if (cols[j] != cols[j - 1] + 1) {
8080def2e27SBarry Smith           counts[2 * (cnt++) + 1] = cnt2;
8090def2e27SBarry Smith           counts[2 * cnt]         = cols[j];
8100def2e27SBarry Smith           cnt2                    = 1;
8110def2e27SBarry Smith         } else cnt2++;
8120def2e27SBarry Smith       }
8130def2e27SBarry Smith       counts[2 * (cnt++) + 1] = cnt2;
8140def2e27SBarry Smith       row += a->inode.size[i];
8150def2e27SBarry Smith     }
8169566063dSJacob Faibussowitsch     PetscCall(PetscIntView(2 * cnt, counts, NULL));
8170def2e27SBarry Smith   }
8183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81938702af4SBarry Smith }
82038702af4SBarry Smith 
821d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A, MatAssemblyType mode)
822d71ae5a4SJacob Faibussowitsch {
82349b5e25fSSatish Balay   Mat_SeqSBAIJ *a      = (Mat_SeqSBAIJ *)A->data;
8248f8f2f0dSBarry Smith   PetscInt      fshift = 0, i, *ai = a->i, *aj = a->j, *imax = a->imax;
825d0f46423SBarry Smith   PetscInt      m = A->rmap->N, *ip, N, *ailen = a->ilen;
82613f74950SBarry Smith   PetscInt      mbs = a->mbs, bs2 = a->bs2, rmax = 0;
82749b5e25fSSatish Balay   MatScalar    *aa = a->a, *ap;
82849b5e25fSSatish Balay 
82949b5e25fSSatish Balay   PetscFunctionBegin;
8303ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
83149b5e25fSSatish Balay 
83249b5e25fSSatish Balay   if (m) rmax = ailen[0];
83349b5e25fSSatish Balay   for (i = 1; i < mbs; i++) {
83449b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
83549b5e25fSSatish Balay     fshift += imax[i - 1] - ailen[i - 1];
83649b5e25fSSatish Balay     rmax = PetscMax(rmax, ailen[i]);
83749b5e25fSSatish Balay     if (fshift) {
838580bdb30SBarry Smith       ip = aj + ai[i];
839580bdb30SBarry Smith       ap = aa + bs2 * ai[i];
84049b5e25fSSatish Balay       N  = ailen[i];
8419566063dSJacob Faibussowitsch       PetscCall(PetscArraymove(ip - fshift, ip, N));
8429566063dSJacob Faibussowitsch       PetscCall(PetscArraymove(ap - bs2 * fshift, ap, bs2 * N));
84349b5e25fSSatish Balay     }
84449b5e25fSSatish Balay     ai[i] = ai[i - 1] + ailen[i - 1];
84549b5e25fSSatish Balay   }
84649b5e25fSSatish Balay   if (mbs) {
84749b5e25fSSatish Balay     fshift += imax[mbs - 1] - ailen[mbs - 1];
84849b5e25fSSatish Balay     ai[mbs] = ai[mbs - 1] + ailen[mbs - 1];
84949b5e25fSSatish Balay   }
85049b5e25fSSatish Balay   /* reset ilen and imax for each row */
851ad540459SPierre Jolivet   for (i = 0; i < mbs; i++) ailen[i] = imax[i] = ai[i + 1] - ai[i];
8526c6c5352SBarry Smith   a->nz = ai[mbs];
85349b5e25fSSatish Balay 
854b424e231SHong Zhang   /* diagonals may have moved, reset it */
8551baa6e33SBarry Smith   if (a->diag) PetscCall(PetscArraycpy(a->diag, ai, mbs));
856aed4548fSBarry Smith   PetscCheck(!fshift || a->nounused != -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unused space detected in matrix: %" PetscInt_FMT " X %" PetscInt_FMT " block size %" PetscInt_FMT ", %" PetscInt_FMT " unneeded", m, A->cmap->n, A->rmap->bs, fshift * bs2);
85726fbe8dcSKarl Rupp 
8589566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT ", block size %" PetscInt_FMT "; storage space: %" PetscInt_FMT " unneeded, %" PetscInt_FMT " used\n", m, A->rmap->N, A->rmap->bs, fshift * bs2, a->nz * bs2));
8599566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues is %" PetscInt_FMT "\n", a->reallocs));
8609566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Most nonzeros blocks in any row is %" PetscInt_FMT "\n", rmax));
86126fbe8dcSKarl Rupp 
8628e58a170SBarry Smith   A->info.mallocs += a->reallocs;
86349b5e25fSSatish Balay   a->reallocs         = 0;
86449b5e25fSSatish Balay   A->info.nz_unneeded = (PetscReal)fshift * bs2;
865061b2667SBarry Smith   a->idiagvalid       = PETSC_FALSE;
8664dcd73b1SHong Zhang   a->rmax             = rmax;
86738702af4SBarry Smith 
86838702af4SBarry Smith   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
86944e1c64aSLisandro Dalcin     if (a->jshort && a->free_jshort) {
87017803ae8SHong Zhang       /* when matrix data structure is changed, previous jshort must be replaced */
8719566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->jshort));
87217803ae8SHong Zhang     }
8739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(a->i[A->rmap->n], &a->jshort));
87438702af4SBarry Smith     for (i = 0; i < a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
87538702af4SBarry Smith     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
87641f059aeSBarry Smith     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
8774da8f245SBarry Smith     a->free_jshort = PETSC_TRUE;
87838702af4SBarry Smith   }
8793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88049b5e25fSSatish Balay }
88149b5e25fSSatish Balay 
88249b5e25fSSatish Balay /*
88349b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
88449b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
885a5b23f4aSJose E. Roman    then the resulting sizes = [3,1,1,3,1] corresponding to sets [(0,1,2),(3),(5),(6,7,8),(9)]
88649b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
88749b5e25fSSatish Balay */
888d71ae5a4SJacob Faibussowitsch PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[], PetscInt n, PetscInt bs, PetscInt sizes[], PetscInt *bs_max)
889d71ae5a4SJacob Faibussowitsch {
89013f74950SBarry Smith   PetscInt  i, j, k, row;
891ace3abfcSBarry Smith   PetscBool flg;
89249b5e25fSSatish Balay 
89349b5e25fSSatish Balay   PetscFunctionBegin;
89449b5e25fSSatish Balay   for (i = 0, j = 0; i < n; j++) {
89549b5e25fSSatish Balay     row = idx[i];
896a5b23f4aSJose E. Roman     if (row % bs != 0) { /* Not the beginning of a block */
89749b5e25fSSatish Balay       sizes[j] = 1;
89849b5e25fSSatish Balay       i++;
89949b5e25fSSatish Balay     } else if (i + bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
90049b5e25fSSatish Balay       sizes[j] = 1;          /* Also makes sure at least 'bs' values exist for next else */
90149b5e25fSSatish Balay       i++;
9026aad120cSJose E. Roman     } else { /* Beginning of the block, so check if the complete block exists */
90349b5e25fSSatish Balay       flg = PETSC_TRUE;
90449b5e25fSSatish Balay       for (k = 1; k < bs; k++) {
90549b5e25fSSatish Balay         if (row + k != idx[i + k]) { /* break in the block */
90649b5e25fSSatish Balay           flg = PETSC_FALSE;
90749b5e25fSSatish Balay           break;
90849b5e25fSSatish Balay         }
90949b5e25fSSatish Balay       }
910abc0a331SBarry Smith       if (flg) { /* No break in the bs */
91149b5e25fSSatish Balay         sizes[j] = bs;
91249b5e25fSSatish Balay         i += bs;
91349b5e25fSSatish Balay       } else {
91449b5e25fSSatish Balay         sizes[j] = 1;
91549b5e25fSSatish Balay         i++;
91649b5e25fSSatish Balay       }
91749b5e25fSSatish Balay     }
91849b5e25fSSatish Balay   }
91949b5e25fSSatish Balay   *bs_max = j;
9203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92149b5e25fSSatish Balay }
92249b5e25fSSatish Balay 
92349b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
924da81f932SPierre Jolivet    Any a(i,j) with i>j input by user is ignored.
92549b5e25fSSatish Balay */
92649b5e25fSSatish Balay 
927d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValues_SeqSBAIJ(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
928d71ae5a4SJacob Faibussowitsch {
92949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
930e2ee6c50SBarry Smith   PetscInt     *rp, k, low, high, t, ii, row, nrow, i, col, l, rmax, N, lastcol = -1;
93113f74950SBarry Smith   PetscInt     *imax = a->imax, *ai = a->i, *ailen = a->ilen, roworiented = a->roworiented;
932d0f46423SBarry Smith   PetscInt     *aj = a->j, nonew = a->nonew, bs = A->rmap->bs, brow, bcol;
93313f74950SBarry Smith   PetscInt      ridx, cidx, bs2                 = a->bs2;
93449b5e25fSSatish Balay   MatScalar    *ap, value, *aa                  = a->a, *bap;
93549b5e25fSSatish Balay 
93649b5e25fSSatish Balay   PetscFunctionBegin;
93749b5e25fSSatish Balay   for (k = 0; k < m; k++) { /* loop over added rows */
93849b5e25fSSatish Balay     row  = im[k];           /* row number */
93949b5e25fSSatish Balay     brow = row / bs;        /* block row number */
94049b5e25fSSatish Balay     if (row < 0) continue;
9416bdcaf15SBarry 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);
94249b5e25fSSatish Balay     rp   = aj + ai[brow];       /*ptr to beginning of column value of the row block*/
94349b5e25fSSatish Balay     ap   = aa + bs2 * ai[brow]; /*ptr to beginning of element value of the row block*/
94449b5e25fSSatish Balay     rmax = imax[brow];          /* maximum space allocated for this row */
94549b5e25fSSatish Balay     nrow = ailen[brow];         /* actual length of this row */
94649b5e25fSSatish Balay     low  = 0;
9478509e838SStefano Zampini     high = nrow;
94849b5e25fSSatish Balay     for (l = 0; l < n; l++) { /* loop over added columns */
94949b5e25fSSatish Balay       if (in[l] < 0) continue;
9506bdcaf15SBarry 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);
95149b5e25fSSatish Balay       col  = in[l];
95249b5e25fSSatish Balay       bcol = col / bs; /* block col number */
95349b5e25fSSatish Balay 
954941593c8SHong Zhang       if (brow > bcol) {
95526fbe8dcSKarl Rupp         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
95626fbe8dcSKarl Rupp         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
957941593c8SHong Zhang       }
958f4989cb3SHong Zhang 
9599371c9d4SSatish Balay       ridx = row % bs;
9609371c9d4SSatish Balay       cidx = col % bs; /*row and col index inside the block */
9618549e402SHong Zhang       if ((brow == bcol && ridx <= cidx) || (brow < bcol)) {
96249b5e25fSSatish Balay         /* element value a(k,l) */
96326fbe8dcSKarl Rupp         if (roworiented) value = v[l + k * n];
96426fbe8dcSKarl Rupp         else value = v[k + l * m];
96549b5e25fSSatish Balay 
96649b5e25fSSatish Balay         /* move pointer bap to a(k,l) quickly and add/insert value */
96726fbe8dcSKarl Rupp         if (col <= lastcol) low = 0;
9688509e838SStefano Zampini         else high = nrow;
9698509e838SStefano Zampini 
970e2ee6c50SBarry Smith         lastcol = col;
97149b5e25fSSatish Balay         while (high - low > 7) {
97249b5e25fSSatish Balay           t = (low + high) / 2;
97349b5e25fSSatish Balay           if (rp[t] > bcol) high = t;
97449b5e25fSSatish Balay           else low = t;
97549b5e25fSSatish Balay         }
97649b5e25fSSatish Balay         for (i = low; i < high; i++) {
97749b5e25fSSatish Balay           if (rp[i] > bcol) break;
97849b5e25fSSatish Balay           if (rp[i] == bcol) {
97949b5e25fSSatish Balay             bap = ap + bs2 * i + bs * cidx + ridx;
98049b5e25fSSatish Balay             if (is == ADD_VALUES) *bap += value;
98149b5e25fSSatish Balay             else *bap = value;
9828549e402SHong Zhang             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9838549e402SHong Zhang             if (brow == bcol && ridx < cidx) {
9848549e402SHong Zhang               bap = ap + bs2 * i + bs * ridx + cidx;
9858549e402SHong Zhang               if (is == ADD_VALUES) *bap += value;
9868549e402SHong Zhang               else *bap = value;
9878549e402SHong Zhang             }
98849b5e25fSSatish Balay             goto noinsert1;
98949b5e25fSSatish Balay           }
99049b5e25fSSatish Balay         }
99149b5e25fSSatish Balay 
99249b5e25fSSatish Balay         if (nonew == 1) goto noinsert1;
99308401ef6SPierre Jolivet         PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
994fef13f97SBarry Smith         MatSeqXAIJReallocateAIJ(A, a->mbs, bs2, nrow, brow, bcol, rmax, aa, ai, aj, rp, ap, imax, nonew, MatScalar);
99549b5e25fSSatish Balay 
9969371c9d4SSatish Balay         N = nrow++ - 1;
9979371c9d4SSatish Balay         high++;
99849b5e25fSSatish Balay         /* shift up all the later entries in this row */
9999566063dSJacob Faibussowitsch         PetscCall(PetscArraymove(rp + i + 1, rp + i, N - i + 1));
10009566063dSJacob Faibussowitsch         PetscCall(PetscArraymove(ap + bs2 * (i + 1), ap + bs2 * i, bs2 * (N - i + 1)));
10019566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(ap + bs2 * i, bs2));
100249b5e25fSSatish Balay         rp[i]                          = bcol;
100349b5e25fSSatish Balay         ap[bs2 * i + bs * cidx + ridx] = value;
10048509e838SStefano Zampini         /* for diag block, add/insert its symmetric element a(cidx,ridx) */
1005ad540459SPierre Jolivet         if (brow == bcol && ridx < cidx) ap[bs2 * i + bs * ridx + cidx] = value;
1006e56f5c9eSBarry Smith         A->nonzerostate++;
100749b5e25fSSatish Balay       noinsert1:;
100849b5e25fSSatish Balay         low = i;
10098549e402SHong Zhang       }
101049b5e25fSSatish Balay     } /* end of loop over added columns */
101149b5e25fSSatish Balay     ailen[brow] = nrow;
101249b5e25fSSatish Balay   } /* end of loop over added rows */
10133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
101449b5e25fSSatish Balay }
101549b5e25fSSatish Balay 
1016d71ae5a4SJacob Faibussowitsch PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA, IS row, const MatFactorInfo *info)
1017d71ae5a4SJacob Faibussowitsch {
10184ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)inA->data;
101949b5e25fSSatish Balay   Mat           outA;
1020ace3abfcSBarry Smith   PetscBool     row_identity;
102149b5e25fSSatish Balay 
102249b5e25fSSatish Balay   PetscFunctionBegin;
102308401ef6SPierre Jolivet   PetscCheck(info->levels == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only levels=0 is supported for in-place icc");
10249566063dSJacob Faibussowitsch   PetscCall(ISIdentity(row, &row_identity));
102528b400f6SJacob Faibussowitsch   PetscCheck(row_identity, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix reordering is not supported");
102608401ef6SPierre Jolivet   PetscCheck(inA->rmap->bs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix block size %" PetscInt_FMT " is not supported", inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */
1027c84f5b01SHong Zhang 
102849b5e25fSSatish Balay   outA            = inA;
1029d5f3da31SBarry Smith   inA->factortype = MAT_FACTOR_ICC;
10309566063dSJacob Faibussowitsch   PetscCall(PetscFree(inA->solvertype));
10319566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &inA->solvertype));
103249b5e25fSSatish Balay 
10339566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSBAIJ(inA));
10349566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetNumericFactorization_inplace(inA, row_identity));
103549b5e25fSSatish Balay 
10369566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)row));
10379566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->row));
1038c84f5b01SHong Zhang   a->row = row;
10399566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)row));
10409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->col));
1041c84f5b01SHong Zhang   a->col = row;
1042c84f5b01SHong Zhang 
1043c84f5b01SHong Zhang   /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
10449566063dSJacob Faibussowitsch   if (a->icol) PetscCall(ISInvertPermutation(row, PETSC_DECIDE, &a->icol));
104549b5e25fSSatish Balay 
1046aa624791SPierre Jolivet   if (!a->solve_work) PetscCall(PetscMalloc1(inA->rmap->N + inA->rmap->bs, &a->solve_work));
104749b5e25fSSatish Balay 
10489566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(outA, inA, info));
10493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105049b5e25fSSatish Balay }
1051950f1e5bSHong Zhang 
1052d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat, PetscInt *indices)
1053d71ae5a4SJacob Faibussowitsch {
1054045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
105513f74950SBarry Smith   PetscInt      i, nz, n;
105649b5e25fSSatish Balay 
105749b5e25fSSatish Balay   PetscFunctionBegin;
10586c6c5352SBarry Smith   nz = baij->maxnz;
1059d0f46423SBarry Smith   n  = mat->cmap->n;
106026fbe8dcSKarl Rupp   for (i = 0; i < nz; i++) baij->j[i] = indices[i];
106126fbe8dcSKarl Rupp 
10626c6c5352SBarry Smith   baij->nz = nz;
106326fbe8dcSKarl Rupp   for (i = 0; i < n; i++) baij->ilen[i] = baij->imax[i];
106426fbe8dcSKarl Rupp 
10659566063dSJacob Faibussowitsch   PetscCall(MatSetOption(mat, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
10663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106749b5e25fSSatish Balay }
106849b5e25fSSatish Balay 
106949b5e25fSSatish Balay /*@
107019585528SSatish Balay   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
107111a5261eSBarry Smith   in a `MATSEQSBAIJ` matrix.
107249b5e25fSSatish Balay 
107349b5e25fSSatish Balay   Input Parameters:
107411a5261eSBarry Smith +  mat     - the `MATSEQSBAIJ` matrix
107549b5e25fSSatish Balay -  indices - the column indices
107649b5e25fSSatish Balay 
107749b5e25fSSatish Balay   Level: advanced
107849b5e25fSSatish Balay 
107949b5e25fSSatish Balay   Notes:
108049b5e25fSSatish Balay   This can be called if you have precomputed the nonzero structure of the
108149b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
108211a5261eSBarry Smith   of the `MatSetValues()` operation.
108349b5e25fSSatish Balay 
108449b5e25fSSatish Balay   You MUST have set the correct numbers of nonzeros per row in the call to
108511a5261eSBarry Smith   `MatCreateSeqSBAIJ()`, and the columns indices MUST be sorted.
108649b5e25fSSatish Balay 
10872ef1f0ffSBarry Smith   MUST be called before any calls to `MatSetValues()`
108849b5e25fSSatish Balay 
10891cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreateSeqSBAIJ`
109049b5e25fSSatish Balay @*/
1091d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat, PetscInt *indices)
1092d71ae5a4SJacob Faibussowitsch {
109349b5e25fSSatish Balay   PetscFunctionBegin;
10940700a824SBarry Smith   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
1095dadcf809SJacob Faibussowitsch   PetscValidIntPointer(indices, 2);
1096cac4c232SBarry Smith   PetscUseMethod(mat, "MatSeqSBAIJSetColumnIndices_C", (Mat, PetscInt *), (mat, indices));
10973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109849b5e25fSSatish Balay }
109949b5e25fSSatish Balay 
1100d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_SeqSBAIJ(Mat A, Mat B, MatStructure str)
1101d71ae5a4SJacob Faibussowitsch {
11024c7a3774SStefano Zampini   PetscBool isbaij;
11033c896bc6SHong Zhang 
11043c896bc6SHong Zhang   PetscFunctionBegin;
11059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &isbaij, MATSEQSBAIJ, MATMPISBAIJ, ""));
110628b400f6SJacob Faibussowitsch   PetscCheck(isbaij, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "Not for matrix type %s", ((PetscObject)B)->type_name);
11074c7a3774SStefano Zampini   /* If the two matrices have the same copy implementation and nonzero pattern, use fast copy. */
11084c7a3774SStefano Zampini   if (str == SAME_NONZERO_PATTERN && A->ops->copy == B->ops->copy) {
11093c896bc6SHong Zhang     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
11103c896bc6SHong Zhang     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)B->data;
11113c896bc6SHong Zhang 
111208401ef6SPierre Jolivet     PetscCheck(a->i[a->mbs] == b->i[b->mbs], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of nonzeros in two matrices are different");
111308401ef6SPierre Jolivet     PetscCheck(a->mbs == b->mbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of rows in two matrices are different");
111408401ef6SPierre Jolivet     PetscCheck(a->bs2 == b->bs2, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Different block size");
11159566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->a, a->a, a->bs2 * a->i[a->mbs]));
11169566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateIncrease((PetscObject)B));
11173c896bc6SHong Zhang   } else {
11189566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
11199566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
11209566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
11213c896bc6SHong Zhang   }
11223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11233c896bc6SHong Zhang }
11243c896bc6SHong Zhang 
1125d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A, PetscScalar *array[])
1126d71ae5a4SJacob Faibussowitsch {
1127a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
11285fd66863SKarl Rupp 
1129a6ece127SHong Zhang   PetscFunctionBegin;
1130a6ece127SHong Zhang   *array = a->a;
11313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1132a6ece127SHong Zhang }
1133a6ece127SHong Zhang 
1134d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A, PetscScalar *array[])
1135d71ae5a4SJacob Faibussowitsch {
1136a6ece127SHong Zhang   PetscFunctionBegin;
1137cda14afcSprj-   *array = NULL;
11383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1139a6ece127SHong Zhang }
1140a6ece127SHong Zhang 
1141d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y, Mat X, PetscInt *nnz)
1142d71ae5a4SJacob Faibussowitsch {
1143b264fe52SHong Zhang   PetscInt      bs = Y->rmap->bs, mbs = Y->rmap->N / bs;
114452768537SHong Zhang   Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ *)X->data;
114552768537SHong Zhang   Mat_SeqSBAIJ *y = (Mat_SeqSBAIJ *)Y->data;
114652768537SHong Zhang 
114752768537SHong Zhang   PetscFunctionBegin;
114852768537SHong Zhang   /* Set the number of nonzeros in the new matrix */
11499566063dSJacob Faibussowitsch   PetscCall(MatAXPYGetPreallocation_SeqX_private(mbs, x->i, x->j, y->i, y->j, nnz));
11503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115152768537SHong Zhang }
115252768537SHong Zhang 
1153d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y, PetscScalar a, Mat X, MatStructure str)
1154d71ae5a4SJacob Faibussowitsch {
115542ee4b1aSHong Zhang   Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ *)X->data, *y = (Mat_SeqSBAIJ *)Y->data;
115631ce2d13SHong Zhang   PetscInt      bs = Y->rmap->bs, bs2 = bs * bs;
1157e838b9e7SJed Brown   PetscBLASInt  one = 1;
115842ee4b1aSHong Zhang 
115942ee4b1aSHong Zhang   PetscFunctionBegin;
1160134adf20SPierre Jolivet   if (str == UNKNOWN_NONZERO_PATTERN || (PetscDefined(USE_DEBUG) && str == SAME_NONZERO_PATTERN)) {
1161134adf20SPierre Jolivet     PetscBool e = x->nz == y->nz && x->mbs == y->mbs ? PETSC_TRUE : PETSC_FALSE;
1162134adf20SPierre Jolivet     if (e) {
11639566063dSJacob Faibussowitsch       PetscCall(PetscArraycmp(x->i, y->i, x->mbs + 1, &e));
1164134adf20SPierre Jolivet       if (e) {
11659566063dSJacob Faibussowitsch         PetscCall(PetscArraycmp(x->j, y->j, x->i[x->mbs], &e));
1166134adf20SPierre Jolivet         if (e) str = SAME_NONZERO_PATTERN;
1167134adf20SPierre Jolivet       }
1168134adf20SPierre Jolivet     }
116954c59aa7SJacob Faibussowitsch     if (!e) PetscCheck(str != SAME_NONZERO_PATTERN, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatStructure is not SAME_NONZERO_PATTERN");
1170134adf20SPierre Jolivet   }
117142ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1172f4df32b1SMatthew Knepley     PetscScalar  alpha = a;
1173c5df96a5SBarry Smith     PetscBLASInt bnz;
11749566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(x->nz * bs2, &bnz));
1175792fecdfSBarry Smith     PetscCallBLAS("BLASaxpy", BLASaxpy_(&bnz, &alpha, x->a, &one, y->a, &one));
11769566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
1177ab784542SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
11789566063dSJacob Faibussowitsch     PetscCall(MatSetOption(X, MAT_GETROW_UPPERTRIANGULAR, PETSC_TRUE));
11799566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic(Y, a, X, str));
11809566063dSJacob Faibussowitsch     PetscCall(MatSetOption(X, MAT_GETROW_UPPERTRIANGULAR, PETSC_FALSE));
118142ee4b1aSHong Zhang   } else {
118252768537SHong Zhang     Mat       B;
118352768537SHong Zhang     PetscInt *nnz;
118454c59aa7SJacob Faibussowitsch     PetscCheck(bs == X->rmap->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Matrices must have same block size");
11859566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
11869566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
11879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Y->rmap->N, &nnz));
11889566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y), &B));
11899566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)B, ((PetscObject)Y)->name));
11909566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, Y->rmap->n, Y->cmap->n, Y->rmap->N, Y->cmap->N));
11919566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(B, Y, Y));
11929566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, ((PetscObject)Y)->type_name));
11939566063dSJacob Faibussowitsch     PetscCall(MatAXPYGetPreallocation_SeqSBAIJ(Y, X, nnz));
11949566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, nnz));
119552768537SHong Zhang 
11969566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B, Y, a, X, str));
119752768537SHong Zhang 
11989566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y, &B));
11999566063dSJacob Faibussowitsch     PetscCall(PetscFree(nnz));
12009566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
12019566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
120242ee4b1aSHong Zhang   }
12033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120442ee4b1aSHong Zhang }
120542ee4b1aSHong Zhang 
1206d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A, PetscReal tol, PetscBool *flg)
1207d71ae5a4SJacob Faibussowitsch {
1208efcf0fc3SBarry Smith   PetscFunctionBegin;
1209efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
12103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1211efcf0fc3SBarry Smith }
1212efcf0fc3SBarry Smith 
1213d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A, PetscBool *flg)
1214d71ae5a4SJacob Faibussowitsch {
1215efcf0fc3SBarry Smith   PetscFunctionBegin;
1216efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
12173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1218efcf0fc3SBarry Smith }
1219efcf0fc3SBarry Smith 
1220d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A, PetscReal tol, PetscBool *flg)
1221d71ae5a4SJacob Faibussowitsch {
1222efcf0fc3SBarry Smith   PetscFunctionBegin;
1223efcf0fc3SBarry Smith   *flg = PETSC_FALSE;
12243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1225efcf0fc3SBarry Smith }
1226efcf0fc3SBarry Smith 
1227d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConjugate_SeqSBAIJ(Mat A)
1228d71ae5a4SJacob Faibussowitsch {
12292726fb6dSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
12302726fb6dSPierre Jolivet   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
12312726fb6dSPierre Jolivet   PetscInt      i, nz = a->bs2 * a->i[a->mbs];
12322726fb6dSPierre Jolivet   MatScalar    *aa = a->a;
12332726fb6dSPierre Jolivet 
12342726fb6dSPierre Jolivet   PetscFunctionBegin;
12352726fb6dSPierre Jolivet   for (i = 0; i < nz; i++) aa[i] = PetscConj(aa[i]);
12362726fb6dSPierre Jolivet #else
12372726fb6dSPierre Jolivet   PetscFunctionBegin;
12382726fb6dSPierre Jolivet #endif
12393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12402726fb6dSPierre Jolivet }
12412726fb6dSPierre Jolivet 
1242d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1243d71ae5a4SJacob Faibussowitsch {
124499cafbc1SBarry Smith   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
124599cafbc1SBarry Smith   PetscInt      i, nz = a->bs2 * a->i[a->mbs];
1246dd6ea824SBarry Smith   MatScalar    *aa = a->a;
124799cafbc1SBarry Smith 
124899cafbc1SBarry Smith   PetscFunctionBegin;
124999cafbc1SBarry Smith   for (i = 0; i < nz; i++) aa[i] = PetscRealPart(aa[i]);
12503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125199cafbc1SBarry Smith }
125299cafbc1SBarry Smith 
1253d71ae5a4SJacob Faibussowitsch PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1254d71ae5a4SJacob Faibussowitsch {
125599cafbc1SBarry Smith   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data;
125699cafbc1SBarry Smith   PetscInt      i, nz = a->bs2 * a->i[a->mbs];
1257dd6ea824SBarry Smith   MatScalar    *aa = a->a;
125899cafbc1SBarry Smith 
125999cafbc1SBarry Smith   PetscFunctionBegin;
126099cafbc1SBarry Smith   for (i = 0; i < nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
12613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126299cafbc1SBarry Smith }
126399cafbc1SBarry Smith 
1264d71ae5a4SJacob Faibussowitsch PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A, PetscInt is_n, const PetscInt is_idx[], PetscScalar diag, Vec x, Vec b)
1265d71ae5a4SJacob Faibussowitsch {
12663bededecSBarry Smith   Mat_SeqSBAIJ      *baij = (Mat_SeqSBAIJ *)A->data;
12673bededecSBarry Smith   PetscInt           i, j, k, count;
12683bededecSBarry Smith   PetscInt           bs = A->rmap->bs, bs2 = baij->bs2, row, col;
12693bededecSBarry Smith   PetscScalar        zero = 0.0;
12703bededecSBarry Smith   MatScalar         *aa;
12713bededecSBarry Smith   const PetscScalar *xx;
12723bededecSBarry Smith   PetscScalar       *bb;
127356777dd2SBarry Smith   PetscBool         *zeroed, vecs = PETSC_FALSE;
12743bededecSBarry Smith 
12753bededecSBarry Smith   PetscFunctionBegin;
12763bededecSBarry Smith   /* fix right hand side if needed */
12773bededecSBarry Smith   if (x && b) {
12789566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x, &xx));
12799566063dSJacob Faibussowitsch     PetscCall(VecGetArray(b, &bb));
128056777dd2SBarry Smith     vecs = PETSC_TRUE;
12813bededecSBarry Smith   }
12823bededecSBarry Smith 
12833bededecSBarry Smith   /* zero the columns */
12849566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->n, &zeroed));
12853bededecSBarry Smith   for (i = 0; i < is_n; i++) {
1286aed4548fSBarry Smith     PetscCheck(is_idx[i] >= 0 && is_idx[i] < A->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "row %" PetscInt_FMT " out of range", is_idx[i]);
12873bededecSBarry Smith     zeroed[is_idx[i]] = PETSC_TRUE;
12883bededecSBarry Smith   }
128956777dd2SBarry Smith   if (vecs) {
129056777dd2SBarry Smith     for (i = 0; i < A->rmap->N; i++) {
129156777dd2SBarry Smith       row = i / bs;
129256777dd2SBarry Smith       for (j = baij->i[row]; j < baij->i[row + 1]; j++) {
129356777dd2SBarry Smith         for (k = 0; k < bs; k++) {
129456777dd2SBarry Smith           col = bs * baij->j[j] + k;
129556777dd2SBarry Smith           if (col <= i) continue;
129656777dd2SBarry Smith           aa = ((MatScalar *)(baij->a)) + j * bs2 + (i % bs) + bs * k;
129726fbe8dcSKarl Rupp           if (!zeroed[i] && zeroed[col]) bb[i] -= aa[0] * xx[col];
129826fbe8dcSKarl Rupp           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0] * xx[i];
129956777dd2SBarry Smith         }
130056777dd2SBarry Smith       }
130156777dd2SBarry Smith     }
130226fbe8dcSKarl Rupp     for (i = 0; i < is_n; i++) bb[is_idx[i]] = diag * xx[is_idx[i]];
130356777dd2SBarry Smith   }
130456777dd2SBarry Smith 
13053bededecSBarry Smith   for (i = 0; i < A->rmap->N; i++) {
13063bededecSBarry Smith     if (!zeroed[i]) {
13073bededecSBarry Smith       row = i / bs;
13083bededecSBarry Smith       for (j = baij->i[row]; j < baij->i[row + 1]; j++) {
13093bededecSBarry Smith         for (k = 0; k < bs; k++) {
13103bededecSBarry Smith           col = bs * baij->j[j] + k;
13113bededecSBarry Smith           if (zeroed[col]) {
13123bededecSBarry Smith             aa    = ((MatScalar *)(baij->a)) + j * bs2 + (i % bs) + bs * k;
13133bededecSBarry Smith             aa[0] = 0.0;
13143bededecSBarry Smith           }
13153bededecSBarry Smith         }
13163bededecSBarry Smith       }
13173bededecSBarry Smith     }
13183bededecSBarry Smith   }
13199566063dSJacob Faibussowitsch   PetscCall(PetscFree(zeroed));
132056777dd2SBarry Smith   if (vecs) {
13219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x, &xx));
13229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(b, &bb));
132356777dd2SBarry Smith   }
13243bededecSBarry Smith 
13253bededecSBarry Smith   /* zero the rows */
13263bededecSBarry Smith   for (i = 0; i < is_n; i++) {
13273bededecSBarry Smith     row   = is_idx[i];
13283bededecSBarry Smith     count = (baij->i[row / bs + 1] - baij->i[row / bs]) * bs;
13293bededecSBarry Smith     aa    = ((MatScalar *)(baij->a)) + baij->i[row / bs] * bs2 + (row % bs);
13303bededecSBarry Smith     for (k = 0; k < count; k++) {
13313bededecSBarry Smith       aa[0] = zero;
13323bededecSBarry Smith       aa += bs;
13333bededecSBarry Smith     }
1334dbbe0bcdSBarry Smith     if (diag != 0.0) PetscUseTypeMethod(A, setvalues, 1, &row, 1, &row, &diag, INSERT_VALUES);
13353bededecSBarry Smith   }
13369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd_SeqSBAIJ(A, MAT_FINAL_ASSEMBLY));
13373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13383bededecSBarry Smith }
13393bededecSBarry Smith 
1340d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_SeqSBAIJ(Mat Y, PetscScalar a)
1341d71ae5a4SJacob Faibussowitsch {
13427d68702bSBarry Smith   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)Y->data;
13437d68702bSBarry Smith 
13447d68702bSBarry Smith   PetscFunctionBegin;
134548a46eb9SPierre Jolivet   if (!Y->preallocated || !aij->nz) PetscCall(MatSeqSBAIJSetPreallocation(Y, Y->rmap->bs, 1, NULL));
13469566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
13473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13487d68702bSBarry Smith }
13497d68702bSBarry Smith 
13503964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
135149b5e25fSSatish Balay                                        MatGetRow_SeqSBAIJ,
135249b5e25fSSatish Balay                                        MatRestoreRow_SeqSBAIJ,
135349b5e25fSSatish Balay                                        MatMult_SeqSBAIJ_N,
135497304618SKris Buschelman                                        /*  4*/ MatMultAdd_SeqSBAIJ_N,
1355431c96f7SBarry Smith                                        MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */
1356e005ede5SBarry Smith                                        MatMultAdd_SeqSBAIJ_N,
1357f4259b30SLisandro Dalcin                                        NULL,
1358f4259b30SLisandro Dalcin                                        NULL,
1359f4259b30SLisandro Dalcin                                        NULL,
1360f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
1361f4259b30SLisandro Dalcin                                        NULL,
1362c078aec8SLisandro Dalcin                                        MatCholeskyFactor_SeqSBAIJ,
136341f059aeSBarry Smith                                        MatSOR_SeqSBAIJ,
136449b5e25fSSatish Balay                                        MatTranspose_SeqSBAIJ,
136597304618SKris Buschelman                                        /* 15*/ MatGetInfo_SeqSBAIJ,
136649b5e25fSSatish Balay                                        MatEqual_SeqSBAIJ,
136749b5e25fSSatish Balay                                        MatGetDiagonal_SeqSBAIJ,
136849b5e25fSSatish Balay                                        MatDiagonalScale_SeqSBAIJ,
136949b5e25fSSatish Balay                                        MatNorm_SeqSBAIJ,
1370f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
137149b5e25fSSatish Balay                                        MatAssemblyEnd_SeqSBAIJ,
137249b5e25fSSatish Balay                                        MatSetOption_SeqSBAIJ,
137349b5e25fSSatish Balay                                        MatZeroEntries_SeqSBAIJ,
1374f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
1375f4259b30SLisandro Dalcin                                        NULL,
1376f4259b30SLisandro Dalcin                                        NULL,
1377f4259b30SLisandro Dalcin                                        NULL,
1378f4259b30SLisandro Dalcin                                        NULL,
137926cec326SBarry Smith                                        /* 29*/ MatSetUp_Seq_Hash,
1380f4259b30SLisandro Dalcin                                        NULL,
1381f4259b30SLisandro Dalcin                                        NULL,
1382f4259b30SLisandro Dalcin                                        NULL,
1383f4259b30SLisandro Dalcin                                        NULL,
1384d519adbfSMatthew Knepley                                        /* 34*/ MatDuplicate_SeqSBAIJ,
1385f4259b30SLisandro Dalcin                                        NULL,
1386f4259b30SLisandro Dalcin                                        NULL,
1387f4259b30SLisandro Dalcin                                        NULL,
1388c84f5b01SHong Zhang                                        MatICCFactor_SeqSBAIJ,
1389d519adbfSMatthew Knepley                                        /* 39*/ MatAXPY_SeqSBAIJ,
13907dae84e0SHong Zhang                                        MatCreateSubMatrices_SeqSBAIJ,
139149b5e25fSSatish Balay                                        MatIncreaseOverlap_SeqSBAIJ,
139249b5e25fSSatish Balay                                        MatGetValues_SeqSBAIJ,
13933c896bc6SHong Zhang                                        MatCopy_SeqSBAIJ,
1394f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
139549b5e25fSSatish Balay                                        MatScale_SeqSBAIJ,
13967d68702bSBarry Smith                                        MatShift_SeqSBAIJ,
1397f4259b30SLisandro Dalcin                                        NULL,
13983bededecSBarry Smith                                        MatZeroRowsColumns_SeqSBAIJ,
1399f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
140049b5e25fSSatish Balay                                        MatGetRowIJ_SeqSBAIJ,
140149b5e25fSSatish Balay                                        MatRestoreRowIJ_SeqSBAIJ,
1402f4259b30SLisandro Dalcin                                        NULL,
1403f4259b30SLisandro Dalcin                                        NULL,
1404f4259b30SLisandro Dalcin                                        /* 54*/ NULL,
1405f4259b30SLisandro Dalcin                                        NULL,
1406f4259b30SLisandro Dalcin                                        NULL,
1407dc29a518SPierre Jolivet                                        MatPermute_SeqSBAIJ,
140849b5e25fSSatish Balay                                        MatSetValuesBlocked_SeqSBAIJ,
14097dae84e0SHong Zhang                                        /* 59*/ MatCreateSubMatrix_SeqSBAIJ,
1410f4259b30SLisandro Dalcin                                        NULL,
1411f4259b30SLisandro Dalcin                                        NULL,
1412f4259b30SLisandro Dalcin                                        NULL,
1413f4259b30SLisandro Dalcin                                        NULL,
1414f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
1415f4259b30SLisandro Dalcin                                        NULL,
1416f4259b30SLisandro Dalcin                                        NULL,
1417f4259b30SLisandro Dalcin                                        NULL,
1418f4259b30SLisandro Dalcin                                        NULL,
1419d519adbfSMatthew Knepley                                        /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1420f4259b30SLisandro Dalcin                                        NULL,
142128d58a37SPierre Jolivet                                        MatConvert_MPISBAIJ_Basic,
1422f4259b30SLisandro Dalcin                                        NULL,
1423f4259b30SLisandro Dalcin                                        NULL,
1424f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
1425f4259b30SLisandro Dalcin                                        NULL,
1426f4259b30SLisandro Dalcin                                        NULL,
1427f4259b30SLisandro Dalcin                                        NULL,
1428f4259b30SLisandro Dalcin                                        NULL,
1429f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
1430f4259b30SLisandro Dalcin                                        NULL,
1431f4259b30SLisandro Dalcin                                        NULL,
143297304618SKris Buschelman                                        MatGetInertia_SeqSBAIJ,
14335bba2384SShri Abhyankar                                        MatLoad_SeqSBAIJ,
1434d519adbfSMatthew Knepley                                        /* 84*/ MatIsSymmetric_SeqSBAIJ,
1435865e5f61SKris Buschelman                                        MatIsHermitian_SeqSBAIJ,
1436efcf0fc3SBarry Smith                                        MatIsStructurallySymmetric_SeqSBAIJ,
1437f4259b30SLisandro Dalcin                                        NULL,
1438f4259b30SLisandro Dalcin                                        NULL,
1439f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
1440f4259b30SLisandro Dalcin                                        NULL,
1441f4259b30SLisandro Dalcin                                        NULL,
1442f4259b30SLisandro Dalcin                                        NULL,
1443f4259b30SLisandro Dalcin                                        NULL,
1444f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
1445f4259b30SLisandro Dalcin                                        NULL,
1446f4259b30SLisandro Dalcin                                        NULL,
1447f4259b30SLisandro Dalcin                                        NULL,
1448f4259b30SLisandro Dalcin                                        NULL,
1449f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
1450f4259b30SLisandro Dalcin                                        NULL,
1451f4259b30SLisandro Dalcin                                        NULL,
14522726fb6dSPierre Jolivet                                        MatConjugate_SeqSBAIJ,
1453f4259b30SLisandro Dalcin                                        NULL,
1454f4259b30SLisandro Dalcin                                        /*104*/ NULL,
145599cafbc1SBarry Smith                                        MatRealPart_SeqSBAIJ,
1456f5edf698SHong Zhang                                        MatImaginaryPart_SeqSBAIJ,
1457f5edf698SHong Zhang                                        MatGetRowUpperTriangular_SeqSBAIJ,
14582af78befSBarry Smith                                        MatRestoreRowUpperTriangular_SeqSBAIJ,
1459f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1460f4259b30SLisandro Dalcin                                        NULL,
1461f4259b30SLisandro Dalcin                                        NULL,
1462f4259b30SLisandro Dalcin                                        NULL,
1463547795f9SHong Zhang                                        MatMissingDiagonal_SeqSBAIJ,
1464f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1465f4259b30SLisandro Dalcin                                        NULL,
1466f4259b30SLisandro Dalcin                                        NULL,
1467f4259b30SLisandro Dalcin                                        NULL,
1468f4259b30SLisandro Dalcin                                        NULL,
1469f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1470f4259b30SLisandro Dalcin                                        NULL,
1471f4259b30SLisandro Dalcin                                        NULL,
1472f4259b30SLisandro Dalcin                                        NULL,
1473f4259b30SLisandro Dalcin                                        NULL,
1474f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1475f4259b30SLisandro Dalcin                                        NULL,
1476f4259b30SLisandro Dalcin                                        NULL,
1477f4259b30SLisandro Dalcin                                        NULL,
1478f4259b30SLisandro Dalcin                                        NULL,
1479f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1480f4259b30SLisandro Dalcin                                        NULL,
1481f4259b30SLisandro Dalcin                                        NULL,
1482f4259b30SLisandro Dalcin                                        NULL,
1483f4259b30SLisandro Dalcin                                        NULL,
1484f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1485f4259b30SLisandro Dalcin                                        NULL,
1486f4259b30SLisandro Dalcin                                        NULL,
1487f4259b30SLisandro Dalcin                                        NULL,
1488f4259b30SLisandro Dalcin                                        NULL,
148946533700Sstefano_zampini                                        /*139*/ MatSetBlockSizes_Default,
1490f4259b30SLisandro Dalcin                                        NULL,
1491f4259b30SLisandro Dalcin                                        NULL,
1492f4259b30SLisandro Dalcin                                        NULL,
1493f4259b30SLisandro Dalcin                                        NULL,
1494d70f29a3SPierre Jolivet                                        /*144*/ MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ,
1495d70f29a3SPierre Jolivet                                        NULL,
1496d70f29a3SPierre Jolivet                                        NULL,
149799a7f59eSMark Adams                                        NULL,
149899a7f59eSMark Adams                                        NULL,
14997fb60732SBarry Smith                                        NULL,
1500dec0b466SHong Zhang                                        /*150*/ NULL,
1501dec0b466SHong Zhang                                        NULL};
1502be1d678aSKris Buschelman 
1503d71ae5a4SJacob Faibussowitsch PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
1504d71ae5a4SJacob Faibussowitsch {
15054afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1506d0f46423SBarry Smith   PetscInt      nz  = aij->i[mat->rmap->N] * mat->rmap->bs * aij->bs2;
150749b5e25fSSatish Balay 
150849b5e25fSSatish Balay   PetscFunctionBegin;
150908401ef6SPierre Jolivet   PetscCheck(aij->nonew == 1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
151049b5e25fSSatish Balay 
151149b5e25fSSatish Balay   /* allocate space for values if not already there */
151248a46eb9SPierre Jolivet   if (!aij->saved_values) PetscCall(PetscMalloc1(nz + 1, &aij->saved_values));
151349b5e25fSSatish Balay 
151449b5e25fSSatish Balay   /* copy values over */
15159566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(aij->saved_values, aij->a, nz));
15163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151749b5e25fSSatish Balay }
151849b5e25fSSatish Balay 
1519d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
1520d71ae5a4SJacob Faibussowitsch {
15214afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1522d0f46423SBarry Smith   PetscInt      nz  = aij->i[mat->rmap->N] * mat->rmap->bs * aij->bs2;
152349b5e25fSSatish Balay 
152449b5e25fSSatish Balay   PetscFunctionBegin;
152508401ef6SPierre Jolivet   PetscCheck(aij->nonew == 1, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
152628b400f6SJacob Faibussowitsch   PetscCheck(aij->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
152749b5e25fSSatish Balay 
152849b5e25fSSatish Balay   /* copy values over */
15299566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(aij->a, aij->saved_values, nz));
15303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153149b5e25fSSatish Balay }
153249b5e25fSSatish Balay 
1533d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B, PetscInt bs, PetscInt nz, PetscInt *nnz)
1534d71ae5a4SJacob Faibussowitsch {
1535c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)B->data;
15364dcd73b1SHong Zhang   PetscInt      i, mbs, nbs, bs2;
15372576faa2SJed Brown   PetscBool     skipallocation = PETSC_FALSE, flg = PETSC_FALSE, realalloc = PETSC_FALSE;
153849b5e25fSSatish Balay 
1539b4e2f619SBarry Smith   PetscFunctionBegin;
1540ad79cf63SBarry Smith   if (B->hash_active) {
1541ad79cf63SBarry Smith     PetscInt bs;
1542aea10558SJacob Faibussowitsch     B->ops[0] = b->cops;
1543ad79cf63SBarry Smith     PetscCall(PetscHMapIJVDestroy(&b->ht));
1544ad79cf63SBarry Smith     PetscCall(MatGetBlockSize(B, &bs));
1545ad79cf63SBarry Smith     if (bs > 1) PetscCall(PetscHSetIJDestroy(&b->bht));
1546ad79cf63SBarry Smith     PetscCall(PetscFree(b->dnz));
1547ad79cf63SBarry Smith     PetscCall(PetscFree(b->bdnz));
1548ad79cf63SBarry Smith     B->hash_active = PETSC_FALSE;
1549ad79cf63SBarry Smith   }
15502576faa2SJed Brown   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1551db4efbfdSBarry Smith 
15529566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(B, PetscAbs(bs)));
15539566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
15549566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
155508401ef6SPierre Jolivet   PetscCheck(B->rmap->N <= B->cmap->N, PETSC_COMM_SELF, PETSC_ERR_SUP, "SEQSBAIJ matrix cannot have more rows %" PetscInt_FMT " than columns %" PetscInt_FMT, B->rmap->N, B->cmap->N);
15569566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(B->rmap, &bs));
1557899cda47SBarry Smith 
155821940c7eSstefano_zampini   B->preallocated = PETSC_TRUE;
155921940c7eSstefano_zampini 
1560d0f46423SBarry Smith   mbs = B->rmap->N / bs;
15614dcd73b1SHong Zhang   nbs = B->cmap->n / bs;
156249b5e25fSSatish Balay   bs2 = bs * bs;
156349b5e25fSSatish Balay 
1564aed4548fSBarry Smith   PetscCheck(mbs * bs == B->rmap->N && nbs * bs == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number rows, cols must be divisible by blocksize");
156549b5e25fSSatish Balay 
1566ab93d7beSBarry Smith   if (nz == MAT_SKIP_ALLOCATION) {
1567ab93d7beSBarry Smith     skipallocation = PETSC_TRUE;
1568ab93d7beSBarry Smith     nz             = 0;
1569ab93d7beSBarry Smith   }
1570ab93d7beSBarry Smith 
1571435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
157208401ef6SPierre Jolivet   PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nz cannot be less than 0: value %" PetscInt_FMT, nz);
157349b5e25fSSatish Balay   if (nnz) {
157449b5e25fSSatish Balay     for (i = 0; i < mbs; i++) {
157508401ef6SPierre 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]);
157608401ef6SPierre Jolivet       PetscCheck(nnz[i] <= nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nnz cannot be greater than block row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " block rowlength %" PetscInt_FMT, i, nnz[i], nbs);
157749b5e25fSSatish Balay     }
157849b5e25fSSatish Balay   }
157949b5e25fSSatish Balay 
1580db4efbfdSBarry Smith   B->ops->mult             = MatMult_SeqSBAIJ_N;
1581db4efbfdSBarry Smith   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1582db4efbfdSBarry Smith   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1583db4efbfdSBarry Smith   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;
158426fbe8dcSKarl Rupp 
15859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_no_unroll", &flg, NULL));
158649b5e25fSSatish Balay   if (!flg) {
158749b5e25fSSatish Balay     switch (bs) {
158849b5e25fSSatish Balay     case 1:
158949b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_1;
159049b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1591431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1592431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
159349b5e25fSSatish Balay       break;
159449b5e25fSSatish Balay     case 2:
159549b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_2;
159649b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1597431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1598431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
159949b5e25fSSatish Balay       break;
160049b5e25fSSatish Balay     case 3:
160149b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_3;
160249b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1603431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1604431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
160549b5e25fSSatish Balay       break;
160649b5e25fSSatish Balay     case 4:
160749b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_4;
160849b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1609431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1610431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
161149b5e25fSSatish Balay       break;
161249b5e25fSSatish Balay     case 5:
161349b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_5;
161449b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1615431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1616431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
161749b5e25fSSatish Balay       break;
161849b5e25fSSatish Balay     case 6:
161949b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_6;
162049b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1621431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1622431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
162349b5e25fSSatish Balay       break;
162449b5e25fSSatish Balay     case 7:
1625de53e5efSHong Zhang       B->ops->mult             = MatMult_SeqSBAIJ_7;
162649b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1627431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1628431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
162949b5e25fSSatish Balay       break;
163049b5e25fSSatish Balay     }
163149b5e25fSSatish Balay   }
163249b5e25fSSatish Balay 
163349b5e25fSSatish Balay   b->mbs = mbs;
16344dcd73b1SHong Zhang   b->nbs = nbs;
1635ab93d7beSBarry Smith   if (!skipallocation) {
16362ee49352SLisandro Dalcin     if (!b->imax) {
16379566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(mbs, &b->imax, mbs, &b->ilen));
163826fbe8dcSKarl Rupp 
1639c760cd28SBarry Smith       b->free_imax_ilen = PETSC_TRUE;
16402ee49352SLisandro Dalcin     }
164149b5e25fSSatish Balay     if (!nnz) {
1642435da068SBarry Smith       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
164349b5e25fSSatish Balay       else if (nz <= 0) nz = 1;
16445d2a9ed1SStefano Zampini       nz = PetscMin(nbs, nz);
164526fbe8dcSKarl Rupp       for (i = 0; i < mbs; i++) b->imax[i] = nz;
16469566063dSJacob Faibussowitsch       PetscCall(PetscIntMultError(nz, mbs, &nz));
164749b5e25fSSatish Balay     } else {
1648c73702f5SBarry Smith       PetscInt64 nz64 = 0;
16499371c9d4SSatish Balay       for (i = 0; i < mbs; i++) {
16509371c9d4SSatish Balay         b->imax[i] = nnz[i];
16519371c9d4SSatish Balay         nz64 += nnz[i];
16529371c9d4SSatish Balay       }
16539566063dSJacob Faibussowitsch       PetscCall(PetscIntCast(nz64, &nz));
165449b5e25fSSatish Balay     }
16552ee49352SLisandro Dalcin     /* b->ilen will count nonzeros in each block row so far. */
165626fbe8dcSKarl Rupp     for (i = 0; i < mbs; i++) b->ilen[i] = 0;
16576c6c5352SBarry Smith     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
165849b5e25fSSatish Balay 
165949b5e25fSSatish Balay     /* allocate the matrix space */
16609566063dSJacob Faibussowitsch     PetscCall(MatSeqXAIJFreeAIJ(B, &b->a, &b->j, &b->i));
16619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(bs2 * nz, &b->a, nz, &b->j, B->rmap->N + 1, &b->i));
16629566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(b->a, nz * bs2));
16639566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(b->j, nz));
166426fbe8dcSKarl Rupp 
166549b5e25fSSatish Balay     b->singlemalloc = PETSC_TRUE;
166649b5e25fSSatish Balay 
166749b5e25fSSatish Balay     /* pointer to beginning of each row */
1668e60cf9a0SBarry Smith     b->i[0] = 0;
166926fbe8dcSKarl Rupp     for (i = 1; i < mbs + 1; i++) b->i[i] = b->i[i - 1] + b->imax[i - 1];
167026fbe8dcSKarl Rupp 
1671e6b907acSBarry Smith     b->free_a  = PETSC_TRUE;
1672e6b907acSBarry Smith     b->free_ij = PETSC_TRUE;
1673e811da20SHong Zhang   } else {
1674e6b907acSBarry Smith     b->free_a  = PETSC_FALSE;
1675e6b907acSBarry Smith     b->free_ij = PETSC_FALSE;
1676ab93d7beSBarry Smith   }
167749b5e25fSSatish Balay 
167849b5e25fSSatish Balay   b->bs2     = bs2;
16796c6c5352SBarry Smith   b->nz      = 0;
1680b32cb4a7SJed Brown   b->maxnz   = nz;
1681f4259b30SLisandro Dalcin   b->inew    = NULL;
1682f4259b30SLisandro Dalcin   b->jnew    = NULL;
1683f4259b30SLisandro Dalcin   b->anew    = NULL;
1684f4259b30SLisandro Dalcin   b->a2anew  = NULL;
16851a3463dfSHong Zhang   b->permute = PETSC_FALSE;
1686cb7b82ddSBarry Smith 
1687cb7b82ddSBarry Smith   B->was_assembled = PETSC_FALSE;
1688cb7b82ddSBarry Smith   B->assembled     = PETSC_FALSE;
16899566063dSJacob Faibussowitsch   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
16903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1691c464158bSHong Zhang }
1692153ea458SHong Zhang 
1693d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B, PetscInt bs, const PetscInt ii[], const PetscInt jj[], const PetscScalar V[])
1694d71ae5a4SJacob Faibussowitsch {
16950cd7f59aSBarry Smith   PetscInt     i, j, m, nz, anz, nz_max = 0, *nnz;
1696f4259b30SLisandro Dalcin   PetscScalar *values      = NULL;
169738f409ebSLisandro Dalcin   PetscBool    roworiented = ((Mat_SeqSBAIJ *)B->data)->roworiented;
16980cd7f59aSBarry Smith 
169938f409ebSLisandro Dalcin   PetscFunctionBegin;
170008401ef6SPierre Jolivet   PetscCheck(bs >= 1, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_OUTOFRANGE, "Invalid block size specified, must be positive but it is %" PetscInt_FMT, bs);
17019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(B->rmap, bs));
17029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(B->cmap, bs));
17039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
17049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
17059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(B->rmap, &bs));
170638f409ebSLisandro Dalcin   m = B->rmap->n / bs;
170738f409ebSLisandro Dalcin 
1708aed4548fSBarry Smith   PetscCheck(!ii[0], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ii[0] must be 0 but it is %" PetscInt_FMT, ii[0]);
17099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m + 1, &nnz));
171038f409ebSLisandro Dalcin   for (i = 0; i < m; i++) {
171138f409ebSLisandro Dalcin     nz = ii[i + 1] - ii[i];
171208401ef6SPierre Jolivet     PetscCheck(nz >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT, i, nz);
17130cd7f59aSBarry Smith     anz = 0;
17140cd7f59aSBarry Smith     for (j = 0; j < nz; j++) {
17150cd7f59aSBarry Smith       /* count only values on the diagonal or above */
17160cd7f59aSBarry Smith       if (jj[ii[i] + j] >= i) {
17170cd7f59aSBarry Smith         anz = nz - j;
17180cd7f59aSBarry Smith         break;
17190cd7f59aSBarry Smith       }
17200cd7f59aSBarry Smith     }
17210cd7f59aSBarry Smith     nz_max = PetscMax(nz_max, anz);
17220cd7f59aSBarry Smith     nnz[i] = anz;
172338f409ebSLisandro Dalcin   }
17249566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(B, bs, 0, nnz));
17259566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
172638f409ebSLisandro Dalcin 
172738f409ebSLisandro Dalcin   values = (PetscScalar *)V;
172848a46eb9SPierre Jolivet   if (!values) PetscCall(PetscCalloc1(bs * bs * nz_max, &values));
172938f409ebSLisandro Dalcin   for (i = 0; i < m; i++) {
173038f409ebSLisandro Dalcin     PetscInt        ncols = ii[i + 1] - ii[i];
173138f409ebSLisandro Dalcin     const PetscInt *icols = jj + ii[i];
173238f409ebSLisandro Dalcin     if (!roworiented || bs == 1) {
173338f409ebSLisandro Dalcin       const PetscScalar *svals = values + (V ? (bs * bs * ii[i]) : 0);
17349566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked_SeqSBAIJ(B, 1, &i, ncols, icols, svals, INSERT_VALUES));
173538f409ebSLisandro Dalcin     } else {
173638f409ebSLisandro Dalcin       for (j = 0; j < ncols; j++) {
173738f409ebSLisandro Dalcin         const PetscScalar *svals = values + (V ? (bs * bs * (ii[i] + j)) : 0);
17389566063dSJacob Faibussowitsch         PetscCall(MatSetValuesBlocked_SeqSBAIJ(B, 1, &i, 1, &icols[j], svals, INSERT_VALUES));
173938f409ebSLisandro Dalcin       }
174038f409ebSLisandro Dalcin     }
174138f409ebSLisandro Dalcin   }
17429566063dSJacob Faibussowitsch   if (!V) PetscCall(PetscFree(values));
17439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
17449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
17459566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
17463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
174738f409ebSLisandro Dalcin }
174838f409ebSLisandro Dalcin 
1749db4efbfdSBarry Smith /*
1750db4efbfdSBarry Smith    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1751db4efbfdSBarry Smith */
1752d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B, PetscBool natural)
1753d71ae5a4SJacob Faibussowitsch {
1754ace3abfcSBarry Smith   PetscBool flg = PETSC_FALSE;
1755db4efbfdSBarry Smith   PetscInt  bs  = B->rmap->bs;
1756db4efbfdSBarry Smith 
1757db4efbfdSBarry Smith   PetscFunctionBegin;
17589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_no_unroll", &flg, NULL));
1759db4efbfdSBarry Smith   if (flg) bs = 8;
1760db4efbfdSBarry Smith 
1761db4efbfdSBarry Smith   if (!natural) {
1762db4efbfdSBarry Smith     switch (bs) {
1763d71ae5a4SJacob Faibussowitsch     case 1:
1764d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1765d71ae5a4SJacob Faibussowitsch       break;
1766d71ae5a4SJacob Faibussowitsch     case 2:
1767d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1768d71ae5a4SJacob Faibussowitsch       break;
1769d71ae5a4SJacob Faibussowitsch     case 3:
1770d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1771d71ae5a4SJacob Faibussowitsch       break;
1772d71ae5a4SJacob Faibussowitsch     case 4:
1773d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1774d71ae5a4SJacob Faibussowitsch       break;
1775d71ae5a4SJacob Faibussowitsch     case 5:
1776d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1777d71ae5a4SJacob Faibussowitsch       break;
1778d71ae5a4SJacob Faibussowitsch     case 6:
1779d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1780d71ae5a4SJacob Faibussowitsch       break;
1781d71ae5a4SJacob Faibussowitsch     case 7:
1782d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1783d71ae5a4SJacob Faibussowitsch       break;
1784d71ae5a4SJacob Faibussowitsch     default:
1785d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1786d71ae5a4SJacob Faibussowitsch       break;
1787db4efbfdSBarry Smith     }
1788db4efbfdSBarry Smith   } else {
1789db4efbfdSBarry Smith     switch (bs) {
1790d71ae5a4SJacob Faibussowitsch     case 1:
1791d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1792d71ae5a4SJacob Faibussowitsch       break;
1793d71ae5a4SJacob Faibussowitsch     case 2:
1794d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1795d71ae5a4SJacob Faibussowitsch       break;
1796d71ae5a4SJacob Faibussowitsch     case 3:
1797d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1798d71ae5a4SJacob Faibussowitsch       break;
1799d71ae5a4SJacob Faibussowitsch     case 4:
1800d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1801d71ae5a4SJacob Faibussowitsch       break;
1802d71ae5a4SJacob Faibussowitsch     case 5:
1803d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1804d71ae5a4SJacob Faibussowitsch       break;
1805d71ae5a4SJacob Faibussowitsch     case 6:
1806d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1807d71ae5a4SJacob Faibussowitsch       break;
1808d71ae5a4SJacob Faibussowitsch     case 7:
1809d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1810d71ae5a4SJacob Faibussowitsch       break;
1811d71ae5a4SJacob Faibussowitsch     default:
1812d71ae5a4SJacob Faibussowitsch       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1813d71ae5a4SJacob Faibussowitsch       break;
1814db4efbfdSBarry Smith     }
1815db4efbfdSBarry Smith   }
18163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1817db4efbfdSBarry Smith }
1818db4efbfdSBarry Smith 
1819cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType, MatReuse, Mat *);
1820cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *);
1821d71ae5a4SJacob Faibussowitsch static PetscErrorCode       MatFactorGetSolverType_petsc(Mat A, MatSolverType *type)
1822d71ae5a4SJacob Faibussowitsch {
18234ac6704cSBarry Smith   PetscFunctionBegin;
18244ac6704cSBarry Smith   *type = MATSOLVERPETSC;
18253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18264ac6704cSBarry Smith }
1827d769727bSBarry Smith 
1828d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A, MatFactorType ftype, Mat *B)
1829d71ae5a4SJacob Faibussowitsch {
1830d0f46423SBarry Smith   PetscInt n = A->rmap->n;
18315c9eb25fSBarry Smith 
18325c9eb25fSBarry Smith   PetscFunctionBegin;
18330e92d65fSHong Zhang #if defined(PETSC_USE_COMPLEX)
183403e5aca4SStefano Zampini   if ((ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) && A->hermitian == PETSC_BOOL3_TRUE && A->symmetric != PETSC_BOOL3_TRUE) {
183503e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Hermitian MAT_FACTOR_CHOLESKY or MAT_FACTOR_ICC are not supported. Use MAT_FACTOR_LU instead.\n"));
183603e5aca4SStefano Zampini     *B = NULL;
183703e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
183803e5aca4SStefano Zampini   }
18390e92d65fSHong Zhang #endif
1840eb1ec7c1SStefano Zampini 
18419566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
18429566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, n, n, n, n));
18435c9eb25fSBarry Smith   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
18449566063dSJacob Faibussowitsch     PetscCall(MatSetType(*B, MATSEQSBAIJ));
18459566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(*B, A->rmap->bs, MAT_SKIP_ALLOCATION, NULL));
184626fbe8dcSKarl Rupp 
18477b056e98SHong Zhang     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1848c6d0d4f0SHong Zhang     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
18499566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
18509566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC]));
1851e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
185200c67f3bSHong Zhang 
1853d5f3da31SBarry Smith   (*B)->factortype     = ftype;
1854f73b0415SBarry Smith   (*B)->canuseordering = PETSC_TRUE;
18559566063dSJacob Faibussowitsch   PetscCall(PetscFree((*B)->solvertype));
18569566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
18579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_petsc));
18583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18595c9eb25fSBarry Smith }
18605c9eb25fSBarry Smith 
18618397e458SBarry Smith /*@C
18622ef1f0ffSBarry Smith    MatSeqSBAIJGetArray - gives access to the array where the numerical data for a `MATSEQSBAIJ` matrix is stored
18638397e458SBarry Smith 
18648397e458SBarry Smith    Not Collective
18658397e458SBarry Smith 
18668397e458SBarry Smith    Input Parameter:
186711a5261eSBarry Smith .  mat - a `MATSEQSBAIJ` matrix
18688397e458SBarry Smith 
18698397e458SBarry Smith    Output Parameter:
18708397e458SBarry Smith .   array - pointer to the data
18718397e458SBarry Smith 
18728397e458SBarry Smith    Level: intermediate
18738397e458SBarry Smith 
18741cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatSeqSBAIJRestoreArray()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`
18758397e458SBarry Smith @*/
1876d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJGetArray(Mat A, PetscScalar **array)
1877d71ae5a4SJacob Faibussowitsch {
18788397e458SBarry Smith   PetscFunctionBegin;
1879cac4c232SBarry Smith   PetscUseMethod(A, "MatSeqSBAIJGetArray_C", (Mat, PetscScalar **), (A, array));
18803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18818397e458SBarry Smith }
18828397e458SBarry Smith 
18838397e458SBarry Smith /*@C
18842ef1f0ffSBarry Smith    MatSeqSBAIJRestoreArray - returns access to the array where the numerical data for a `MATSEQSBAIJ` matrix is stored obtained by `MatSeqSBAIJGetArray()`
18858397e458SBarry Smith 
18868397e458SBarry Smith    Not Collective
18878397e458SBarry Smith 
18888397e458SBarry Smith    Input Parameters:
18892ef1f0ffSBarry Smith +  mat - a `MATSEQSBAIJ` matrix
1890a2b725a8SWilliam Gropp -  array - pointer to the data
18918397e458SBarry Smith 
18928397e458SBarry Smith    Level: intermediate
18938397e458SBarry Smith 
18941cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatSeqSBAIJGetArray()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`
18958397e458SBarry Smith @*/
1896d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJRestoreArray(Mat A, PetscScalar **array)
1897d71ae5a4SJacob Faibussowitsch {
18988397e458SBarry Smith   PetscFunctionBegin;
1899cac4c232SBarry Smith   PetscUseMethod(A, "MatSeqSBAIJRestoreArray_C", (Mat, PetscScalar **), (A, array));
19003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19018397e458SBarry Smith }
19028397e458SBarry Smith 
19030bad9183SKris Buschelman /*MC
1904fafad747SKris Buschelman   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
19050bad9183SKris Buschelman   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
19060bad9183SKris Buschelman 
1907828413b8SBarry Smith   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
190811a5261eSBarry Smith   can call `MatSetOption`(`Mat`, `MAT_HERMITIAN`).
1909828413b8SBarry Smith 
19102ef1f0ffSBarry Smith   Options Database Key:
191111a5261eSBarry Smith   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to `MatSetFromOptions()`
19120bad9183SKris Buschelman 
19132ef1f0ffSBarry Smith   Level: beginner
19142ef1f0ffSBarry Smith 
191595452b02SPatrick Sanan   Notes:
191695452b02SPatrick Sanan     By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
191711a5261eSBarry Smith      stored and it is assumed they symmetric to the upper triangular). If you call `MatSetOption`(`Mat`,`MAT_IGNORE_LOWER_TRIANGULAR`,`PETSC_FALSE`) or use
19182ef1f0ffSBarry Smith      the options database `-mat_ignore_lower_triangular` false it will generate an error if you try to set a value in the lower triangular portion.
191971dad5bbSBarry Smith 
1920476417e5SBarry Smith     The number of rows in the matrix must be less than or equal to the number of columns
192171dad5bbSBarry Smith 
19221cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreateSeqSBAIJ()`, `MatType`, `MATMPISBAIJ`
19230bad9183SKris Buschelman M*/
1924d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1925d71ae5a4SJacob Faibussowitsch {
1926a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
192713f74950SBarry Smith   PetscMPIInt   size;
1928ace3abfcSBarry Smith   PetscBool     no_unroll = PETSC_FALSE, no_inode = PETSC_FALSE;
1929a23d5eceSKris Buschelman 
1930a23d5eceSKris Buschelman   PetscFunctionBegin;
19319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
193208401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Comm must be of size 1");
1933a23d5eceSKris Buschelman 
19344dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
1935a23d5eceSKris Buschelman   B->data   = (void *)b;
1936aea10558SJacob Faibussowitsch   B->ops[0] = MatOps_Values;
193726fbe8dcSKarl Rupp 
1938a23d5eceSKris Buschelman   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1939a23d5eceSKris Buschelman   B->ops->view       = MatView_SeqSBAIJ;
1940f4259b30SLisandro Dalcin   b->row             = NULL;
1941f4259b30SLisandro Dalcin   b->icol            = NULL;
1942a23d5eceSKris Buschelman   b->reallocs        = 0;
1943f4259b30SLisandro Dalcin   b->saved_values    = NULL;
19440def2e27SBarry Smith   b->inode.limit     = 5;
19450def2e27SBarry Smith   b->inode.max_limit = 5;
1946a23d5eceSKris Buschelman 
1947a23d5eceSKris Buschelman   b->roworiented        = PETSC_TRUE;
1948a23d5eceSKris Buschelman   b->nonew              = 0;
1949f4259b30SLisandro Dalcin   b->diag               = NULL;
1950f4259b30SLisandro Dalcin   b->solve_work         = NULL;
1951f4259b30SLisandro Dalcin   b->mult_work          = NULL;
1952f4259b30SLisandro Dalcin   B->spptr              = NULL;
1953f2cbd3d5SJed Brown   B->info.nz_unneeded   = (PetscReal)b->maxnz * b->bs2;
1954a9817697SBarry Smith   b->keepnonzeropattern = PETSC_FALSE;
1955a23d5eceSKris Buschelman 
1956f4259b30SLisandro Dalcin   b->inew    = NULL;
1957f4259b30SLisandro Dalcin   b->jnew    = NULL;
1958f4259b30SLisandro Dalcin   b->anew    = NULL;
1959f4259b30SLisandro Dalcin   b->a2anew  = NULL;
1960a23d5eceSKris Buschelman   b->permute = PETSC_FALSE;
1961a23d5eceSKris Buschelman 
196271dad5bbSBarry Smith   b->ignore_ltriangular = PETSC_TRUE;
196326fbe8dcSKarl Rupp 
19649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_ignore_lower_triangular", &b->ignore_ltriangular, NULL));
1965941593c8SHong Zhang 
1966f5edf698SHong Zhang   b->getrow_utriangular = PETSC_FALSE;
196726fbe8dcSKarl Rupp 
19689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)B)->options, ((PetscObject)B)->prefix, "-mat_getrow_uppertriangular", &b->getrow_utriangular, NULL));
1969f5edf698SHong Zhang 
19709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJGetArray_C", MatSeqSBAIJGetArray_SeqSBAIJ));
19719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJRestoreArray_C", MatSeqSBAIJRestoreArray_SeqSBAIJ));
19729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqSBAIJ));
19739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqSBAIJ));
19749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJSetColumnIndices_C", MatSeqSBAIJSetColumnIndices_SeqSBAIJ));
19759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_seqaij_C", MatConvert_SeqSBAIJ_SeqAIJ));
19769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_seqbaij_C", MatConvert_SeqSBAIJ_SeqBAIJ));
19779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJSetPreallocation_C", MatSeqSBAIJSetPreallocation_SeqSBAIJ));
19789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSBAIJSetPreallocationCSR_C", MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ));
19796214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
19809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_elemental_C", MatConvert_SeqSBAIJ_Elemental));
19816214f412SHong Zhang #endif
1982d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
19839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsbaij_scalapack_C", MatConvert_SBAIJ_ScaLAPACK));
1984d24d4204SJose E. Roman #endif
198523ce1328SBarry Smith 
1986b94d7dedSBarry Smith   B->symmetry_eternal            = PETSC_TRUE;
1987b94d7dedSBarry Smith   B->structural_symmetry_eternal = PETSC_TRUE;
1988b94d7dedSBarry Smith   B->symmetric                   = PETSC_BOOL3_TRUE;
1989b94d7dedSBarry Smith   B->structurally_symmetric      = PETSC_BOOL3_TRUE;
1990eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
1991b94d7dedSBarry Smith   B->hermitian = PETSC_BOOL3_FALSE;
1992eb1ec7c1SStefano Zampini #else
1993b94d7dedSBarry Smith   B->hermitian = PETSC_BOOL3_TRUE;
1994eb1ec7c1SStefano Zampini #endif
199513647f61SHong Zhang 
19969566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSBAIJ));
19970def2e27SBarry Smith 
1998d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "Options for SEQSBAIJ matrix", "Mat");
19999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_unroll", "Do not optimize for inodes (slower)", NULL, no_unroll, &no_unroll, NULL));
200048a46eb9SPierre Jolivet   if (no_unroll) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_unroll\n"));
20019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_inode", "Do not optimize for inodes (slower)", NULL, no_inode, &no_inode, NULL));
20029566063dSJacob Faibussowitsch   if (no_inode) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_inode\n"));
20039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_inode_limit", "Do not use inodes larger then this value", NULL, b->inode.limit, &b->inode.limit, NULL));
2004d0609cedSBarry Smith   PetscOptionsEnd();
2005ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
20060def2e27SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
20073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2008a23d5eceSKris Buschelman }
2009a23d5eceSKris Buschelman 
2010a23d5eceSKris Buschelman /*@C
2011a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
201211a5261eSBarry Smith    compressed row) `MATSEQSBAIJ` format.  For good matrix assembly performance the
201320f4b53cSBarry Smith    user should preallocate the matrix storage by setting the parameter `nz`
201420f4b53cSBarry Smith    (or the array `nnz`).
2015a23d5eceSKris Buschelman 
2016c3339decSBarry Smith    Collective
2017a23d5eceSKris Buschelman 
2018a23d5eceSKris Buschelman    Input Parameters:
20191c4f3114SJed Brown +  B - the symmetric matrix
202011a5261eSBarry Smith .  bs - size of block, the blocks are ALWAYS square. One can use `MatSetBlockSizes()` to set a different row and column blocksize but the row
202111a5261eSBarry Smith           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with `MatCreateVecs()`
2022a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
2023a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
20242ef1f0ffSBarry Smith          diagonal portion of each block (possibly different for each block row) or `NULL`
2025a23d5eceSKris Buschelman 
2026a23d5eceSKris Buschelman    Options Database Keys:
2027a2b725a8SWilliam Gropp +   -mat_no_unroll - uses code that does not unroll the loops in the
2028a23d5eceSKris Buschelman                      block calculations (much slower)
2029a2b725a8SWilliam Gropp -   -mat_block_size - size of the blocks to use (only works if a negative bs is passed in
2030a23d5eceSKris Buschelman 
2031a23d5eceSKris Buschelman    Level: intermediate
2032a23d5eceSKris Buschelman 
2033a23d5eceSKris Buschelman    Notes:
203420f4b53cSBarry Smith    Specify the preallocated storage with either `nz` or `nnz` (not both).
20352ef1f0ffSBarry Smith    Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
2036651615e1SBarry Smith    allocation.  See [Sparse Matrices](sec_matsparse) for details.
2037a23d5eceSKris Buschelman 
203811a5261eSBarry Smith    You can call `MatGetInfo()` to get information on how effective the preallocation was;
2039aa95bbe8SBarry Smith    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
20402ef1f0ffSBarry Smith    You can also run with the option `-info` and look for messages with the string
2041aa95bbe8SBarry Smith    malloc in them to see if additional memory allocation was needed.
2042aa95bbe8SBarry Smith 
20432ef1f0ffSBarry Smith    If the `nnz` parameter is given then the `nz` parameter is ignored
204449a6f317SBarry Smith 
20451cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, [Sparse Matrices](sec_matsparse), `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatCreateSBAIJ()`
2046a23d5eceSKris Buschelman @*/
2047d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B, PetscInt bs, PetscInt nz, const PetscInt nnz[])
2048d71ae5a4SJacob Faibussowitsch {
2049a23d5eceSKris Buschelman   PetscFunctionBegin;
20506ba663aaSJed Brown   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
20516ba663aaSJed Brown   PetscValidType(B, 1);
20526ba663aaSJed Brown   PetscValidLogicalCollectiveInt(B, bs, 2);
2053cac4c232SBarry Smith   PetscTryMethod(B, "MatSeqSBAIJSetPreallocation_C", (Mat, PetscInt, PetscInt, const PetscInt[]), (B, bs, nz, nnz));
20543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2055a23d5eceSKris Buschelman }
205649b5e25fSSatish Balay 
205738f409ebSLisandro Dalcin /*@C
205811a5261eSBarry Smith    MatSeqSBAIJSetPreallocationCSR - Creates a sparse parallel matrix in `MATSEQSBAIJ` format using the given nonzero structure and (optional) numerical values
205938f409ebSLisandro Dalcin 
206038f409ebSLisandro Dalcin    Input Parameters:
20611c4f3114SJed Brown +  B - the matrix
2062eab78319SHong Zhang .  bs - size of block, the blocks are ALWAYS square.
206338f409ebSLisandro Dalcin .  i - the indices into j for the start of each local row (starts with zero)
206438f409ebSLisandro Dalcin .  j - the column indices for each local row (starts with zero) these must be sorted for each row
206538f409ebSLisandro Dalcin -  v - optional values in the matrix
206638f409ebSLisandro Dalcin 
2067664954b6SBarry Smith    Level: advanced
206838f409ebSLisandro Dalcin 
206938f409ebSLisandro Dalcin    Notes:
207011a5261eSBarry Smith    The order of the entries in values is specified by the `MatOption` `MAT_ROW_ORIENTED`.  For example, C programs
207111a5261eSBarry Smith    may want to use the default `MAT_ROW_ORIENTED` = `PETSC_TRUE` and use an array v[nnz][bs][bs] where the second index is
207238f409ebSLisandro Dalcin    over rows within a block and the last index is over columns within a block row.  Fortran programs will likely set
207311a5261eSBarry Smith    `MAT_ROW_ORIENTED` = `PETSC_FALSE` and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
207438f409ebSLisandro Dalcin    block column and the second index is over columns within a block.
207538f409ebSLisandro Dalcin 
207650c5228eSBarry Smith    Any entries below the diagonal are ignored
20770cd7f59aSBarry Smith 
20780cd7f59aSBarry Smith    Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
20790cd7f59aSBarry Smith    and usually the numerical values as well
2080664954b6SBarry Smith 
20811cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqSBAIJ()`, `MatSetValuesBlocked()`, `MatSeqSBAIJSetPreallocation()`, `MATSEQSBAIJ`
208238f409ebSLisandro Dalcin @*/
2083d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B, PetscInt bs, const PetscInt i[], const PetscInt j[], const PetscScalar v[])
2084d71ae5a4SJacob Faibussowitsch {
208538f409ebSLisandro Dalcin   PetscFunctionBegin;
208638f409ebSLisandro Dalcin   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
208738f409ebSLisandro Dalcin   PetscValidType(B, 1);
208838f409ebSLisandro Dalcin   PetscValidLogicalCollectiveInt(B, bs, 2);
2089cac4c232SBarry Smith   PetscTryMethod(B, "MatSeqSBAIJSetPreallocationCSR_C", (Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]), (B, bs, i, j, v));
20903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
209138f409ebSLisandro Dalcin }
209238f409ebSLisandro Dalcin 
2093c464158bSHong Zhang /*@C
20942ef1f0ffSBarry Smith    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in (block
209511a5261eSBarry Smith    compressed row) `MATSEQSBAIJ` format.  For good matrix assembly performance the
20962ef1f0ffSBarry Smith    user should preallocate the matrix storage by setting the parameter `nz`
20972ef1f0ffSBarry Smith    (or the array `nnz`).
209849b5e25fSSatish Balay 
2099d083f849SBarry Smith    Collective
2100c464158bSHong Zhang 
2101c464158bSHong Zhang    Input Parameters:
210211a5261eSBarry Smith +  comm - MPI communicator, set to `PETSC_COMM_SELF`
210311a5261eSBarry Smith .  bs - size of block, the blocks are ALWAYS square. One can use `MatSetBlockSizes()` to set a different row and column blocksize but the row
2104bb7ae925SBarry Smith           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
210520f4b53cSBarry Smith .  m - number of rows
210620f4b53cSBarry Smith .  n - number of columns
2107c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
2108744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
21092ef1f0ffSBarry Smith          diagonal portion of each block (possibly different for each block row) or `NULL`
2110c464158bSHong Zhang 
2111c464158bSHong Zhang    Output Parameter:
2112c464158bSHong Zhang .  A - the symmetric matrix
2113c464158bSHong Zhang 
2114c464158bSHong Zhang    Options Database Keys:
2115a2b725a8SWilliam Gropp +   -mat_no_unroll - uses code that does not unroll the loops in the
2116c464158bSHong Zhang                      block calculations (much slower)
2117a2b725a8SWilliam Gropp -   -mat_block_size - size of the blocks to use
2118c464158bSHong Zhang 
2119c464158bSHong Zhang    Level: intermediate
2120c464158bSHong Zhang 
21212ef1f0ffSBarry Smith    Notes:
212211a5261eSBarry Smith    It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
2123f6f02116SRichard Tran Mills    MatXXXXSetPreallocation() paradigm instead of this routine directly.
212411a5261eSBarry Smith    [MatXXXXSetPreallocation() is, for example, `MatSeqAIJSetPreallocation()`]
2125175b88e8SBarry Smith 
21266d6d819aSHong Zhang    The number of rows and columns must be divisible by blocksize.
21276d6d819aSHong Zhang    This matrix type does not support complex Hermitian operation.
2128c464158bSHong Zhang 
21292ef1f0ffSBarry Smith    Specify the preallocated storage with either `nz` or `nnz` (not both).
21302ef1f0ffSBarry Smith    Set `nz` = `PETSC_DEFAULT` and `nnz` = `NULL` for PETSc to control dynamic memory
2131651615e1SBarry Smith    allocation.  See [Sparse Matrices](sec_matsparse) for details.
2132c464158bSHong Zhang 
21332ef1f0ffSBarry Smith    If the `nnz` parameter is given then the `nz` parameter is ignored
213449a6f317SBarry Smith 
21351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, [Sparse Matrices](sec_matsparse), `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatCreateSBAIJ()`
2136c464158bSHong Zhang @*/
2137d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
2138d71ae5a4SJacob Faibussowitsch {
2139c464158bSHong Zhang   PetscFunctionBegin;
21409566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
21419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, m, n));
21429566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATSEQSBAIJ));
21439566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(*A, bs, nz, (PetscInt *)nnz));
21443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
214549b5e25fSSatish Balay }
214649b5e25fSSatish Balay 
2147d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A, MatDuplicateOption cpvalues, Mat *B)
2148d71ae5a4SJacob Faibussowitsch {
214949b5e25fSSatish Balay   Mat           C;
215049b5e25fSSatish Balay   Mat_SeqSBAIJ *c, *a  = (Mat_SeqSBAIJ *)A->data;
2151b40805acSSatish Balay   PetscInt      i, mbs = a->mbs, nz = a->nz, bs2 = a->bs2;
215249b5e25fSSatish Balay 
215349b5e25fSSatish Balay   PetscFunctionBegin;
215431fe6a7dSBarry Smith   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Cannot duplicate unassembled matrix");
215508401ef6SPierre Jolivet   PetscCheck(a->i[mbs] == nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Corrupt matrix");
215649b5e25fSSatish Balay 
2157f4259b30SLisandro Dalcin   *B = NULL;
21589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
21599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, A->rmap->N, A->cmap->n, A->rmap->N, A->cmap->n));
21609566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(C, A, A));
21619566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATSEQSBAIJ));
2162692f9cbeSHong Zhang   c = (Mat_SeqSBAIJ *)C->data;
2163692f9cbeSHong Zhang 
2164273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
2165d5f3da31SBarry Smith   C->factortype         = A->factortype;
2166f4259b30SLisandro Dalcin   c->row                = NULL;
2167f4259b30SLisandro Dalcin   c->icol               = NULL;
2168f4259b30SLisandro Dalcin   c->saved_values       = NULL;
2169a9817697SBarry Smith   c->keepnonzeropattern = a->keepnonzeropattern;
217049b5e25fSSatish Balay   C->assembled          = PETSC_TRUE;
217149b5e25fSSatish Balay 
21729566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
21739566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &C->cmap));
217449b5e25fSSatish Balay   c->bs2 = a->bs2;
217549b5e25fSSatish Balay   c->mbs = a->mbs;
217649b5e25fSSatish Balay   c->nbs = a->nbs;
217749b5e25fSSatish Balay 
2178c760cd28SBarry Smith   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2179c760cd28SBarry Smith     c->imax           = a->imax;
2180c760cd28SBarry Smith     c->ilen           = a->ilen;
2181c760cd28SBarry Smith     c->free_imax_ilen = PETSC_FALSE;
2182c760cd28SBarry Smith   } else {
21839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2((mbs + 1), &c->imax, (mbs + 1), &c->ilen));
218449b5e25fSSatish Balay     for (i = 0; i < mbs; i++) {
218549b5e25fSSatish Balay       c->imax[i] = a->imax[i];
218649b5e25fSSatish Balay       c->ilen[i] = a->ilen[i];
218749b5e25fSSatish Balay     }
2188c760cd28SBarry Smith     c->free_imax_ilen = PETSC_TRUE;
2189c760cd28SBarry Smith   }
219049b5e25fSSatish Balay 
219149b5e25fSSatish Balay   /* allocate the matrix space */
21924da8f245SBarry Smith   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
21939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs2 * nz, &c->a));
219444e1c64aSLisandro Dalcin     c->i            = a->i;
219544e1c64aSLisandro Dalcin     c->j            = a->j;
21964da8f245SBarry Smith     c->singlemalloc = PETSC_FALSE;
219744e1c64aSLisandro Dalcin     c->free_a       = PETSC_TRUE;
21984da8f245SBarry Smith     c->free_ij      = PETSC_FALSE;
21994da8f245SBarry Smith     c->parent       = A;
22009566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
22019566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
22029566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
22034da8f245SBarry Smith   } else {
22049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(bs2 * nz, &c->a, nz, &c->j, mbs + 1, &c->i));
22059566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->i, a->i, mbs + 1));
22064da8f245SBarry Smith     c->singlemalloc = PETSC_TRUE;
220744e1c64aSLisandro Dalcin     c->free_a       = PETSC_TRUE;
22084da8f245SBarry Smith     c->free_ij      = PETSC_TRUE;
22094da8f245SBarry Smith   }
221049b5e25fSSatish Balay   if (mbs > 0) {
221148a46eb9SPierre Jolivet     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) PetscCall(PetscArraycpy(c->j, a->j, nz));
221249b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
22139566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->a, a->a, bs2 * nz));
221449b5e25fSSatish Balay     } else {
22159566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(c->a, bs2 * nz));
221649b5e25fSSatish Balay     }
2217a1c3900fSBarry Smith     if (a->jshort) {
221844e1c64aSLisandro Dalcin       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
221944e1c64aSLisandro Dalcin       /* if the parent matrix is reassembled, this child matrix will never notice */
22209566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nz, &c->jshort));
22219566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->jshort, a->jshort, nz));
222226fbe8dcSKarl Rupp 
22234da8f245SBarry Smith       c->free_jshort = PETSC_TRUE;
22244da8f245SBarry Smith     }
2225a1c3900fSBarry Smith   }
222649b5e25fSSatish Balay 
222749b5e25fSSatish Balay   c->roworiented = a->roworiented;
222849b5e25fSSatish Balay   c->nonew       = a->nonew;
222949b5e25fSSatish Balay 
223049b5e25fSSatish Balay   if (a->diag) {
2231c760cd28SBarry Smith     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2232c760cd28SBarry Smith       c->diag      = a->diag;
2233c760cd28SBarry Smith       c->free_diag = PETSC_FALSE;
2234c760cd28SBarry Smith     } else {
22359566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mbs, &c->diag));
223626fbe8dcSKarl Rupp       for (i = 0; i < mbs; i++) c->diag[i] = a->diag[i];
2237c760cd28SBarry Smith       c->free_diag = PETSC_TRUE;
2238c760cd28SBarry Smith     }
223944e1c64aSLisandro Dalcin   }
22406c6c5352SBarry Smith   c->nz         = a->nz;
2241f2cbd3d5SJed Brown   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
2242f4259b30SLisandro Dalcin   c->solve_work = NULL;
2243f4259b30SLisandro Dalcin   c->mult_work  = NULL;
224426fbe8dcSKarl Rupp 
224549b5e25fSSatish Balay   *B = C;
22469566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
22473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224849b5e25fSSatish Balay }
224949b5e25fSSatish Balay 
2250618cc2edSLisandro Dalcin /* Used for both SeqBAIJ and SeqSBAIJ matrices */
2251618cc2edSLisandro Dalcin #define MatLoad_SeqSBAIJ_Binary MatLoad_SeqBAIJ_Binary
2252618cc2edSLisandro Dalcin 
2253d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLoad_SeqSBAIJ(Mat mat, PetscViewer viewer)
2254d71ae5a4SJacob Faibussowitsch {
22557f489da9SVaclav Hapla   PetscBool isbinary;
22562f480046SShri Abhyankar 
22572f480046SShri Abhyankar   PetscFunctionBegin;
22589566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
225928b400f6SJacob Faibussowitsch   PetscCheck(isbinary, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer type %s not yet supported for reading %s matrices", ((PetscObject)viewer)->type_name, ((PetscObject)mat)->type_name);
22609566063dSJacob Faibussowitsch   PetscCall(MatLoad_SeqSBAIJ_Binary(mat, viewer));
22613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22622f480046SShri Abhyankar }
22632f480046SShri Abhyankar 
2264c75a6043SHong Zhang /*@
226511a5261eSBarry Smith      MatCreateSeqSBAIJWithArrays - Creates an sequential `MATSEQSBAIJ` matrix using matrix elements
2266c75a6043SHong Zhang               (upper triangular entries in CSR format) provided by the user.
2267c75a6043SHong Zhang 
2268d083f849SBarry Smith      Collective
2269c75a6043SHong Zhang 
2270c75a6043SHong Zhang    Input Parameters:
2271c75a6043SHong Zhang +  comm - must be an MPI communicator of size 1
2272c75a6043SHong Zhang .  bs - size of block
2273c75a6043SHong Zhang .  m - number of rows
2274c75a6043SHong Zhang .  n - number of columns
2275483a2f95SBarry Smith .  i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix
2276c75a6043SHong Zhang .  j - column indices
2277c75a6043SHong Zhang -  a - matrix values
2278c75a6043SHong Zhang 
2279c75a6043SHong Zhang    Output Parameter:
2280c75a6043SHong Zhang .  mat - the matrix
2281c75a6043SHong Zhang 
2282dfb205c3SBarry Smith    Level: advanced
2283c75a6043SHong Zhang 
2284c75a6043SHong Zhang    Notes:
22852ef1f0ffSBarry Smith        The `i`, `j`, and `a` arrays are not copied by this routine, the user must free these arrays
2286c75a6043SHong Zhang     once the matrix is destroyed
2287c75a6043SHong Zhang 
2288c75a6043SHong Zhang        You cannot set new nonzero locations into this matrix, that will generate an error.
2289c75a6043SHong Zhang 
22902ef1f0ffSBarry Smith        The `i` and `j` indices are 0 based
2291c75a6043SHong Zhang 
22922ef1f0ffSBarry Smith        When block size is greater than 1 the matrix values must be stored using the `MATSBAIJ` storage format. For block size of 1
2293dfb205c3SBarry Smith        it is the regular CSR format excluding the lower triangular elements.
2294dfb205c3SBarry Smith 
22951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSEQSBAIJ`, `MatCreate()`, `MatCreateSBAIJ()`, `MatCreateSeqSBAIJ()`
2296c75a6043SHong Zhang @*/
2297d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm comm, PetscInt bs, PetscInt m, PetscInt n, PetscInt i[], PetscInt j[], PetscScalar a[], Mat *mat)
2298d71ae5a4SJacob Faibussowitsch {
2299c75a6043SHong Zhang   PetscInt      ii;
2300c75a6043SHong Zhang   Mat_SeqSBAIJ *sbaij;
2301c75a6043SHong Zhang 
2302c75a6043SHong Zhang   PetscFunctionBegin;
230308401ef6SPierre Jolivet   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "block size %" PetscInt_FMT " > 1 is not supported yet", bs);
2304aed4548fSBarry Smith   PetscCheck(m == 0 || i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i (row indices) must start with 0");
2305c75a6043SHong Zhang 
23069566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, mat));
23079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mat, m, n, m, n));
23089566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mat, MATSEQSBAIJ));
23099566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(*mat, bs, MAT_SKIP_ALLOCATION, NULL));
2310c75a6043SHong Zhang   sbaij = (Mat_SeqSBAIJ *)(*mat)->data;
23119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m, &sbaij->imax, m, &sbaij->ilen));
2312c75a6043SHong Zhang 
2313c75a6043SHong Zhang   sbaij->i = i;
2314c75a6043SHong Zhang   sbaij->j = j;
2315c75a6043SHong Zhang   sbaij->a = a;
231626fbe8dcSKarl Rupp 
2317c75a6043SHong Zhang   sbaij->singlemalloc   = PETSC_FALSE;
2318c75a6043SHong Zhang   sbaij->nonew          = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2319e6b907acSBarry Smith   sbaij->free_a         = PETSC_FALSE;
2320e6b907acSBarry Smith   sbaij->free_ij        = PETSC_FALSE;
2321ddf7884eSMatthew Knepley   sbaij->free_imax_ilen = PETSC_TRUE;
2322c75a6043SHong Zhang 
2323c75a6043SHong Zhang   for (ii = 0; ii < m; ii++) {
2324c75a6043SHong Zhang     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii + 1] - i[ii];
23256bdcaf15SBarry Smith     PetscCheck(i[ii + 1] >= i[ii], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative row length in i (row indices) row = %" PetscInt_FMT " length = %" PetscInt_FMT, ii, i[ii + 1] - i[ii]);
2326c75a6043SHong Zhang   }
232776bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
2328c75a6043SHong Zhang     for (ii = 0; ii < sbaij->i[m]; ii++) {
23296bdcaf15SBarry Smith       PetscCheck(j[ii] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative column index at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
23306bdcaf15SBarry Smith       PetscCheck(j[ii] < n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index too large at location = %" PetscInt_FMT " index = %" PetscInt_FMT, ii, j[ii]);
2331c75a6043SHong Zhang     }
233276bd3646SJed Brown   }
2333c75a6043SHong Zhang 
23349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY));
23359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY));
23363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2337c75a6043SHong Zhang }
2338d06b337dSHong Zhang 
2339d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm comm, Mat inmat, PetscInt n, MatReuse scall, Mat *outmat)
2340d71ae5a4SJacob Faibussowitsch {
234159f5e6ceSHong Zhang   PetscFunctionBegin;
23429566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm, inmat, n, scall, outmat));
23433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234459f5e6ceSHong Zhang }
2345