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