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 PetscTruth 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 PetscTruth 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 PetscTruth 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 PetscTruth keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/\ 23 PetscTruth ignorezeroentries; \ 24 PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 25 Mat XtoY; /* used by MatAXPY() */\ 26 PetscTruth free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 27 PetscTruth 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 PetscTruth 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 Mat parent /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); 38 means that this shares some data structures with the parent including diag, ilen, imax, i, j */ 39 40 /* 41 MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 42 format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 43 j[i[k]+p] is the pth column in row k. Note that the diagonal 44 matrix elements are stored with the rest of the nonzeros (not separately). 45 */ 46 47 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 48 typedef struct { 49 MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatRelax_Inode() */ 50 PetscInt bdiagsize; /* length of bdiag and ibdiag */ 51 PetscTruth ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 52 53 PetscTruth use; 54 PetscInt node_count; /* number of inodes */ 55 PetscInt *size; /* size of each inode */ 56 PetscInt limit; /* inode limit */ 57 PetscInt max_limit; /* maximum supported inode limit */ 58 PetscTruth checked; /* if inodes have been checked for */ 59 } Mat_Inode; 60 61 EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer); 62 EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType); 63 EXTERN PetscErrorCode MatDestroy_Inode(Mat); 64 EXTERN PetscErrorCode MatCreate_Inode(Mat); 65 EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth); 66 EXTERN PetscErrorCode MatDuplicate_Inode(Mat,MatDuplicateOption,Mat*); 67 EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat,IS,IS,const MatFactorInfo*,Mat*); 68 EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 69 EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 70 71 typedef struct { 72 SEQAIJHEADER(MatScalar); 73 Mat_Inode inode; 74 MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 75 76 PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 77 PetscTruth idiagvalid; /* current idiag[] and mdiag[] are valid */ 78 PetscScalar fshift,omega; /* last used omega and fshift */ 79 80 ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 81 } Mat_SeqAIJ; 82 83 /* 84 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 85 */ 86 #undef __FUNCT__ 87 #define __FUNCT__ "MatSeqXAIJFreeAIJ" 88 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 89 { 90 PetscErrorCode ierr; 91 Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 92 if (A->singlemalloc) { 93 ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 94 } else { 95 if (A->free_a && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 96 if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);} 97 if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);} 98 } 99 *a = 0; *j = 0; *i = 0; 100 return 0; 101 } 102 103 /* 104 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 105 */ 106 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 107 if (NROW >= RMAX) {\ 108 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 109 /* there is no extra room in row, therefore enlarge */ \ 110 PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 111 datatype *new_a; \ 112 \ 113 if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 114 /* malloc new storage space */ \ 115 ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 116 \ 117 /* copy over old data into new slots */ \ 118 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 119 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 120 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 121 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 122 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 123 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 124 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 125 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 126 /* free up old matrix storage */ \ 127 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 128 AA = new_a; \ 129 Ain->a = (MatScalar*) new_a; \ 130 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 131 Ain->singlemalloc = PETSC_TRUE; \ 132 \ 133 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 134 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 135 Ain->maxnz += CHUNKSIZE; \ 136 Ain->reallocs++; \ 137 } \ 138 139 EXTERN_C_BEGIN 140 EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*); 141 EXTERN_C_END 142 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 143 EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 144 EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 145 EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 146 EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 147 EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*); 148 EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 149 150 EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 151 EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 152 EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 153 EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 154 EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 155 156 EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 157 EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 158 EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 159 160 EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 161 EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 162 EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 163 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 164 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 165 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_newdatastruct(Mat,Mat,IS,IS,const MatFactorInfo*); 166 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 167 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_newdatastruct(Mat,Mat,const MatFactorInfo*); 168 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 169 EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 170 EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 171 EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 172 EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 173 EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 174 EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 175 EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 176 EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 177 EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 178 EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 179 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*,Mat*); 180 EXTERN PetscErrorCode MatILUDTFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 181 EXTERN PetscErrorCode MatILUDTFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 182 EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*); 183 EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 184 EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 185 EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 186 EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 187 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 188 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 189 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 190 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 191 EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 192 EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 193 EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 194 EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 195 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 196 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 197 EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 198 EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 199 EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 200 EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 201 EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 202 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 203 EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 204 205 EXTERN_C_BEGIN 206 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 207 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 208 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*); 209 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 210 EXTERN_C_END 211 212 /* 213 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 214 215 Input Parameters: 216 + nnz - the number of entries 217 . r - the array of vector values 218 . xv - the matrix values for the row 219 - xi - the column indices of the nonzeros in the row 220 221 Output Parameter: 222 . sum - negative the sum of results 223 224 PETSc compile flags: 225 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 226 - PETSC_KERNEL_USE_UNROLL_2 - 227 228 .seealso: PetscSparseDensePlusDot() 229 230 */ 231 #ifdef PETSC_KERNEL_USE_UNROLL_4 232 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 233 if (nnz > 0) {\ 234 switch (nnz & 0x3) {\ 235 case 3: sum -= *xv++ * r[*xi++];\ 236 case 2: sum -= *xv++ * r[*xi++];\ 237 case 1: sum -= *xv++ * r[*xi++];\ 238 nnz -= 4;}\ 239 while (nnz > 0) {\ 240 sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 241 xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 242 xv += 4; xi += 4; nnz -= 4; }}} 243 244 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 245 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 246 PetscInt __i,__i1,__i2;\ 247 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 248 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 249 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 250 251 #else 252 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 253 PetscInt __i;\ 254 for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 255 #endif 256 257 258 259 /* 260 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 261 262 Input Parameters: 263 + nnz - the number of entries 264 . r - the array of vector values 265 . xv - the matrix values for the row 266 - xi - the column indices of the nonzeros in the row 267 268 Output Parameter: 269 . sum - the sum of results 270 271 PETSc compile flags: 272 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 273 - PETSC_KERNEL_USE_UNROLL_2 - 274 275 .seealso: PetscSparseDenseMinusDot() 276 277 */ 278 #ifdef PETSC_KERNEL_USE_UNROLL_4 279 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 280 if (nnz > 0) {\ 281 switch (nnz & 0x3) {\ 282 case 3: sum += *xv++ * r[*xi++];\ 283 case 2: sum += *xv++ * r[*xi++];\ 284 case 1: sum += *xv++ * r[*xi++];\ 285 nnz -= 4;}\ 286 while (nnz > 0) {\ 287 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 288 xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 289 xv += 4; xi += 4; nnz -= 4; }}} 290 291 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 292 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 293 PetscInt __i,__i1,__i2;\ 294 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 295 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 296 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 297 298 #else 299 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 300 PetscInt __i;\ 301 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 302 #endif 303 304 #endif 305