#if !defined(__BAIJ_H) #define __BAIJ_H #include "include/private/matimpl.h" #include "src/mat/impls/aij/seq/aij.h" /* MATSEQBAIJ format - Block compressed row storage. The i[] and j[] arrays start at 0. */ /* This header is shared by the SeqSBAIJ matrix */ #define SEQBAIJHEADER \ PetscInt bs2; /* square of block size */ \ PetscInt mbs,nbs; /* rows/bs, columns/bs */ \ PetscScalar *mult_work; /* work array for matrix vector product*/ \ PetscScalar *saved_values; \ \ Mat sbaijMat; /* mat in sbaij format */ \ \ PetscInt setvalueslen; /* only used for single precision */ \ MatScalar *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied \ before calling MatSetValuesXXX_SeqBAIJ_MatScalar() */ \ \ PetscTruth pivotinblocks; /* pivot inside factorization of each diagonal block */ \ \ PetscScalar *idiag; /* inverse of block diagonal */ \ PetscTruth idiagvalid /* if above has correct/current values */ typedef struct { SEQAIJHEADER(PetscScalar); SEQBAIJHEADER; } Mat_SeqBAIJ; EXTERN_C_BEGIN EXTERN PetscErrorCode MatSeqBAIJSetPreallocation_SeqBAIJ(Mat,PetscInt,PetscInt,PetscInt*); EXTERN_C_END EXTERN PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *); EXTERN PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat,IS,MatFactorInfo*,Mat *); EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat,IS,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*); EXTERN PetscErrorCode MatMissingDiagonal_SeqBAIJ(Mat); EXTERN PetscErrorCode MatMarkDiagonal_SeqBAIJ(Mat); EXTERN PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*); EXTERN PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat,PetscInt,IS*,PetscInt); EXTERN PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); EXTERN PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]); EXTERN PetscErrorCode MatMultTranspose_SeqBAIJ(Mat,Vec,Vec); EXTERN PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatScale_SeqBAIJ(Mat,PetscScalar); EXTERN PetscErrorCode MatNorm_SeqBAIJ(Mat,NormType,PetscReal *); EXTERN PetscErrorCode MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*); EXTERN PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat,Vec); EXTERN PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); EXTERN PetscErrorCode MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); EXTERN PetscErrorCode MatZeroEntries_SeqBAIJ(Mat); EXTERN PetscErrorCode MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat); EXTERN PetscErrorCode MatSeqBAIJ_UpdateSolvers(Mat); EXTERN PetscErrorCode MatSolve_SeqBAIJ_Update(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_1(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_2(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_3(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_4(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec); #if defined(PETSC_HAVE_SSE) EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec); #endif EXTERN PetscErrorCode MatSolve_SeqBAIJ_5(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_6(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_7(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolve_SeqBAIJ_N(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec); EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,MatFactorInfo*,Mat*); #if defined(PETSC_HAVE_SSE) EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,MatFactorInfo*,Mat*); #else #endif EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat,MatFactorInfo*,Mat*); EXTERN PetscErrorCode MatMult_SeqBAIJ_1(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_2(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_3(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_4(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_5(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_6(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_7(Mat,Vec,Vec); EXTERN PetscErrorCode MatMult_SeqBAIJ_N(Mat,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); EXTERN PetscErrorCode MatLoad_SeqBAIJ(PetscViewer, MatType,Mat*); #endif