xref: /petsc/src/mat/impls/sbaij/seq/sbaij.h (revision 421e10b8f79bbed49dcc4e803c884835d979c6ea)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay #if !defined(__SBAIJ_H)
349b5e25fSSatish Balay #define __SBAIJ_H
4c833c1fbSBarry Smith #include "src/mat/matimpl.h"
56c6c5352SBarry 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 {
13*421e10b8SBarry Smith   SEQAIJHEADER(PetscScalar);
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 */
20e2ee2a47SBarry Smith   PetscScalar      *relax_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 */
2549b5e25fSSatish Balay } Mat_SeqSBAIJ;
2649b5e25fSSatish Balay 
27f349cc81SSatish Balay EXTERN_C_BEGIN
28ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
29f349cc81SSatish Balay EXTERN_C_END
30dfbe8321SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat *);
31dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*);
32dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat);
3349b5e25fSSatish Balay 
34af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
35dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
36266135f8SHong Zhang 
37af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
38dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
39266135f8SHong Zhang 
40af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
41dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
42266135f8SHong Zhang 
43af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
44dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
45266135f8SHong Zhang 
46af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
47dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
48266135f8SHong Zhang 
49af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
50dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
51266135f8SHong Zhang 
52af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
53dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
5449b5e25fSSatish Balay 
55af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
56dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
57671cb588SHong Zhang 
5813f74950SBarry Smith EXTERN PetscErrorCode MatRelax_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
59f69a0ea3SMatthew Knepley EXTERN PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer, MatType,Mat*);
60d06b337dSHong Zhang 
6149b5e25fSSatish Balay #endif
62