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 146c6c5352SBarry Smith 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 */ 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 */ 23b00f7748SHong Zhang 24b45a75daSHong Zhang /* carry MatFactorInfo from symbolic factor to numeric factor */ 2513f74950SBarry Smith PetscInt factor_levels; 26b45a75daSHong Zhang PetscReal factor_damping; 27f4cce38bSHong Zhang PetscReal factor_shift; 28b45a75daSHong Zhang PetscReal factor_zeropivot; 2949b5e25fSSatish Balay } Mat_SeqSBAIJ; 3049b5e25fSSatish Balay 31*f349cc81SSatish Balay EXTERN_C_BEGIN 32ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat,PetscInt,PetscInt,PetscInt*); 33*f349cc81SSatish Balay EXTERN_C_END 34dfbe8321SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat *); 35dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*); 36dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat); 3749b5e25fSSatish Balay 38af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 39dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 40266135f8SHong Zhang 41af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 42dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 43266135f8SHong Zhang 44af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 45dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 46266135f8SHong Zhang 47af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 48dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 49266135f8SHong Zhang 50af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 51dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 52266135f8SHong Zhang 53af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 54dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 55266135f8SHong Zhang 56af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 57dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 5849b5e25fSSatish Balay 59af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat,MatFactorInfo*,Mat*); 60dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 61671cb588SHong Zhang 6213f74950SBarry Smith EXTERN PetscErrorCode MatRelax_SeqSBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 63dfbe8321SBarry Smith EXTERN PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer,const MatType,Mat*); 64d06b337dSHong Zhang 6549b5e25fSSatish Balay #endif 66