1 2 #if !defined(__AIJ_H) 3 #define __AIJ_H 4 5 #include <private/matimpl.h> 6 7 /* 8 Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 9 */ 10 #define SEQAIJHEADER(datatype) \ 11 PetscBool roworiented; /* if true, row-oriented input, default */\ 12 PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */\ 13 PetscInt nounused; /* -1 generate error on unused space */\ 14 PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\ 15 PetscInt maxnz; /* allocated nonzeros */\ 16 PetscInt *imax; /* maximum space allocated for each row */\ 17 PetscInt *ilen; /* actual length of each row */\ 18 PetscBool free_imax_ilen; \ 19 PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 20 as more values are set than were prealloced */\ 21 PetscInt rmax; /* max nonzeros in any row */\ 22 PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/\ 23 PetscBool ignorezeroentries; \ 24 PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 25 Mat XtoY; /* used by MatAXPY() */\ 26 PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 27 PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 28 Mat_CompressedRow compressedrow; /* use compressed row format */ \ 29 PetscInt nz; /* nonzeros */ \ 30 PetscInt *i; /* pointer to beginning of each row */ \ 31 PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 32 PetscInt *diag; /* pointers to diagonal elements */ \ 33 PetscBool free_diag; \ 34 datatype *a; /* nonzero elements */ \ 35 PetscScalar *solve_work; /* work space used in MatSolve */ \ 36 IS row, col, icol; /* index sets, used for reorderings */ \ 37 PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 38 Mat parent /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); 39 means that this shares some data structures with the parent including diag, ilen, imax, i, j */ 40 41 /* 42 MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 43 format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 44 j[i[k]+p] is the pth column in row k. Note that the diagonal 45 matrix elements are stored with the rest of the nonzeros (not separately). 46 */ 47 48 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 49 typedef struct { 50 MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 51 PetscInt bdiagsize; /* length of bdiag and ibdiag */ 52 PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 53 54 PetscBool use; 55 PetscInt node_count; /* number of inodes */ 56 PetscInt *size; /* size of each inode */ 57 PetscInt limit; /* inode limit */ 58 PetscInt max_limit; /* maximum supported inode limit */ 59 PetscBool checked; /* if inodes have been checked for */ 60 } Mat_SeqAIJ_Inode; 61 62 extern PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 63 extern PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 64 extern PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 65 extern PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 66 extern PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool ); 67 extern PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 68 extern PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool ); 69 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 70 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 71 72 typedef struct { 73 SEQAIJHEADER(MatScalar); 74 Mat_SeqAIJ_Inode inode; 75 MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 76 77 PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 78 PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 79 PetscScalar *ibdiag; /* inverses of block diagonals */ 80 PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 81 PetscScalar fshift,omega; /* last used omega and fshift */ 82 83 ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 84 } Mat_SeqAIJ; 85 86 /* 87 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 88 */ 89 #undef __FUNCT__ 90 #define __FUNCT__ "MatSeqXAIJFreeAIJ" 91 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 92 { 93 PetscErrorCode ierr; 94 Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 95 if (A->singlemalloc) { 96 ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 97 } else { 98 if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 99 if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 100 if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 101 } 102 return 0; 103 } 104 /* 105 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 106 This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 107 */ 108 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 109 if (NROW >= RMAX) {\ 110 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 111 /* there is no extra room in row, therefore enlarge */ \ 112 PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 113 datatype *new_a; \ 114 \ 115 if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 116 /* malloc new storage space */ \ 117 ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 118 \ 119 /* copy over old data into new slots */ \ 120 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 121 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 122 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 123 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 124 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 125 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 126 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 127 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 128 /* free up old matrix storage */ \ 129 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 130 AA = new_a; \ 131 Ain->a = (MatScalar*) new_a; \ 132 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 133 Ain->singlemalloc = PETSC_TRUE; \ 134 \ 135 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 136 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 137 Ain->maxnz += BS2*CHUNKSIZE; \ 138 Ain->reallocs++; \ 139 } \ 140 141 142 EXTERN_C_BEGIN 143 extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 144 EXTERN_C_END 145 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 146 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 147 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 148 149 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 150 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 151 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 152 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 153 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 154 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 155 extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 156 extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 157 extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*); 158 extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 159 160 extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 161 extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 162 extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 163 extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 164 extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 165 166 extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 167 extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 168 extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 169 170 extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 171 extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 172 extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 173 extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 174 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 175 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 176 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 177 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 178 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 179 extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 180 extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 181 extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 182 extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 183 extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 184 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 185 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 186 extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 187 extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 188 extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 189 extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 190 extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 191 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 192 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 193 extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 194 extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 195 extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg); 196 extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 197 extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 198 extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 199 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 200 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 201 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 202 extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 203 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 204 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 205 extern PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 206 extern PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 207 extern PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 208 extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 209 extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 210 extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 211 extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 212 extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 213 extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 214 extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 215 extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 216 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 217 extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 218 219 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool ); 220 extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool ); 221 222 extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 223 224 EXTERN_C_BEGIN 225 extern PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 226 extern PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 227 extern PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*); 228 extern PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 229 extern PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 230 extern PetscErrorCode MatCreate_SeqAIJ(Mat); 231 EXTERN_C_END 232 extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 233 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 234 235 236 /* 237 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 238 239 Input Parameters: 240 + nnz - the number of entries 241 . r - the array of vector values 242 . xv - the matrix values for the row 243 - xi - the column indices of the nonzeros in the row 244 245 Output Parameter: 246 . sum - negative the sum of results 247 248 PETSc compile flags: 249 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 250 - PETSC_KERNEL_USE_UNROLL_2 - 251 252 .seealso: PetscSparseDensePlusDot() 253 254 */ 255 #ifdef PETSC_KERNEL_USE_UNROLL_4 256 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 257 if (nnz > 0) {\ 258 switch (nnz & 0x3) {\ 259 case 3: sum -= *xv++ * r[*xi++];\ 260 case 2: sum -= *xv++ * r[*xi++];\ 261 case 1: sum -= *xv++ * r[*xi++];\ 262 nnz -= 4;}\ 263 while (nnz > 0) {\ 264 sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 265 xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 266 xv += 4; xi += 4; nnz -= 4; }}} 267 268 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 269 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 270 PetscInt __i,__i1,__i2;\ 271 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 272 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 273 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 274 275 #else 276 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 277 PetscInt __i;\ 278 for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 279 #endif 280 281 282 283 /* 284 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 285 286 Input Parameters: 287 + nnz - the number of entries 288 . r - the array of vector values 289 . xv - the matrix values for the row 290 - xi - the column indices of the nonzeros in the row 291 292 Output Parameter: 293 . sum - the sum of results 294 295 PETSc compile flags: 296 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 297 - PETSC_KERNEL_USE_UNROLL_2 - 298 299 .seealso: PetscSparseDenseMinusDot() 300 301 */ 302 #ifdef PETSC_KERNEL_USE_UNROLL_4 303 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 304 if (nnz > 0) {\ 305 switch (nnz & 0x3) {\ 306 case 3: sum += *xv++ * r[*xi++];\ 307 case 2: sum += *xv++ * r[*xi++];\ 308 case 1: sum += *xv++ * r[*xi++];\ 309 nnz -= 4;}\ 310 while (nnz > 0) {\ 311 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 312 xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 313 xv += 4; xi += 4; nnz -= 4; }}} 314 315 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 316 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 317 PetscInt __i,__i1,__i2;\ 318 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 319 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 320 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 321 322 #else 323 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 324 PetscInt __i;\ 325 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 326 #endif 327 328 #endif 329