xref: /petsc/src/mat/impls/sbaij/seq/sbaij.h (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
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 */
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 
29f349cc81SSatish Balay EXTERN_C_BEGIN
30*09573ac7SBarry Smith extern PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
31f349cc81SSatish Balay EXTERN_C_END
32*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
33*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
34*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,const MatFactorInfo*);
35*09573ac7SBarry Smith extern PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,Mat,IS,const MatFactorInfo*);
36*09573ac7SBarry Smith extern PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
37*09573ac7SBarry Smith extern PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*);
38*09573ac7SBarry Smith extern PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat);
39*09573ac7SBarry Smith extern PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt);
40*09573ac7SBarry Smith extern PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,MatReuse,Mat*);
41*09573ac7SBarry Smith extern PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
42*09573ac7SBarry Smith extern PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar);
43*09573ac7SBarry Smith extern PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
44*09573ac7SBarry Smith extern PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscBool *);
45*09573ac7SBarry Smith extern PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
46*09573ac7SBarry Smith extern PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
47*09573ac7SBarry Smith extern PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
48*09573ac7SBarry Smith extern PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
49*09573ac7SBarry Smith extern PetscErrorCode MatGetRowMaxAbs_SeqSBAIJ(Mat,Vec,PetscInt[]);
50*09573ac7SBarry Smith extern PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*);
51*09573ac7SBarry Smith extern PetscErrorCode MatDestroy_SeqSBAIJ(Mat);
52*09573ac7SBarry Smith extern PetscErrorCode MatView_SeqSBAIJ(Mat,PetscViewer);
5349b5e25fSSatish Balay 
54*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
55*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Mat,const MatFactorInfo*);
56*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
57*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
580a3c351aSHong Zhang 
59*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
60*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
61*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
62*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
630a3c351aSHong Zhang 
64*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
65*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
66*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_1(Mat,Vec,Vec);
67*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_1(Mat,Vec,Vec);
68266135f8SHong Zhang 
69*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
70*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
71*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
72*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
73266135f8SHong Zhang 
74*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
75*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
76*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
77*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
78266135f8SHong Zhang 
79*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
80*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
81*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
82*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
83266135f8SHong Zhang 
84*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
85*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
86*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
87*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
88266135f8SHong Zhang 
89*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
90*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
91*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
92*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
93266135f8SHong Zhang 
94*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
95*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
96*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
97*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
9849b5e25fSSatish Balay 
99*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,Mat,const MatFactorInfo*);
100*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
101*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
102*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
103*09573ac7SBarry Smith extern PetscErrorCode MatForwardSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
104*09573ac7SBarry Smith extern PetscErrorCode MatBackwardSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
1059495ac64SHong Zhang 
106*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat,const MatFactorInfo*);
107*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat,Mat,const MatFactorInfo*);
108*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat,const MatFactorInfo*);
109*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat,const MatFactorInfo*);
110*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat,const MatFactorInfo*);
111*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat,const MatFactorInfo*);
112*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat,const MatFactorInfo*);
113*09573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat,const MatFactorInfo*);
1149495ac64SHong Zhang 
115*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_N_inplace(Mat,Vec,Vec);
116*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_1_inplace(Mat,Vec,Vec);
117*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
118*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_2_inplace(Mat,Vec,Vec);
119*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_3_inplace(Mat,Vec,Vec);
120*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_4_inplace(Mat,Vec,Vec);
121*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_5_inplace(Mat,Vec,Vec);
122*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_6_inplace(Mat,Vec,Vec);
123*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_7_inplace(Mat,Vec,Vec);
1249495ac64SHong Zhang 
125*09573ac7SBarry Smith extern PetscErrorCode MatSolves_SeqSBAIJ_1_inplace(Mat,Vecs,Vecs);
126*09573ac7SBarry Smith extern PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1279495ac64SHong Zhang 
128*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace(Mat,Vec,Vec);
129*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
130*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace(Mat,Vec,Vec);
131*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace(Mat,Vec,Vec);
132*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace(Mat,Vec,Vec);
133*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering_inplace(Mat,Vec,Vec);
134*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace(Mat,Vec,Vec);
135*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering_inplace(Mat,Vec,Vec);
136*09573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace(Mat,Vec,Vec);
1379495ac64SHong Zhang 
138*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
139*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
140*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
141*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
142*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
143*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
144*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
145*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
1469495ac64SHong Zhang 
147*09573ac7SBarry Smith extern PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat,Vec,Vec);
148eeffb40dSHong Zhang 
149*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
150*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
151*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
152*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
153*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
154*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
155*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
156*09573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
157671cb588SHong Zhang 
158*09573ac7SBarry Smith extern PetscErrorCode MatSOR_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
159*09573ac7SBarry Smith extern PetscErrorCode MatLoad_SeqSBAIJ(Mat,PetscViewer);
160*09573ac7SBarry Smith extern PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat,PetscBool );
161d595f711SHong Zhang 
16249b5e25fSSatish Balay #endif
163