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 typedef struct { 87 MatTransposeColoring matcoloring; 88 Mat Bt_den; /* dense matrix of B^T */ 89 Mat ABt_den; /* dense matrix of A*B^T */ 90 PetscBool usecoloring; 91 PetscErrorCode (*destroy)(Mat); 92 } Mat_MatMatMultTrans; 93 94 typedef struct { 95 MatTransposeColoring matcoloring; 96 Mat Rt; /* dense matrix of R^T */ 97 Mat ARt; /* dense matrix of A*R^T */ 98 Mat RARt; /* dense matrix of R*A*R^T */ 99 PetscErrorCode (*destroy)(Mat); 100 } Mat_RARt; 101 102 /* 103 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 104 */ 105 #undef __FUNCT__ 106 #define __FUNCT__ "MatSeqXAIJFreeAIJ" 107 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 108 { 109 PetscErrorCode ierr; 110 Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 111 if (A->singlemalloc) { 112 ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 113 } else { 114 if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 115 if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 116 if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 117 } 118 return 0; 119 } 120 /* 121 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 122 This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 123 */ 124 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 125 if (NROW >= RMAX) {\ 126 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 127 /* there is no extra room in row, therefore enlarge */ \ 128 PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 129 datatype *new_a; \ 130 \ 131 if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 132 /* malloc new storage space */ \ 133 ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 134 \ 135 /* copy over old data into new slots */ \ 136 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 137 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 138 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 139 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 140 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 141 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 142 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 143 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 144 /* free up old matrix storage */ \ 145 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 146 AA = new_a; \ 147 Ain->a = (MatScalar*) new_a; \ 148 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 149 Ain->singlemalloc = PETSC_TRUE; \ 150 \ 151 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 152 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 153 Ain->maxnz += BS2*CHUNKSIZE; \ 154 Ain->reallocs++; \ 155 } \ 156 157 158 EXTERN_C_BEGIN 159 extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 160 EXTERN_C_END 161 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 162 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 163 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 164 165 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 166 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 167 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 168 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 169 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 170 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 171 extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 172 extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 173 extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*); 174 extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 175 176 extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 177 extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 178 extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 179 extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 180 extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 181 182 extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 183 extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 184 extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 185 186 extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 187 extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 188 extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 189 extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 190 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 191 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 192 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 193 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 194 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 195 extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 196 extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 197 extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 198 extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 199 extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 200 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 201 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 202 extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 203 extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 204 extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 205 extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 206 extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 207 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 208 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 209 extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 210 extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 211 extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg); 212 extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 213 extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 214 extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 215 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 216 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 217 218 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 219 extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 220 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 221 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 222 extern PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 223 extern PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 224 225 extern PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 226 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 227 extern PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 228 extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 229 extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 230 extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 231 extern PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 232 extern PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 233 extern PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 234 235 extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 236 extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 237 extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 238 extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 239 extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 240 extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 241 extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 242 extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 243 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 244 extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 245 246 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool ); 247 extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool ); 248 249 extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 250 251 EXTERN_C_BEGIN 252 extern PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 253 extern PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 254 extern PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*); 255 extern PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 256 extern PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 257 extern PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 258 extern PetscErrorCode MatCreate_SeqAIJ(Mat); 259 EXTERN_C_END 260 extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 261 extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 262 263 264 /* 265 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 266 267 Input Parameters: 268 + nnz - the number of entries 269 . r - the array of vector values 270 . xv - the matrix values for the row 271 - xi - the column indices of the nonzeros in the row 272 273 Output Parameter: 274 . sum - negative the sum of results 275 276 PETSc compile flags: 277 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 278 - PETSC_KERNEL_USE_UNROLL_2 - 279 280 .seealso: PetscSparseDensePlusDot() 281 282 */ 283 #ifdef PETSC_KERNEL_USE_UNROLL_4 284 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 285 if (nnz > 0) {\ 286 switch (nnz & 0x3) {\ 287 case 3: sum -= *xv++ * r[*xi++];\ 288 case 2: sum -= *xv++ * r[*xi++];\ 289 case 1: sum -= *xv++ * r[*xi++];\ 290 nnz -= 4;}\ 291 while (nnz > 0) {\ 292 sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 293 xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 294 xv += 4; xi += 4; nnz -= 4; }}} 295 296 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 297 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 298 PetscInt __i,__i1,__i2;\ 299 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 300 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 301 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 302 303 #else 304 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 305 PetscInt __i;\ 306 for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 307 #endif 308 309 310 311 /* 312 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 313 314 Input Parameters: 315 + nnz - the number of entries 316 . r - the array of vector values 317 . xv - the matrix values for the row 318 - xi - the column indices of the nonzeros in the row 319 320 Output Parameter: 321 . sum - the sum of results 322 323 PETSc compile flags: 324 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 325 - PETSC_KERNEL_USE_UNROLL_2 - 326 327 .seealso: PetscSparseDenseMinusDot() 328 329 */ 330 #ifdef PETSC_KERNEL_USE_UNROLL_4 331 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 332 if (nnz > 0) {\ 333 switch (nnz & 0x3) {\ 334 case 3: sum += *xv++ * r[*xi++];\ 335 case 2: sum += *xv++ * r[*xi++];\ 336 case 1: sum += *xv++ * r[*xi++];\ 337 nnz -= 4;}\ 338 while (nnz > 0) {\ 339 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 340 xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 341 xv += 4; xi += 4; nnz -= 4; }}} 342 343 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 344 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 345 PetscInt __i,__i1,__i2;\ 346 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 347 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 348 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 349 350 #else 351 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 352 PetscInt __i;\ 353 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 354 #endif 355 356 #endif 357