xref: /petsc/src/mat/impls/sbaij/seq/sbaij.h (revision e8c0849ab8fe171bed529bea27238c9b402db591)
1a4963045SJacob Faibussowitsch #pragma once
2af0996ceSBarry Smith #include <petsc/private/matimpl.h>
3c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
449b5e25fSSatish Balay 
549b5e25fSSatish Balay /*
649b5e25fSSatish Balay   MATSEQSBAIJ format - Block compressed row storage. The i[] and j[]
749b5e25fSSatish Balay   arrays start at 0.
849b5e25fSSatish Balay */
949b5e25fSSatish Balay 
1049b5e25fSSatish Balay typedef struct {
11dd6ea824SBarry Smith   SEQAIJHEADER(MatScalar);
12e6b907acSBarry Smith   SEQBAIJHEADER;
1313f74950SBarry Smith   PetscInt        *inew;               /* pointer to beginning of each row of reordered matrix */
1413f74950SBarry Smith   PetscInt        *jnew;               /* column values: jnew + i[k] is start of row k */
15a6175056SHong Zhang   MatScalar       *anew;               /* nonzero diagonal and superdiagonal elements of reordered matrix */
16d59c15a7SBarry Smith   PetscScalar     *solves_work;        /* work space used in MatSolves */
1713f74950SBarry Smith   PetscInt         solves_work_n;      /* size of solves_work */
1813f74950SBarry Smith   PetscInt        *a2anew;             /* map used for symm permutation */
19ace3abfcSBarry Smith   PetscBool        permute;            /* if true, a non-trivial permutation is used for factorization */
20ace3abfcSBarry Smith   PetscBool        ignore_ltriangular; /* if true, ignore the lower triangular values inserted by users */
21ace3abfcSBarry Smith   PetscBool        getrow_utriangular; /* if true, MatGetRow_SeqSBAIJ() is enabled to get the upper part of the row */
224108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
2338702af4SBarry Smith   unsigned short  *jshort;
24ace3abfcSBarry Smith   PetscBool        free_jshort;
2549b5e25fSSatish Balay } Mat_SeqSBAIJ;
2649b5e25fSSatish Balay 
275a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat, Mat, IS, const MatFactorInfo *);
285a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
295a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat, IS, const MatFactorInfo *);
305a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat, Mat, IS, const MatFactorInfo *);
315a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
325a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat, MatDuplicateOption, Mat *);
33*421480d9SBarry Smith PETSC_INTERN PetscErrorCode MatGetDiagonalMarkers_SeqSBAIJ(Mat, const PetscInt **, PetscBool *);
345a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat, PetscInt, IS[], PetscInt);
35b2fa50c1SHong Zhang PETSC_INTERN PetscErrorCode MatSeqSBAIJZeroOps_Private(Mat);
367dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqSBAIJ(Mat, IS, IS, MatReuse, Mat *);
377dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_SeqSBAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
385a576424SJed Brown PETSC_INTERN PetscErrorCode MatScale_SeqSBAIJ(Mat, PetscScalar);
395a576424SJed Brown PETSC_INTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat, NormType, PetscReal *);
405a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat, Mat, PetscBool *);
415a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat, Vec);
425a576424SJed Brown PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat, Vec, Vec);
435a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat, MatInfoType, MatInfo *);
445a576424SJed Brown PETSC_INTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
455a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowMaxAbs_SeqSBAIJ(Mat, Vec, PetscInt[]);
465a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat, PetscInt *, PetscInt *, PetscInt *);
475a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqSBAIJ(Mat);
485a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqSBAIJ(Mat, PetscViewer);
4949b5e25fSSatish Balay 
505a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
515a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat, Mat, const MatFactorInfo *);
525a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat, Vec, Vec);
535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat, Vec, Vec);
54910cf402Sprj- PETSC_INTERN PetscErrorCode MatMatSolve_SeqSBAIJ_1_NaturalOrdering(Mat, Mat, Mat);
550a3c351aSHong Zhang 
565a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat, Vec, Vec);
575a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat, Vec, Vec);
585a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_inplace(Mat, Vec, Vec);
595a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_inplace(Mat, Vec, Vec);
600a3c351aSHong Zhang 
615a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat, Vec, Vec);
625a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat, Vec, Vec);
635a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1(Mat, Vec, Vec);
645a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1(Mat, Vec, Vec);
65266135f8SHong Zhang 
665a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
675a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat, Vec, Vec);
685a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat, Vec, Vec);
695a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat, Vec, Vec);
70266135f8SHong Zhang 
715a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
725a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat, Vec, Vec);
735a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat, Vec, Vec);
745a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat, Vec, Vec);
75266135f8SHong Zhang 
765a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
775a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat, Vec, Vec);
785a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat, Vec, Vec);
795a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat, Vec, Vec);
80266135f8SHong Zhang 
815a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
825a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat, Vec, Vec);
835a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat, Vec, Vec);
845a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat, Vec, Vec);
85266135f8SHong Zhang 
865a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
875a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat, Vec, Vec);
885a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat, Vec, Vec);
895a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat, Vec, Vec);
90266135f8SHong Zhang 
915a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
925a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat, Vec, Vec);
935a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat, Vec, Vec);
945a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat, Vec, Vec);
9549b5e25fSSatish Balay 
965a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat, Mat, const MatFactorInfo *);
975a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat, Vec, Vec);
985a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat, Vec, Vec);
995a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat, Vec, Vec);
1005a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_inplace(Mat, Vec, Vec);
1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_inplace(Mat, Vec, Vec);
1029495ac64SHong Zhang 
1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat, Mat, const MatFactorInfo *);
1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat, Mat, const MatFactorInfo *);
1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat, Mat, const MatFactorInfo *);
1065a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat, Mat, const MatFactorInfo *);
1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat, Mat, const MatFactorInfo *);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat, Mat, const MatFactorInfo *);
1095a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat, Mat, const MatFactorInfo *);
1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat, Mat, const MatFactorInfo *);
1119495ac64SHong Zhang 
1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_inplace(Mat, Vec, Vec);
1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_inplace(Mat, Vec, Vec);
1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat, Vec, Vec);
1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_inplace(Mat, Vec, Vec);
1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_inplace(Mat, Vec, Vec);
1175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_inplace(Mat, Vec, Vec);
1185a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_inplace(Mat, Vec, Vec);
1195a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_inplace(Mat, Vec, Vec);
1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_inplace(Mat, Vec, Vec);
1219495ac64SHong Zhang 
1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolves_SeqSBAIJ_1_inplace(Mat, Vecs, Vecs);
1235a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat, Vecs, Vecs);
1249495ac64SHong Zhang 
1255a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat, Vec, Vec);
1265a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat, Vec, Vec);
1275a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat, Vec, Vec);
1285a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat, Vec, Vec);
1295a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat, Vec, Vec);
1305a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat, Vec, Vec);
1315a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat, Vec, Vec);
1325a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat, Vec, Vec);
1339495ac64SHong Zhang 
1345a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat, Vec, Vec, Vec);
1355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat, Vec, Vec, Vec);
1365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat, Vec, Vec, Vec);
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat, Vec, Vec, Vec);
1385a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat, Vec, Vec, Vec);
1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat, Vec, Vec, Vec);
1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat, Vec, Vec, Vec);
1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat, Vec, Vec, Vec);
142671cb588SHong Zhang 
1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqSBAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqSBAIJ(Mat, PetscViewer);
1455a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat, PetscBool);
146d595f711SHong Zhang 
1474de5dceeSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat, Mat, PetscInt *);
1484de5dceeSHong Zhang 
14959f5e6ceSHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
15059f5e6ceSHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
151d79853d5SSatish Balay /* required by mpisbaij.c */
152d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatGetValues_SeqSBAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
153d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
154d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
155d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatGetRow_SeqSBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
156d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
157d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatZeroRows_SeqSBAIJ(Mat, IS, PetscScalar *, Vec, Vec);
158d79853d5SSatish Balay 
15917ea310bSPierre Jolivet PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqSBAIJ(Mat, PetscBool);
160