xref: /petsc/src/mat/impls/sbaij/seq/sbaij.h (revision b2fa50c1850a415a4262ea91c382596971d9ac8b)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay #if !defined(__SBAIJ_H)
349b5e25fSSatish Balay #define __SBAIJ_H
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
649b5e25fSSatish Balay 
749b5e25fSSatish Balay /*
849b5e25fSSatish Balay   MATSEQSBAIJ format - Block compressed row storage. The i[] and j[]
949b5e25fSSatish Balay   arrays start at 0.
1049b5e25fSSatish Balay */
1149b5e25fSSatish Balay 
1249b5e25fSSatish Balay typedef struct {
13dd6ea824SBarry Smith   SEQAIJHEADER(MatScalar);
14e6b907acSBarry Smith   SEQBAIJHEADER;
1513f74950SBarry Smith   PetscInt         *inew;        /* pointer to beginning of each row of reordered matrix */
1613f74950SBarry Smith   PetscInt         *jnew;        /* column values: jnew + i[k] is start of row k */
17a6175056SHong Zhang   MatScalar        *anew;        /* nonzero diagonal and superdiagonal elements of reordered matrix */
18d59c15a7SBarry Smith   PetscScalar      *solves_work; /* work space used in MatSolves */
1913f74950SBarry Smith   PetscInt         solves_work_n; /* size of solves_work */
2013f74950SBarry Smith   PetscInt         *a2anew;        /* map used for symm permutation */
21ace3abfcSBarry Smith   PetscBool        permute;        /* if true, a non-trivial permutation is used for factorization */
22ace3abfcSBarry Smith   PetscBool        ignore_ltriangular; /* if true, ignore the lower triangular values inserted by users */
23ace3abfcSBarry Smith   PetscBool        getrow_utriangular; /* if true, MatGetRow_SeqSBAIJ() is enabled to get the upper part of the row */
244108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
2538702af4SBarry Smith   unsigned short   *jshort;
26ace3abfcSBarry Smith   PetscBool        free_jshort;
2749b5e25fSSatish Balay } Mat_SeqSBAIJ;
2849b5e25fSSatish Balay 
298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
305a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
315a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
325a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,const MatFactorInfo*);
335a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
345a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
355a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*);
365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat);
375a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt);
38*b2fa50c1SHong Zhang PETSC_INTERN PetscErrorCode MatSeqSBAIJZeroOps_Private(Mat);
395a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,MatReuse,Mat*);
405a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
415a576424SJed Brown PETSC_INTERN PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar);
425a576424SJed Brown PETSC_INTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal*);
435a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscBool*);
445a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
455a576424SJed Brown PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
465a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo*);
475a576424SJed Brown PETSC_INTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
485a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowMaxAbs_SeqSBAIJ(Mat,Vec,PetscInt[]);
495a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*);
505a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqSBAIJ(Mat);
515a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqSBAIJ(Mat,PetscViewer);
5249b5e25fSSatish Balay 
535a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
545a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Mat,const MatFactorInfo*);
555a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
570a3c351aSHong Zhang 
585a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
595a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
605a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
615a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
620a3c351aSHong Zhang 
635a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
645a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
655a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1(Mat,Vec,Vec);
665a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1(Mat,Vec,Vec);
67266135f8SHong Zhang 
685a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
695a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
705a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
715a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
72266135f8SHong Zhang 
735a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
745a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
755a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
765a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
77266135f8SHong Zhang 
785a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
795a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
805a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
815a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
82266135f8SHong Zhang 
835a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
845a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
855a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
865a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
87266135f8SHong Zhang 
885a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
895a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
905a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
915a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
92266135f8SHong Zhang 
935a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
945a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
955a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
965a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
9749b5e25fSSatish Balay 
985a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
995a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
1005a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
1049495ac64SHong Zhang 
1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat,const MatFactorInfo*);
1065a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat,Mat,const MatFactorInfo*);
1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat,const MatFactorInfo*);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat,const MatFactorInfo*);
1095a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat,const MatFactorInfo*);
1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat,const MatFactorInfo*);
1115a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat,const MatFactorInfo*);
1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat,const MatFactorInfo*);
1139495ac64SHong Zhang 
1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
1175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_inplace(Mat,Vec,Vec);
1185a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_inplace(Mat,Vec,Vec);
1195a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_inplace(Mat,Vec,Vec);
1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_inplace(Mat,Vec,Vec);
1215a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_inplace(Mat,Vec,Vec);
1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_inplace(Mat,Vec,Vec);
1239495ac64SHong Zhang 
1245a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolves_SeqSBAIJ_1_inplace(Mat,Vecs,Vecs);
1255a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1269495ac64SHong Zhang 
1275a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
1285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
1295a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
1305a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
1315a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
1325a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
1335a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
1345a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
1355a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
1369495ac64SHong Zhang 
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
1385a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
1425a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
1459495ac64SHong Zhang 
1465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat,Vec,Vec);
147eeffb40dSHong Zhang 
1485a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
1495a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
1505a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
1515a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
1525a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
1535a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
1545a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
1555a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
156671cb588SHong Zhang 
1575a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1585a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqSBAIJ(Mat,PetscViewer);
1595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat,PetscBool);
160d595f711SHong Zhang 
1614de5dceeSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat,Mat,PetscInt*);
1624de5dceeSHong Zhang 
16359f5e6ceSHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
16459f5e6ceSHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
165d79853d5SSatish Balay /* required by mpisbaij.c */
166d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatGetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
167d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
168d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
169d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatGetRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
170d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
171d79853d5SSatish Balay PETSC_INTERN PetscErrorCode MatZeroRows_SeqSBAIJ(Mat,IS,PetscScalar*,Vec,Vec);
172d79853d5SSatish Balay 
17349b5e25fSSatish Balay #endif
174