1*49b5e25fSSatish Balay /* $Id: sbaij.h,v 1.21 2000/01/11 21:00:52 bsmith Exp $ */ 2*49b5e25fSSatish Balay 3*49b5e25fSSatish Balay #include "src/mat/matimpl.h" 4*49b5e25fSSatish Balay 5*49b5e25fSSatish Balay #if !defined(__SBAIJ_H) 6*49b5e25fSSatish Balay #define __SBAIJ_H 7*49b5e25fSSatish Balay 8*49b5e25fSSatish Balay /* 9*49b5e25fSSatish Balay MATSEQSBAIJ format - Block compressed row storage. The i[] and j[] 10*49b5e25fSSatish Balay arrays start at 0. 11*49b5e25fSSatish Balay */ 12*49b5e25fSSatish Balay 13*49b5e25fSSatish Balay typedef struct { 14*49b5e25fSSatish Balay PetscTruth sorted; /* if true, rows are sorted by increasing columns */ 15*49b5e25fSSatish Balay PetscTruth roworiented; /* if true, row-oriented input, default */ 16*49b5e25fSSatish Balay int nonew; /* 1 don't add new nonzeros, -1 generate error on new */ 17*49b5e25fSSatish Balay PetscTruth singlemalloc; /* if true a, i, and j have been obtained with 18*49b5e25fSSatish Balay one big malloc */ 19*49b5e25fSSatish Balay int m,n; /* rows or columns */ 20*49b5e25fSSatish Balay int bs,bs2; /* block size, square of block size */ 21*49b5e25fSSatish Balay int mbs,nbs; /* rows/bs or columns/bs */ 22*49b5e25fSSatish Balay int s_nz,s_maxnz; /* total diagonal and superdiagonal nonzero blocks, 23*49b5e25fSSatish Balay total allocated diagonal and superdiagonal nonzero blocks*/ 24*49b5e25fSSatish Balay int *diag; /* pointers to diagonal elements */ 25*49b5e25fSSatish Balay int *i; /* pointer to beginning of each row */ 26*49b5e25fSSatish Balay int *imax; /* maximum space allocated for each row */ 27*49b5e25fSSatish Balay int *ilen; /* actual length of each row */ 28*49b5e25fSSatish Balay int *j; /* column values: j + i[k] is start of row k */ 29*49b5e25fSSatish Balay MatScalar *a; /* nonzero diagonal and subdiagonal elements */ 30*49b5e25fSSatish Balay IS row,col,icol; /* index sets, used for reorderings */ 31*49b5e25fSSatish Balay Scalar *solve_work; /* work space used in MatSolve */ 32*49b5e25fSSatish Balay void *spptr; /* pointer for special library like SuperLU */ 33*49b5e25fSSatish Balay int reallocs; /* number of mallocs done during MatSetValues() 34*49b5e25fSSatish Balay as more values are set then were preallocated */ 35*49b5e25fSSatish Balay Scalar *mult_work; /* work array for matrix vector product*/ 36*49b5e25fSSatish Balay Scalar *saved_values; 37*49b5e25fSSatish Balay 38*49b5e25fSSatish Balay PetscTruth keepzeroedrows; /* if true, MatZeroRows() will not change nonzero structure */ 39*49b5e25fSSatish Balay } Mat_SeqSBAIJ; 40*49b5e25fSSatish Balay 41*49b5e25fSSatish Balay extern int MatILUFactorSymbolic_SeqSBAIJ(Mat,IS,IS,MatILUInfo*,Mat *); 42*49b5e25fSSatish Balay extern int MatConvert_SeqSBAIJ(Mat,MatType,Mat *); 43*49b5e25fSSatish Balay extern int MatDuplicate_SeqSBAIJ(Mat,MatDuplicateOption,Mat*); 44*49b5e25fSSatish Balay extern int MatMarkDiagonal_SeqSBAIJ(Mat); 45*49b5e25fSSatish Balay 46*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 47*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat,Mat*); 48*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 49*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 50*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat,Mat*); 51*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 52*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 53*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat,Mat*); 54*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 55*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 56*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_5_NaturalOrdering(Mat,Mat*); 57*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 58*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 59*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat,Mat*); 60*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 61*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 62*49b5e25fSSatish Balay extern int MatLUFactorNumeric_SeqSBAIJ_7_NaturalOrdering(Mat,Mat*); 63*49b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 64*49b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 65*49b5e25fSSatish Balay 66*49b5e25fSSatish Balay extern int MatReordering_SeqSBAIJ(Mat,IS); 67*49b5e25fSSatish Balay /* 68*49b5e25fSSatish Balay extern int MatCreateSeqSBAIJ(MPI_Comm,int,int,int,int,int*,Mat*); 69*49b5e25fSSatish Balay extern int MatSetValuesBlocked_SeqSBAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode); 70*49b5e25fSSatish Balay */ 71*49b5e25fSSatish Balay #endif 72