xref: /petsc/src/mat/impls/sbaij/seq/sbaij.h (revision 0a3c351a51b8afeba0e28ebe856db6cb3362db0e)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay #if !defined(__SBAIJ_H)
349b5e25fSSatish Balay #define __SBAIJ_H
47c4f633dSBarry Smith #include "private/matimpl.h"
57c4f633dSBarry Smith #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 */
2041f059aeSBarry Smith   PetscScalar      *sor_work;
2113f74950SBarry Smith   PetscInt         *a2anew;        /* map used for symm permutation */
22d36791b2SHong Zhang   PetscTruth       permute;        /* if true, a non-trivial permutation is used for factorization */
23941593c8SHong Zhang   PetscTruth       ignore_ltriangular; /* if true, ignore the lower triangular values inserted by users */
24f5edf698SHong Zhang   PetscTruth       getrow_utriangular; /* if true, MatGetRow_SeqSBAIJ() is enabled to get the upper part of the row */
254108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
2638702af4SBarry Smith   unsigned short   *jshort;
274da8f245SBarry Smith   PetscTruth       free_jshort;
2849b5e25fSSatish Balay } Mat_SeqSBAIJ;
2949b5e25fSSatish Balay 
30f349cc81SSatish Balay EXTERN_C_BEGIN
31ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
32f349cc81SSatish Balay EXTERN_C_END
330481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
340481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,const MatFactorInfo*);
350481f469SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
36dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*);
37dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat);
389495ac64SHong Zhang EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt);
394aa3045dSJed Brown EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,MatReuse,Mat*);
409495ac64SHong Zhang EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
419495ac64SHong Zhang EXTERN PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar);
429495ac64SHong Zhang EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
439495ac64SHong Zhang EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
449495ac64SHong Zhang EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
459495ac64SHong Zhang EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
469495ac64SHong Zhang EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
479495ac64SHong Zhang EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
48985db425SBarry Smith EXTERN PetscErrorCode MatGetRowMaxAbs_SeqSBAIJ(Mat,Vec,PetscInt[]);
499495ac64SHong Zhang EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*);
5049b5e25fSSatish Balay 
510481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
52*0a3c351aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
535f44c854SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_newdatastruct(Mat,Vec,Vec);
54*0a3c351aSHong Zhang 
55*0a3c351aSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
56*0a3c351aSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
57*0a3c351aSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
58*0a3c351aSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
59*0a3c351aSHong Zhang 
60*0a3c351aSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_newdatastruct(Mat,Vec,Vec);
61*0a3c351aSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_newdatastruct(Mat,Vec,Vec);
62*0a3c351aSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_1_newdatastruct(Mat,Vec,Vec);
63*0a3c351aSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_newdatastruct(Mat,Vec,Vec);
64266135f8SHong Zhang 
650481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
66dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
67759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
68759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
69266135f8SHong Zhang 
700481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
71dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
72759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
73759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
74266135f8SHong Zhang 
750481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
76dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
77759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
78759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
79266135f8SHong Zhang 
800481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
81dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
82759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
83759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
84266135f8SHong Zhang 
850481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
86dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
87759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
88759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
89266135f8SHong Zhang 
900481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
91dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
92759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
93759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
9449b5e25fSSatish Balay 
950481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
96dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
97759322afSHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
98759322afSHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
999495ac64SHong Zhang EXTERN PetscErrorCode MatForwardSolve_SeqSBAIJ_N(Mat,Vec,Vec);
1009495ac64SHong Zhang EXTERN PetscErrorCode MatBackwardSolve_SeqSBAIJ_N(Mat,Vec,Vec);
1019495ac64SHong Zhang 
1020481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat,const MatFactorInfo*);
1030481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat,const MatFactorInfo*);
1040481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat,const MatFactorInfo*);
1050481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat,const MatFactorInfo*);
1060481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat,const MatFactorInfo*);
1070481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat,const MatFactorInfo*);
1080481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat,const MatFactorInfo*);
1090481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat,const MatFactorInfo*);
1109495ac64SHong Zhang 
1119495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
112*0a3c351aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
113783ef271SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_newdatastruct(Mat,Vec,Vec);
1149495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
1159495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
1169495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
1179495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
1189495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
1199495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
1209495ac64SHong Zhang 
1219495ac64SHong Zhang EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1229495ac64SHong Zhang 
123*0a3c351aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
124*0a3c351aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_newdatastruct(Mat,Vec,Vec);
1259495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
1269495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
1279495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
1289495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
1299495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
1309495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
1319495ac64SHong Zhang EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
1329495ac64SHong Zhang 
1339495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
1349495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
1359495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
1369495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
1379495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
1389495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
1399495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
1409495ac64SHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
1419495ac64SHong Zhang 
142eeffb40dSHong Zhang EXTERN PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat,Vec,Vec);
143eeffb40dSHong Zhang 
1449495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
1459495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
1469495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
1479495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
1489495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
1499495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
1509495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
1519495ac64SHong Zhang EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
152671cb588SHong Zhang 
15341f059aeSBarry Smith EXTERN PetscErrorCode MatSOR_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
154a313700dSBarry Smith EXTERN PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer, const MatType,Mat*);
155d06b337dSHong Zhang 
15649b5e25fSSatish Balay #endif
157