1 2 #if !defined(__AIJ_H) 3 #define __AIJ_H 4 #include "src/mat/matimpl.h" 5 6 /* Info about i-nodes (identical nodes) */ 7 typedef struct { 8 PetscTruth use; 9 int node_count; /* number of inodes */ 10 int *size; /* size of each inode */ 11 int limit; /* inode limit */ 12 int max_limit; /* maximum supported inode limit */ 13 PetscTruth checked; /* if inodes have been checked for */ 14 } Mat_SeqAIJ_Inode; 15 16 /* 17 MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 18 format). The i[] and j[] arrays start at 0. For example, 19 j[i[k]+p] is the pth column in row k. Note that the diagonal 20 matrix elements are stored with the rest of the nonzeros (not separately). 21 */ 22 23 typedef struct { 24 PetscTruth sorted; /* if true, rows are sorted by increasing columns */ 25 PetscTruth roworiented; /* if true, row-oriented input, default */ 26 int nonew; /* 1 don't add new nonzeros, -1 generate error on new */ 27 PetscTruth singlemalloc; /* if true a, i, and j have been obtained with 28 one big malloc */ 29 PetscTruth freedata; /* free the i,j,a data when the matrix is destroyed; true by default */ 30 int nz,maxnz; /* nonzeros, allocated nonzeros */ 31 int *diag; /* pointers to diagonal elements */ 32 int *i; /* pointer to beginning of each row */ 33 int *imax; /* maximum space allocated for each row */ 34 int *ilen; /* actual length of each row */ 35 int *j; /* column values: j + i[k] - 1 is start of row k */ 36 PetscScalar *a; /* nonzero elements */ 37 IS row,col,icol; /* index sets, used for reorderings */ 38 PetscScalar *solve_work; /* work space used in MatSolve */ 39 Mat_SeqAIJ_Inode inode; /* identical node informaton */ 40 int reallocs; /* number of mallocs done during MatSetValues() 41 as more values are set than were prealloced */ 42 int rmax; /* max nonzeros in any row */ 43 PetscTruth ilu_preserve_row_sums; 44 PetscReal lu_dtcol; 45 PetscReal lu_damping; 46 PetscReal lu_shift; /* Manteuffel shift switch, fraction */ 47 PetscReal lu_shift_fraction; 48 PetscReal lu_zeropivot; 49 PetscScalar *saved_values; /* location for stashing nonzero values of matrix */ 50 PetscScalar *idiag,*ssor; /* inverse of diagonal entries; space for eisen */ 51 52 PetscTruth keepzeroedrows; /* keeps matrix structure same in calls to MatZeroRows()*/ 53 PetscTruth ignorezeroentries; 54 ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 55 Mat sbaijMat; /* mat in sbaij format */ 56 57 int *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */ 58 Mat XtoY; /* used by MatAXPY() */ 59 } Mat_SeqAIJ; 60 61 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat *); 62 EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat *); 63 EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat*); 64 EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat *); 65 EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 66 EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat); 67 EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 68 69 EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 70 EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 71 EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 72 EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 73 EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,int,int,Vec); 74 75 EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 76 EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 77 EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,int,void*); 78 79 EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,int *[],int *[]); 80 EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,int *[],int *[]); 81 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 82 EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth); 83 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*); 84 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 85 EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*); 86 EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 87 EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 88 EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 89 EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 90 EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 91 EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 92 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,MatFactorInfo*,IS,IS,Mat*); 93 EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer,const MatType,Mat*); 94 EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 95 EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 96 EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 97 EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 98 EXTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 99 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 100 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 101 EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 102 EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 103 EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 104 EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,int,const int[],int,const int[],const PetscScalar[],InsertMode); 105 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,int,int*,int**,PetscScalar**); 106 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,int,int*,int**,PetscScalar**); 107 EXTERN PetscErrorCode MatPrintHelp_SeqAIJ(Mat); 108 EXTERN PetscErrorCode MatAXPY_SeqAIJ(const PetscScalar[],Mat,Mat,MatStructure); 109 110 EXTERN_C_BEGIN 111 EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,Mat*); 112 EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,Mat*); 113 EXTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 114 EXTERN PetscErrorCode MatAdjustForInodes_SeqAIJ(Mat,IS*,IS*); 115 EXTERN PetscErrorCode MatSeqAIJGetInodeSizes_SeqAIJ(Mat,int*,int*[],int*); 116 EXTERN_C_END 117 118 #endif 119