1b8a66259SBarry Smith 22d40f771SBarry Smith #if !defined(__AIJ_H) 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 57c4f633dSBarry Smith #include "private/matimpl.h" 68f690400SShri Abhyankar 7e6b907acSBarry Smith /* 8e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 9e6b907acSBarry Smith */ 10421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 11e6b907acSBarry Smith PetscTruth roworiented; /* if true, row-oriented input, default */\ 12e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */\ 1328b2fa4aSMatthew Knepley PetscInt nounused; /* -1 generate error on unused space */\ 14e6b907acSBarry Smith PetscTruth singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\ 15e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */\ 16e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */\ 17e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */\ 18c760cd28SBarry Smith PetscTruth free_imax_ilen; \ 19e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 20e6b907acSBarry Smith as more values are set than were prealloced */\ 21e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */\ 22a9817697SBarry Smith PetscTruth keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/\ 23e6b907acSBarry Smith PetscTruth ignorezeroentries; \ 24e6b907acSBarry Smith PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 25e6b907acSBarry Smith Mat XtoY; /* used by MatAXPY() */\ 26e6b907acSBarry Smith PetscTruth free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 27e6b907acSBarry Smith PetscTruth free_a; /* free the numerical values when matrix is destroy */ \ 28e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 29e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 30e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 31e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 32e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 33c760cd28SBarry Smith PetscTruth free_diag; \ 34421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 35e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 36*4fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 37*4fd072dbSBarry Smith Mat parent /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); 38*4fd072dbSBarry Smith means that this shares some data structures with the parent including diag, ilen, imax, i, j */ 39e6b907acSBarry Smith 402d40f771SBarry Smith /* 41ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 42e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 43dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 445768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 452d40f771SBarry Smith */ 46d35516d3SLois Curfman McInnes 47e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 48b8a66259SBarry Smith typedef struct { 4989c6957cSBarry Smith MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatRelax_Inode() */ 50f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 5171f1c65dSBarry Smith PetscTruth ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 52f0d39aaaSBarry Smith 53e6b907acSBarry Smith PetscTruth use; 54e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 55e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 56e6b907acSBarry Smith PetscInt limit; /* inode limit */ 57e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 58e6b907acSBarry Smith PetscTruth checked; /* if inodes have been checked for */ 59e6b907acSBarry Smith } Mat_Inode; 60e6b907acSBarry Smith 61e6b907acSBarry Smith EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer); 62e6b907acSBarry Smith EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType); 63e6b907acSBarry Smith EXTERN PetscErrorCode MatDestroy_Inode(Mat); 64e6b907acSBarry Smith EXTERN PetscErrorCode MatCreate_Inode(Mat); 654e0d8c25SBarry Smith EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth); 66719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicate_Inode(Mat,MatDuplicateOption,Mat*); 670481f469SBarry Smith EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat,IS,IS,const MatFactorInfo*,Mat*); 680481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 690481f469SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 70e6b907acSBarry Smith 71e6b907acSBarry Smith typedef struct { 7254f21887SBarry Smith SEQAIJHEADER(MatScalar); 73e6b907acSBarry Smith Mat_Inode inode; 7454f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 7571f1c65dSBarry Smith 7671f1c65dSBarry Smith PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 7771f1c65dSBarry Smith PetscTruth idiagvalid; /* current idiag[] and mdiag[] are valid */ 7871f1c65dSBarry Smith PetscScalar fshift,omega; /* last used omega and fshift */ 7971f1c65dSBarry Smith 803a7fca6bSBarry Smith ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 81ec8511deSBarry Smith } Mat_SeqAIJ; 82b8a66259SBarry Smith 83e6b907acSBarry Smith /* 84e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 85e6b907acSBarry Smith */ 86e6b907acSBarry Smith #undef __FUNCT__ 87e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ" 88d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 89d0f46423SBarry Smith { 90e6b907acSBarry Smith PetscErrorCode ierr; 91e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 92e6b907acSBarry Smith if (A->singlemalloc) { 93e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 94e6b907acSBarry Smith } else { 95e6b907acSBarry Smith if (A->free_a && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 96e6b907acSBarry Smith if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);} 97e6b907acSBarry Smith if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);} 98e6b907acSBarry Smith } 99e6b907acSBarry Smith *a = 0; *j = 0; *i = 0; 100e6b907acSBarry Smith return 0; 101e6b907acSBarry Smith } 102e6b907acSBarry Smith 103e6b907acSBarry Smith /* 104e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 105e6b907acSBarry Smith */ 106421e10b8SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 107e6b907acSBarry Smith if (NROW >= RMAX) {\ 108e6b907acSBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 109e6b907acSBarry Smith /* there is no extra room in row, therefore enlarge */ \ 110e6b907acSBarry Smith PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 111421e10b8SBarry Smith datatype *new_a; \ 112e6b907acSBarry Smith \ 113e6b907acSBarry Smith if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 114e6b907acSBarry Smith /* malloc new storage space */ \ 115421e10b8SBarry Smith ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 116e6b907acSBarry Smith \ 117e6b907acSBarry Smith /* copy over old data into new slots */ \ 118e6b907acSBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 119e6b907acSBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 120e6b907acSBarry Smith ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 121e6b907acSBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 122e6b907acSBarry Smith ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 123421e10b8SBarry Smith ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 124421e10b8SBarry Smith ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 125421e10b8SBarry Smith ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 126e6b907acSBarry Smith /* free up old matrix storage */ \ 127e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 128421e10b8SBarry Smith AA = new_a; \ 129421e10b8SBarry Smith Ain->a = (MatScalar*) new_a; \ 130421e10b8SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 131e6b907acSBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 132e6b907acSBarry Smith \ 133e6b907acSBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 134e6b907acSBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 135e6b907acSBarry Smith Ain->maxnz += CHUNKSIZE; \ 136e6b907acSBarry Smith Ain->reallocs++; \ 137e6b907acSBarry Smith } \ 138e6b907acSBarry Smith 139f349cc81SSatish Balay EXTERN_C_BEGIN 140ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*); 141f349cc81SSatish Balay EXTERN_C_END 1420481f469SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 1430481f469SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 1440481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 1450481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 146dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 14709f38230SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*); 148dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 14908480c60SBarry Smith 150dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 151dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 152dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 153dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 15438baddfdSBarry Smith EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 15508480c60SBarry Smith 156dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 157dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 15838baddfdSBarry Smith EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 1593a7fca6bSBarry Smith 16038baddfdSBarry Smith EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 1610390132cSHong Zhang EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 16238baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 16338baddfdSBarry Smith EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 1640481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 165ce3d78c0SShri Abhyankar EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_newdatastruct(Mat,Mat,IS,IS,const MatFactorInfo*); 1660481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 167ce3d78c0SShri Abhyankar EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_newdatastruct(Mat,Mat,const MatFactorInfo*); 1680481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 1690481f469SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 170dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 1718d8c7c43SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 172137fb511SHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 173dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 174dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 175dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 1762bebee5dSHong Zhang EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 177dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 178dfbe8321SBarry Smith EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1790481f469SBarry Smith EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*,Mat*); 1803c2a7987SHong Zhang EXTERN PetscErrorCode MatILUDTFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 1813c2a7987SHong Zhang EXTERN PetscErrorCode MatILUDTFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 182a313700dSBarry Smith EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*); 183dfbe8321SBarry Smith EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 184dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 185dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 186dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 1877ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 1887ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 189dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 190dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 191bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 192bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 193bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 19438baddfdSBarry Smith EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 195a77337e4SBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 196a77337e4SBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 197f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 1988f7157efSSatish Balay EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 1998f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 2008f7157efSSatish Balay EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 2018f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 202b24902e0SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 203b24902e0SBarry Smith EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 2049af31e4aSHong Zhang 20597304618SKris Buschelman EXTERN_C_BEGIN 2068cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 2078cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 2088cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*); 209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 21097304618SKris Buschelman EXTERN_C_END 21170f19b1fSKris Buschelman 212003131ecSBarry Smith /* 213003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 214003131ecSBarry Smith 215003131ecSBarry Smith Input Parameters: 216003131ecSBarry Smith + nnz - the number of entries 217003131ecSBarry Smith . r - the array of vector values 218003131ecSBarry Smith . xv - the matrix values for the row 219003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 220003131ecSBarry Smith 221003131ecSBarry Smith Output Parameter: 222003131ecSBarry Smith . sum - negative the sum of results 223003131ecSBarry Smith 224003131ecSBarry Smith PETSc compile flags: 225b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 226e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 227003131ecSBarry Smith 228003131ecSBarry Smith .seealso: PetscSparseDensePlusDot() 229003131ecSBarry Smith 230003131ecSBarry Smith */ 231003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4 232003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 233003131ecSBarry Smith if (nnz > 0) {\ 234003131ecSBarry Smith switch (nnz & 0x3) {\ 235003131ecSBarry Smith case 3: sum -= *xv++ * r[*xi++];\ 236003131ecSBarry Smith case 2: sum -= *xv++ * r[*xi++];\ 237003131ecSBarry Smith case 1: sum -= *xv++ * r[*xi++];\ 238003131ecSBarry Smith nnz -= 4;}\ 239003131ecSBarry Smith while (nnz > 0) {\ 240003131ecSBarry Smith sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 241003131ecSBarry Smith xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 242003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 243003131ecSBarry Smith 244003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 245003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 246003131ecSBarry Smith PetscInt __i,__i1,__i2;\ 247003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 248003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 249003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 250003131ecSBarry Smith 251003131ecSBarry Smith #else 252003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 253003131ecSBarry Smith PetscInt __i;\ 254003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 255003131ecSBarry Smith #endif 256003131ecSBarry Smith 257003131ecSBarry Smith 258003131ecSBarry Smith 259003131ecSBarry Smith /* 260003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 261003131ecSBarry Smith 262003131ecSBarry Smith Input Parameters: 263003131ecSBarry Smith + nnz - the number of entries 264003131ecSBarry Smith . r - the array of vector values 265003131ecSBarry Smith . xv - the matrix values for the row 266003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 267003131ecSBarry Smith 268003131ecSBarry Smith Output Parameter: 269003131ecSBarry Smith . sum - the sum of results 270003131ecSBarry Smith 271003131ecSBarry Smith PETSc compile flags: 272b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 273e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 274003131ecSBarry Smith 275003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot() 276003131ecSBarry Smith 277003131ecSBarry Smith */ 278003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4 279003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 280003131ecSBarry Smith if (nnz > 0) {\ 281003131ecSBarry Smith switch (nnz & 0x3) {\ 282003131ecSBarry Smith case 3: sum += *xv++ * r[*xi++];\ 283003131ecSBarry Smith case 2: sum += *xv++ * r[*xi++];\ 284003131ecSBarry Smith case 1: sum += *xv++ * r[*xi++];\ 285003131ecSBarry Smith nnz -= 4;}\ 286003131ecSBarry Smith while (nnz > 0) {\ 287003131ecSBarry Smith sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 288003131ecSBarry Smith xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 289003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 290003131ecSBarry Smith 291003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 292003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 293003131ecSBarry Smith PetscInt __i,__i1,__i2;\ 294003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 295003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 296003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 297003131ecSBarry Smith 298003131ecSBarry Smith #else 299003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 300003131ecSBarry Smith PetscInt __i;\ 301003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 302003131ecSBarry Smith #endif 303003131ecSBarry Smith 3042d40f771SBarry Smith #endif 305