1 2 #if !defined(__AIJ_H) 3 #define __AIJ_H 4 5 #include "private/matimpl.h" 6 /* 7 Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 8 */ 9 #define SEQAIJHEADER(datatype) \ 10 PetscTruth roworiented; /* if true, row-oriented input, default */\ 11 PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */\ 12 PetscInt nounused; /* -1 generate error on unused space */\ 13 PetscTruth singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\ 14 PetscInt maxnz; /* allocated nonzeros */\ 15 PetscInt *imax; /* maximum space allocated for each row */\ 16 PetscInt *ilen; /* actual length of each row */\ 17 PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 18 as more values are set than were prealloced */\ 19 PetscInt rmax; /* max nonzeros in any row */\ 20 PetscTruth keepzeroedrows; /* keeps matrix structure same in calls to MatZeroRows()*/\ 21 PetscTruth ignorezeroentries; \ 22 PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 23 Mat XtoY; /* used by MatAXPY() */\ 24 PetscTruth free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 25 PetscTruth free_a; /* free the numerical values when matrix is destroy */ \ 26 Mat_CompressedRow compressedrow; /* use compressed row format */ \ 27 PetscInt nz; /* nonzeros */ \ 28 PetscInt *i; /* pointer to beginning of each row */ \ 29 PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 30 PetscInt *diag; /* pointers to diagonal elements */ \ 31 datatype *a; /* nonzero elements */ \ 32 PetscScalar *solve_work; /* work space used in MatSolve */ \ 33 IS row, col, icol /* index sets, used for reorderings */ 34 35 /* 36 MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 37 format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 38 j[i[k]+p] is the pth column in row k. Note that the diagonal 39 matrix elements are stored with the rest of the nonzeros (not separately). 40 */ 41 42 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 43 typedef struct { 44 MatScalar *bdiag,*ibdiag; /* diagonal blocks of matrix used for MatRelax_Inode() */ 45 PetscInt bdiagsize; /* length of bdiag and ibdiag */ 46 PetscTruth ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 47 48 PetscTruth use; 49 PetscInt node_count; /* number of inodes */ 50 PetscInt *size; /* size of each inode */ 51 PetscInt limit; /* inode limit */ 52 PetscInt max_limit; /* maximum supported inode limit */ 53 PetscTruth checked; /* if inodes have been checked for */ 54 } Mat_Inode; 55 56 EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer); 57 EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType); 58 EXTERN PetscErrorCode MatDestroy_Inode(Mat); 59 EXTERN PetscErrorCode MatCreate_Inode(Mat); 60 EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth); 61 EXTERN PetscErrorCode MatDuplicate_Inode(Mat,MatDuplicateOption,Mat*); 62 EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat,IS,IS,const MatFactorInfo*,Mat*); 63 EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 64 EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 65 66 67 typedef struct { 68 SEQAIJHEADER(MatScalar); 69 Mat_Inode inode; 70 MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 71 72 PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 73 PetscTruth idiagvalid; /* current idiag[] and mdiag[] are valid */ 74 PetscScalar fshift,omega; /* last used omega and fshift */ 75 76 ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 77 } Mat_SeqAIJ; 78 79 /* 80 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 81 */ 82 #undef __FUNCT__ 83 #define __FUNCT__ "MatSeqXAIJFreeAIJ" 84 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 85 { 86 PetscErrorCode ierr; 87 Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 88 if (A->singlemalloc) { 89 ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 90 } else { 91 if (A->free_a && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 92 if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);} 93 if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);} 94 } 95 *a = 0; *j = 0; *i = 0; 96 return 0; 97 } 98 99 /* 100 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 101 */ 102 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 103 if (NROW >= RMAX) {\ 104 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 105 /* there is no extra room in row, therefore enlarge */ \ 106 PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 107 datatype *new_a; \ 108 \ 109 if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 110 /* malloc new storage space */ \ 111 ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 112 \ 113 /* copy over old data into new slots */ \ 114 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 115 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 116 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 117 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 118 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 119 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 120 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 121 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 122 /* free up old matrix storage */ \ 123 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 124 AA = new_a; \ 125 Ain->a = (MatScalar*) new_a; \ 126 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 127 Ain->singlemalloc = PETSC_TRUE; \ 128 \ 129 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 130 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 131 Ain->maxnz += CHUNKSIZE; \ 132 Ain->reallocs++; \ 133 } \ 134 135 EXTERN_C_BEGIN 136 EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*); 137 EXTERN_C_END 138 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 139 EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 140 EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 141 EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 142 EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 143 EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*); 144 EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 145 146 EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 147 EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 148 EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 149 EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 150 EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 151 152 EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 153 EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 154 EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 155 156 EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 157 EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 158 EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 159 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 160 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 161 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 162 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 163 EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 164 EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 165 EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 166 EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 167 EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 168 EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 169 EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 170 EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 171 EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 172 EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 173 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*,Mat*); 174 EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*); 175 EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 176 EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 177 EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 178 EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 179 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 180 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 181 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 182 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 183 EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 184 EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 185 EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 186 EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 187 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 188 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 189 EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 190 EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 191 EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 192 EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 193 EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 194 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 195 EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 196 197 EXTERN_C_BEGIN 198 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 199 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 200 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*); 201 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 202 EXTERN_C_END 203 204 #endif 205