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 */ \ 364fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 37baabb860SHong Zhang PetscTruth pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 384fd072dbSBarry Smith Mat parent /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); 394fd072dbSBarry Smith means that this shares some data structures with the parent including diag, ilen, imax, i, j */ 40e6b907acSBarry Smith 412d40f771SBarry Smith /* 42ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 43e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 44dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 455768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 462d40f771SBarry Smith */ 47d35516d3SLois Curfman McInnes 48e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 49b8a66259SBarry Smith typedef struct { 504108e4d5SBarry Smith MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 51f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 5271f1c65dSBarry Smith PetscTruth ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 53f0d39aaaSBarry Smith 54e6b907acSBarry Smith PetscTruth use; 55e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 56e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 57e6b907acSBarry Smith PetscInt limit; /* inode limit */ 58e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 59e6b907acSBarry Smith PetscTruth checked; /* if inodes have been checked for */ 604108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 61e6b907acSBarry Smith 624108e4d5SBarry Smith EXTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 634108e4d5SBarry Smith EXTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 644108e4d5SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 654108e4d5SBarry Smith EXTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 664108e4d5SBarry Smith EXTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscTruth); 674108e4d5SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 68628f99d7SShri Abhyankar EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 6928f1b45aSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 70e6b907acSBarry Smith 71e6b907acSBarry Smith typedef struct { 7254f21887SBarry Smith SEQAIJHEADER(MatScalar); 734108e4d5SBarry Smith Mat_SeqAIJ_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 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 104357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 105e6b907acSBarry Smith */ 106*fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 107*fef13f97SBarry Smith if (NROW >= RMAX) {\ 108*fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 109*fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 110*fef13f97SBarry Smith PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 111*fef13f97SBarry Smith datatype *new_a; \ 112*fef13f97SBarry Smith \ 113*fef13f97SBarry Smith if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 114*fef13f97SBarry Smith /* malloc new storage space */ \ 115*fef13f97SBarry Smith ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 116*fef13f97SBarry Smith \ 117*fef13f97SBarry Smith /* copy over old data into new slots */ \ 118*fef13f97SBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 119*fef13f97SBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 120*fef13f97SBarry Smith ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 121*fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 122*fef13f97SBarry Smith ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 123*fef13f97SBarry Smith ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 124*fef13f97SBarry Smith ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 125*fef13f97SBarry Smith ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 126*fef13f97SBarry Smith /* free up old matrix storage */ \ 127*fef13f97SBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 128*fef13f97SBarry Smith AA = new_a; \ 129*fef13f97SBarry Smith Ain->a = (MatScalar*) new_a; \ 130*fef13f97SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 131*fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 132*fef13f97SBarry Smith \ 133*fef13f97SBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 134*fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 135*fef13f97SBarry Smith Ain->maxnz += BS2*CHUNKSIZE; \ 136*fef13f97SBarry Smith Ain->reallocs++; \ 137*fef13f97SBarry Smith } \ 13817454e89SShri Abhyankar 139e6b907acSBarry Smith 140f349cc81SSatish Balay EXTERN_C_BEGIN 141ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*); 142f349cc81SSatish Balay EXTERN_C_END 143ad04f41aSHong Zhang EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 14435233d99SShri Abhyankar EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 14535233d99SShri Abhyankar EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 1461df811f5SHong Zhang 147ad04f41aSHong Zhang EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 14835233d99SShri Abhyankar EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 149ad04f41aSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 15035233d99SShri Abhyankar EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 151ad04f41aSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 15235233d99SShri Abhyankar EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 153dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 154cae5a23dSHong Zhang EXTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 15509f38230SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*); 156dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 15708480c60SBarry Smith 158dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 159dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 160dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 161dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 16241f059aeSBarry Smith EXTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 16308480c60SBarry Smith 164dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 165dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 16638baddfdSBarry Smith EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 1673a7fca6bSBarry Smith 16838baddfdSBarry Smith EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 1690390132cSHong Zhang EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 17038baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 17138baddfdSBarry Smith EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 172ad04f41aSHong Zhang EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 17335233d99SShri Abhyankar EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 174ad04f41aSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 17535233d99SShri Abhyankar EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 1760481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 1770481f469SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 178ad04f41aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 17935233d99SShri Abhyankar EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 180ad04f41aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 18135233d99SShri Abhyankar EXTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 182ad04f41aSHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 18335233d99SShri Abhyankar EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 184137fb511SHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 185ad04f41aSHong Zhang EXTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 18635233d99SShri Abhyankar EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 187ad04f41aSHong Zhang EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 18835233d99SShri Abhyankar EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 189ad04f41aSHong Zhang EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 19035233d99SShri Abhyankar EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 191ad04f41aSHong Zhang EXTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 19235233d99SShri Abhyankar EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 193dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 194dfbe8321SBarry Smith EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 195a313700dSBarry Smith EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*); 196dfbe8321SBarry Smith EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 197dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 198dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 1997ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 2007ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 201dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 202dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 203bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 204bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 205bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 20638baddfdSBarry Smith EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 207a77337e4SBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 208a77337e4SBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 209f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 2108f7157efSSatish Balay EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 2118f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 2128f7157efSSatish Balay EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 2138f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 214b24902e0SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 215b24902e0SBarry Smith EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 2169af31e4aSHong Zhang 217019b515eSShri Abhyankar EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth); 218019b515eSShri Abhyankar EXTERN PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscTruth); 219019b515eSShri Abhyankar 22097304618SKris Buschelman EXTERN_C_BEGIN 2218cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 2228cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 2238cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*); 224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 225ec01deb9SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 22697304618SKris Buschelman EXTERN_C_END 22770f19b1fSKris Buschelman 228003131ecSBarry Smith /* 229003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 230003131ecSBarry Smith 231003131ecSBarry Smith Input Parameters: 232003131ecSBarry Smith + nnz - the number of entries 233003131ecSBarry Smith . r - the array of vector values 234003131ecSBarry Smith . xv - the matrix values for the row 235003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 236003131ecSBarry Smith 237003131ecSBarry Smith Output Parameter: 238003131ecSBarry Smith . sum - negative the sum of results 239003131ecSBarry Smith 240003131ecSBarry Smith PETSc compile flags: 241b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 242e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 243003131ecSBarry Smith 244003131ecSBarry Smith .seealso: PetscSparseDensePlusDot() 245003131ecSBarry Smith 246003131ecSBarry Smith */ 247003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4 248003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 249003131ecSBarry Smith if (nnz > 0) {\ 250003131ecSBarry Smith switch (nnz & 0x3) {\ 251003131ecSBarry Smith case 3: sum -= *xv++ * r[*xi++];\ 252003131ecSBarry Smith case 2: sum -= *xv++ * r[*xi++];\ 253003131ecSBarry Smith case 1: sum -= *xv++ * r[*xi++];\ 254003131ecSBarry Smith nnz -= 4;}\ 255003131ecSBarry Smith while (nnz > 0) {\ 256003131ecSBarry Smith sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 257003131ecSBarry Smith xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 258003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 259003131ecSBarry Smith 260003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 261003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 262003131ecSBarry Smith PetscInt __i,__i1,__i2;\ 263003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 264003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 265003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 266003131ecSBarry Smith 267003131ecSBarry Smith #else 268003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 269003131ecSBarry Smith PetscInt __i;\ 270003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 271003131ecSBarry Smith #endif 272003131ecSBarry Smith 273003131ecSBarry Smith 274003131ecSBarry Smith 275003131ecSBarry Smith /* 276003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 277003131ecSBarry Smith 278003131ecSBarry Smith Input Parameters: 279003131ecSBarry Smith + nnz - the number of entries 280003131ecSBarry Smith . r - the array of vector values 281003131ecSBarry Smith . xv - the matrix values for the row 282003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 283003131ecSBarry Smith 284003131ecSBarry Smith Output Parameter: 285003131ecSBarry Smith . sum - the sum of results 286003131ecSBarry Smith 287003131ecSBarry Smith PETSc compile flags: 288b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 289e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 290003131ecSBarry Smith 291003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot() 292003131ecSBarry Smith 293003131ecSBarry Smith */ 294003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4 295003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 296003131ecSBarry Smith if (nnz > 0) {\ 297003131ecSBarry Smith switch (nnz & 0x3) {\ 298003131ecSBarry Smith case 3: sum += *xv++ * r[*xi++];\ 299003131ecSBarry Smith case 2: sum += *xv++ * r[*xi++];\ 300003131ecSBarry Smith case 1: sum += *xv++ * r[*xi++];\ 301003131ecSBarry Smith nnz -= 4;}\ 302003131ecSBarry Smith while (nnz > 0) {\ 303003131ecSBarry Smith sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 304003131ecSBarry Smith xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 305003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 306003131ecSBarry Smith 307003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 308003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 309003131ecSBarry Smith PetscInt __i,__i1,__i2;\ 310003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 311003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 312003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 313003131ecSBarry Smith 314003131ecSBarry Smith #else 315003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 316003131ecSBarry Smith PetscInt __i;\ 317003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 318003131ecSBarry Smith #endif 319003131ecSBarry Smith 3202d40f771SBarry Smith #endif 321