xref: /petsc/src/mat/impls/sbaij/seq/sbaij.h (revision e2ee2a47cae8248a2ea96b9d6f35ff3f37ecfc00)
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 {
136c6c5352SBarry Smith   SEQBAIJHEADER
1413f74950SBarry Smith   PetscInt         *inew;        /* pointer to beginning of each row of reordered matrix */
1513f74950SBarry Smith   PetscInt         *jnew;        /* column values: jnew + i[k] is start of row k */
16a6175056SHong Zhang   MatScalar        *anew;        /* nonzero diagonal and superdiagonal elements of reordered matrix */
17d59c15a7SBarry Smith   PetscScalar      *solves_work; /* work space used in MatSolves */
1813f74950SBarry Smith   PetscInt         solves_work_n;/* size of solves_work */
19*e2ee2a47SBarry Smith   PetscScalar      *relax_work;
2013f74950SBarry Smith   PetscInt         *a2anew;        /* map used for symm permutation */
21d36791b2SHong Zhang   PetscTruth       permute;        /* if true, a non-trivial permutation is used for factorization */
22941593c8SHong Zhang   PetscTruth       ignore_ltriangular; /* if true, ignore the lower triangular values inserted by users */
23f5edf698SHong Zhang   PetscTruth       getrow_utriangular; /* if true, MatGetRow_SeqSBAIJ() is enabled to get the upper part of the row */
2449b5e25fSSatish Balay } Mat_SeqSBAIJ;
2549b5e25fSSatish Balay 
26f349cc81SSatish Balay EXTERN_C_BEGIN
27ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
28f349cc81SSatish Balay EXTERN_C_END
29dfbe8321SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat *);
30dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*);
31dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat);
3249b5e25fSSatish Balay 
33af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
34dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
35266135f8SHong Zhang 
36af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
37dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
38266135f8SHong Zhang 
39af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
40dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
41266135f8SHong Zhang 
42af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
43dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
44266135f8SHong Zhang 
45af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
46dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
47266135f8SHong Zhang 
48af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
49dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
50266135f8SHong Zhang 
51af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
52dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
5349b5e25fSSatish Balay 
54af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
55dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
56671cb588SHong Zhang 
5713f74950SBarry Smith EXTERN PetscErrorCode MatRelax_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
58f69a0ea3SMatthew Knepley EXTERN PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer, MatType,Mat*);
59d06b337dSHong Zhang 
6049b5e25fSSatish Balay #endif
61