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