1 2 #if !defined(__AIJ_H) 3 #define __AIJ_H 4 5 #include <petsc-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 typedef struct { 42 MatTransposeColoring matcoloring; 43 Mat Bt_den; /* dense matrix of B^T */ 44 Mat ABt_den; /* dense matrix of A*B^T */ 45 PetscBool usecoloring; 46 PetscErrorCode (*destroy)(Mat); 47 } Mat_MatMatTransMult; 48 49 typedef struct { 50 PetscInt *api,*apj; /* symbolic structure of A*P */ 51 PetscScalar *apa; /* temporary array for storing one row of A*P */ 52 PetscErrorCode (*destroy)(Mat); 53 } Mat_PtAP; 54 55 typedef struct { 56 MatTransposeColoring matcoloring; 57 Mat Rt; /* dense matrix of R^T */ 58 Mat RARt; /* dense matrix of R*A*R^T */ 59 MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 60 PetscErrorCode (*destroy)(Mat); 61 } Mat_RARt; 62 63 /* 64 MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 65 format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 66 j[i[k]+p] is the pth column in row k. Note that the diagonal 67 matrix elements are stored with the rest of the nonzeros (not separately). 68 */ 69 70 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 71 typedef struct { 72 MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 73 PetscInt bdiagsize; /* length of bdiag and ibdiag */ 74 PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 75 76 PetscBool use; 77 PetscInt node_count; /* number of inodes */ 78 PetscInt *size; /* size of each inode */ 79 PetscInt limit; /* inode limit */ 80 PetscInt max_limit; /* maximum supported inode limit */ 81 PetscBool checked; /* if inodes have been checked for */ 82 } Mat_SeqAIJ_Inode; 83 84 extern PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 85 extern PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 86 extern PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 87 extern PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 88 extern PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool ); 89 extern PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 90 extern PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool ); 91 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 92 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 93 94 typedef struct { 95 SEQAIJHEADER(MatScalar); 96 Mat_SeqAIJ_Inode inode; 97 MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 98 99 PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 100 PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 101 PetscScalar *ibdiag; /* inverses of block diagonals */ 102 PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 103 PetscScalar fshift,omega; /* last used omega and fshift */ 104 105 ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 106 107 PetscScalar *matmult_abdense; /* used by MatMatMult() */ 108 Mat_PtAP *ptap; /* used by MatPtAP() */ 109 } Mat_SeqAIJ; 110 111 /* 112 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 113 */ 114 #undef __FUNCT__ 115 #define __FUNCT__ "MatSeqXAIJFreeAIJ" 116 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 117 { 118 PetscErrorCode ierr; 119 Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 120 if (A->singlemalloc) { 121 ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 122 } else { 123 if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 124 if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 125 if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 126 } 127 return 0; 128 } 129 /* 130 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 131 This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 132 */ 133 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 134 if (NROW >= RMAX) {\ 135 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 136 /* there is no extra room in row, therefore enlarge */ \ 137 PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 138 datatype *new_a; \ 139 \ 140 if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 141 /* malloc new storage space */ \ 142 ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 143 \ 144 /* copy over old data into new slots */ \ 145 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 146 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 147 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 148 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 149 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 150 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 151 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 152 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 153 /* free up old matrix storage */ \ 154 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 155 AA = new_a; \ 156 Ain->a = (MatScalar*) new_a; \ 157 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 158 Ain->singlemalloc = PETSC_TRUE; \ 159 \ 160 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 161 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 162 Ain->maxnz += BS2*CHUNKSIZE; \ 163 Ain->reallocs++; \ 164 } \ 165 166 167 EXTERN_C_BEGIN 168 extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 169 EXTERN_C_END 170 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 171 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 172 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 173 174 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 175 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 176 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 177 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 178 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 179 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 180 extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 181 extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 182 extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*); 183 extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 184 185 extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 186 extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 187 extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 188 extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 189 extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 190 191 extern PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool); 192 extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 193 extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 194 extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 195 196 extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 197 extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 198 extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 199 extern PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*); 200 extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 201 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 202 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 203 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 204 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 205 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 206 extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 207 extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 208 extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 209 extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 210 extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 211 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 212 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 213 extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 214 extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 215 extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 216 extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 217 extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 218 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 219 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 220 extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 221 extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 222 extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg); 223 extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 224 extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 225 extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 226 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 227 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*); 228 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*); 229 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*); 230 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*); 231 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 232 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat); 233 234 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 235 extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 236 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 237 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(Mat,Mat,PetscReal,Mat*); 238 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 239 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat); 240 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2(Mat,Mat,Mat); 241 242 extern PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 243 extern PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 244 245 extern PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 246 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 247 extern PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 248 extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 249 extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 250 extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 251 extern PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 252 extern PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 253 extern PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 254 255 extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 256 extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 257 extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 258 extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 259 extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 260 extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 261 extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 262 extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 263 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 264 extern PetscErrorCode MatSetUp_SeqAIJ(Mat); 265 extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 266 267 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool ); 268 extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool ); 269 270 extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 271 272 EXTERN_C_BEGIN 273 extern PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 274 extern PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 275 extern PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*); 276 extern PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 277 extern PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 278 extern PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 279 extern PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 280 extern PetscErrorCode MatCreate_SeqAIJ(Mat); 281 EXTERN_C_END 282 extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 283 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 284 285 286 /* 287 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 288 289 Input Parameters: 290 + nnz - the number of entries 291 . r - the array of vector values 292 . xv - the matrix values for the row 293 - xi - the column indices of the nonzeros in the row 294 295 Output Parameter: 296 . sum - negative the sum of results 297 298 PETSc compile flags: 299 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 300 - PETSC_KERNEL_USE_UNROLL_2 - 301 302 .seealso: PetscSparseDensePlusDot() 303 304 */ 305 #ifdef PETSC_KERNEL_USE_UNROLL_4 306 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 307 if (nnz > 0) {\ 308 switch (nnz & 0x3) {\ 309 case 3: sum -= *xv++ * r[*xi++];\ 310 case 2: sum -= *xv++ * r[*xi++];\ 311 case 1: sum -= *xv++ * r[*xi++];\ 312 nnz -= 4;}\ 313 while (nnz > 0) {\ 314 sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 315 xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 316 xv += 4; xi += 4; nnz -= 4; }}} 317 318 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 319 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 320 PetscInt __i,__i1,__i2;\ 321 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 322 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 323 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 324 325 #else 326 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 327 PetscInt __i;\ 328 for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 329 #endif 330 331 332 333 /* 334 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 335 336 Input Parameters: 337 + nnz - the number of entries 338 . r - the array of vector values 339 . xv - the matrix values for the row 340 - xi - the column indices of the nonzeros in the row 341 342 Output Parameter: 343 . sum - the sum of results 344 345 PETSc compile flags: 346 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 347 - PETSC_KERNEL_USE_UNROLL_2 - 348 349 .seealso: PetscSparseDenseMinusDot() 350 351 */ 352 #ifdef PETSC_KERNEL_USE_UNROLL_4 353 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 354 if (nnz > 0) {\ 355 switch (nnz & 0x3) {\ 356 case 3: sum += *xv++ * r[*xi++];\ 357 case 2: sum += *xv++ * r[*xi++];\ 358 case 1: sum += *xv++ * r[*xi++];\ 359 nnz -= 4;}\ 360 while (nnz > 0) {\ 361 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 362 xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 363 xv += 4; xi += 4; nnz -= 4; }}} 364 365 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 366 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 367 PetscInt __i,__i1,__i2;\ 368 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 369 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 370 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 371 372 #else 373 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 374 PetscInt __i;\ 375 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 376 #endif 377 378 #endif 379