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 85 PetscInt matmult_denseaxpy; /* used by MatMatMult(): <=0: use sparse axpy; otherwise: num of dense rows used in MatMatMultNumeric() */ 86 PetscScalar *matmult_abdense; /* used by MatMatMult() */ 87 } Mat_SeqAIJ; 88 89 typedef struct { 90 MatTransposeColoring matcoloring; 91 Mat Bt_den; /* dense matrix of B^T */ 92 Mat ABt_den; /* dense matrix of A*B^T */ 93 PetscBool usecoloring; 94 PetscErrorCode (*destroy)(Mat); 95 } Mat_MatMatTransMult; 96 97 typedef struct { 98 PetscInt *api,*apj; /* symbolic structure of A*P */ 99 PetscScalar *apa; /* temporary array for storing one row of A*P */ 100 PetscErrorCode (*destroy)(Mat); 101 } Mat_PtAP; 102 103 typedef struct { 104 MatTransposeColoring matcoloring; 105 Mat Rt; /* dense matrix of R^T */ 106 Mat RARt; /* dense matrix of R*A*R^T */ 107 MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 108 PetscErrorCode (*destroy)(Mat); 109 } Mat_RARt; 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 MatSetColoring_SeqAIJ(Mat,ISColoring); 192 extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 193 extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 194 195 extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 196 extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 197 extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 198 extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 199 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 200 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 201 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 202 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 203 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 204 extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 205 extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 206 extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 207 extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 208 extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 209 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 210 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 211 extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 212 extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 213 extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 214 extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 215 extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 216 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 217 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 218 extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 219 extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 220 extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg); 221 extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 222 extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 223 extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 224 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 225 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 226 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat); 227 extern PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscReal,PetscInt*[],PetscInt*[],PetscInt*); 228 229 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 230 extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 231 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 232 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 233 extern PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 234 extern PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 235 236 extern PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 237 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 238 extern PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 239 extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 240 extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 241 extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 242 extern PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 243 extern PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 244 extern PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 245 246 extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 247 extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 248 extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 249 extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 250 extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 251 extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 252 extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 253 extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 254 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 255 extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 256 257 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool ); 258 extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool ); 259 260 extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 261 262 EXTERN_C_BEGIN 263 extern PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 264 extern PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 265 extern PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*); 266 extern PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 267 extern PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 268 extern PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 269 extern PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 270 extern PetscErrorCode MatCreate_SeqAIJ(Mat); 271 EXTERN_C_END 272 extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 273 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 274 275 276 /* 277 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 278 279 Input Parameters: 280 + nnz - the number of entries 281 . r - the array of vector values 282 . xv - the matrix values for the row 283 - xi - the column indices of the nonzeros in the row 284 285 Output Parameter: 286 . sum - negative the sum of results 287 288 PETSc compile flags: 289 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 290 - PETSC_KERNEL_USE_UNROLL_2 - 291 292 .seealso: PetscSparseDensePlusDot() 293 294 */ 295 #ifdef PETSC_KERNEL_USE_UNROLL_4 296 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 297 if (nnz > 0) {\ 298 switch (nnz & 0x3) {\ 299 case 3: sum -= *xv++ * r[*xi++];\ 300 case 2: sum -= *xv++ * r[*xi++];\ 301 case 1: sum -= *xv++ * r[*xi++];\ 302 nnz -= 4;}\ 303 while (nnz > 0) {\ 304 sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 305 xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 306 xv += 4; xi += 4; nnz -= 4; }}} 307 308 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 309 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 310 PetscInt __i,__i1,__i2;\ 311 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 312 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 313 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 314 315 #else 316 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 317 PetscInt __i;\ 318 for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 319 #endif 320 321 322 323 /* 324 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 325 326 Input Parameters: 327 + nnz - the number of entries 328 . r - the array of vector values 329 . xv - the matrix values for the row 330 - xi - the column indices of the nonzeros in the row 331 332 Output Parameter: 333 . sum - the sum of results 334 335 PETSc compile flags: 336 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 337 - PETSC_KERNEL_USE_UNROLL_2 - 338 339 .seealso: PetscSparseDenseMinusDot() 340 341 */ 342 #ifdef PETSC_KERNEL_USE_UNROLL_4 343 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 344 if (nnz > 0) {\ 345 switch (nnz & 0x3) {\ 346 case 3: sum += *xv++ * r[*xi++];\ 347 case 2: sum += *xv++ * r[*xi++];\ 348 case 1: sum += *xv++ * r[*xi++];\ 349 nnz -= 4;}\ 350 while (nnz > 0) {\ 351 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 352 xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 353 xv += 4; xi += 4; nnz -= 4; }}} 354 355 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 356 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 357 PetscInt __i,__i1,__i2;\ 358 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 359 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 360 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 361 362 #else 363 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 364 PetscInt __i;\ 365 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 366 #endif 367 368 #endif 369