1b8a66259SBarry Smith 22d40f771SBarry Smith #if !defined(__AIJ_H) 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5c6db04a5SJed Brown #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) \ 11ace3abfcSBarry Smith PetscBool 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 */\ 14ace3abfcSBarry Smith PetscBool 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 */\ 18ace3abfcSBarry Smith PetscBool 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 */\ 22ace3abfcSBarry Smith PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/\ 23ace3abfcSBarry Smith PetscBool ignorezeroentries; \ 24e6b907acSBarry Smith PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 25e6b907acSBarry Smith Mat XtoY; /* used by MatAXPY() */\ 26ace3abfcSBarry Smith PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 27ace3abfcSBarry Smith PetscBool 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 */ \ 33ace3abfcSBarry Smith PetscBool 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 */ \ 37ace3abfcSBarry Smith PetscBool 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 */ 52ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 53f0d39aaaSBarry Smith 54ace3abfcSBarry Smith PetscBool 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 */ 59ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 604108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 61e6b907acSBarry Smith 6209573ac7SBarry Smith extern PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 6309573ac7SBarry Smith extern PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 6409573ac7SBarry Smith extern PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 6509573ac7SBarry Smith extern PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 6609573ac7SBarry Smith extern PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool ); 6709573ac7SBarry Smith extern PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 6809573ac7SBarry Smith extern PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool ); 6909573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 7009573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 71e6b907acSBarry Smith 72e6b907acSBarry Smith typedef struct { 7354f21887SBarry Smith SEQAIJHEADER(MatScalar); 744108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 7554f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 7671f1c65dSBarry Smith 7771f1c65dSBarry Smith PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 78ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 79bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 802291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 8171f1c65dSBarry Smith PetscScalar fshift,omega; /* last used omega and fshift */ 8271f1c65dSBarry Smith 833a7fca6bSBarry Smith ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 84ec8511deSBarry Smith } Mat_SeqAIJ; 85b8a66259SBarry Smith 86*2128a86cSHong Zhang typedef struct { 87*2128a86cSHong Zhang MatMultTransposeColoring matcoloring; 88*2128a86cSHong Zhang Mat Bt_den; /* dense matrix of B^T */ 89*2128a86cSHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 90*2128a86cSHong Zhang PetscBool usecoloring; 91*2128a86cSHong Zhang PetscErrorCode (*destroy)(Mat); 92*2128a86cSHong Zhang } Mat_MatMatMultTrans; 93*2128a86cSHong Zhang 94e6b907acSBarry Smith /* 95e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 96e6b907acSBarry Smith */ 97e6b907acSBarry Smith #undef __FUNCT__ 98e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ" 99d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 100d0f46423SBarry Smith { 101e6b907acSBarry Smith PetscErrorCode ierr; 102e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 103e6b907acSBarry Smith if (A->singlemalloc) { 104e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 105e6b907acSBarry Smith } else { 106c31cb41cSBarry Smith if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 107c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 108c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 109e6b907acSBarry Smith } 110e6b907acSBarry Smith return 0; 111e6b907acSBarry Smith } 112e6b907acSBarry Smith /* 113e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 114357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 115e6b907acSBarry Smith */ 116fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 117fef13f97SBarry Smith if (NROW >= RMAX) {\ 118fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 119fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 120fef13f97SBarry Smith PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 121fef13f97SBarry Smith datatype *new_a; \ 122fef13f97SBarry Smith \ 123fef13f97SBarry Smith if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 124fef13f97SBarry Smith /* malloc new storage space */ \ 125fef13f97SBarry Smith ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 126fef13f97SBarry Smith \ 127fef13f97SBarry Smith /* copy over old data into new slots */ \ 128fef13f97SBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 129fef13f97SBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 130fef13f97SBarry Smith ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 131fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 132fef13f97SBarry Smith ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 133fef13f97SBarry Smith ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 134fef13f97SBarry Smith ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 135fef13f97SBarry Smith ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 136fef13f97SBarry Smith /* free up old matrix storage */ \ 137fef13f97SBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 138fef13f97SBarry Smith AA = new_a; \ 139fef13f97SBarry Smith Ain->a = (MatScalar*) new_a; \ 140fef13f97SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 141fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 142fef13f97SBarry Smith \ 143fef13f97SBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 144fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 145fef13f97SBarry Smith Ain->maxnz += BS2*CHUNKSIZE; \ 146fef13f97SBarry Smith Ain->reallocs++; \ 147fef13f97SBarry Smith } \ 14817454e89SShri Abhyankar 149e6b907acSBarry Smith 150f349cc81SSatish Balay EXTERN_C_BEGIN 15109573ac7SBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 152f349cc81SSatish Balay EXTERN_C_END 15309573ac7SBarry Smith extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 15409573ac7SBarry Smith extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 15509573ac7SBarry Smith extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 1561df811f5SHong Zhang 15709573ac7SBarry Smith extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 15809573ac7SBarry Smith extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 15909573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 16009573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 16109573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 16209573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 16309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 16409573ac7SBarry Smith extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 16509573ac7SBarry Smith extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*); 16609573ac7SBarry Smith extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 16708480c60SBarry Smith 16809573ac7SBarry Smith extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 16909573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 17009573ac7SBarry Smith extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 17109573ac7SBarry Smith extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 17209573ac7SBarry Smith extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 17308480c60SBarry Smith 17409573ac7SBarry Smith extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 17509573ac7SBarry Smith extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 17609573ac7SBarry Smith extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 1773a7fca6bSBarry Smith 17809573ac7SBarry Smith extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 17909573ac7SBarry Smith extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 18009573ac7SBarry Smith extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 18109573ac7SBarry Smith extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 18209573ac7SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 18309573ac7SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 18409573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 18509573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 18609573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 18709573ac7SBarry Smith extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 18809573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 18909573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 19009573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 19109573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 19209573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 19309573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 19409573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 19509573ac7SBarry Smith extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 19609573ac7SBarry Smith extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 19709573ac7SBarry Smith extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 19809573ac7SBarry Smith extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 19909573ac7SBarry Smith extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 20009573ac7SBarry Smith extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 20109573ac7SBarry Smith extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 20209573ac7SBarry Smith extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 20309573ac7SBarry Smith extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg); 20409573ac7SBarry Smith extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 20509573ac7SBarry Smith extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 20609573ac7SBarry Smith extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 20709573ac7SBarry Smith extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 20809573ac7SBarry Smith extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 20909573ac7SBarry Smith extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 21009573ac7SBarry Smith extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 21109573ac7SBarry Smith extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 21209573ac7SBarry Smith extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 2135df89d91SHong Zhang 2145df89d91SHong Zhang extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2155df89d91SHong Zhang extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 2165df89d91SHong Zhang extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 21709573ac7SBarry Smith extern PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 21809573ac7SBarry Smith extern PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 21909573ac7SBarry Smith extern PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 220b1683b59SHong Zhang extern PetscErrorCode MatMultTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatMultTransposeColoring); 221c8db22f6SHong Zhang extern PetscErrorCode MatMultTransposeColoringApply_SeqAIJ(Mat,Mat,MatMultTransposeColoring); 2228972f759SHong Zhang extern PetscErrorCode MatMultTransColoringApplyDenToSp_SeqAIJ(MatMultTransposeColoring,Mat,Mat); 2235df89d91SHong Zhang 22409573ac7SBarry Smith extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 22509573ac7SBarry Smith extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 22609573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 22709573ac7SBarry Smith extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 22809573ac7SBarry Smith extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 22909573ac7SBarry Smith extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 23009573ac7SBarry Smith extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 23109573ac7SBarry Smith extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 23209573ac7SBarry Smith extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 23309573ac7SBarry Smith extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 2349af31e4aSHong Zhang 23509573ac7SBarry Smith extern PetscErrorCode Mat_CheckInode(Mat,PetscBool ); 23609573ac7SBarry Smith extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool ); 237019b515eSShri Abhyankar 23809573ac7SBarry Smith extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 2399f5f6813SShri Abhyankar 24097304618SKris Buschelman EXTERN_C_BEGIN 2417087cfbeSBarry Smith extern PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 2427087cfbeSBarry Smith extern PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 2437087cfbeSBarry Smith extern PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*); 2447087cfbeSBarry Smith extern PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 2457087cfbeSBarry Smith extern PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2467087cfbeSBarry Smith extern PetscErrorCode MatCreate_SeqAIJ(Mat); 24797304618SKris Buschelman EXTERN_C_END 2487087cfbeSBarry Smith extern PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 2497087cfbeSBarry Smith extern PetscErrorCode MatDestroy_SeqAIJ(Mat); 2502f947c57SVictor Minden 25170f19b1fSKris Buschelman 252003131ecSBarry Smith /* 253003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 254003131ecSBarry Smith 255003131ecSBarry Smith Input Parameters: 256003131ecSBarry Smith + nnz - the number of entries 257003131ecSBarry Smith . r - the array of vector values 258003131ecSBarry Smith . xv - the matrix values for the row 259003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 260003131ecSBarry Smith 261003131ecSBarry Smith Output Parameter: 262003131ecSBarry Smith . sum - negative the sum of results 263003131ecSBarry Smith 264003131ecSBarry Smith PETSc compile flags: 265b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 266e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 267003131ecSBarry Smith 268003131ecSBarry Smith .seealso: PetscSparseDensePlusDot() 269003131ecSBarry Smith 270003131ecSBarry Smith */ 271003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4 272003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 273003131ecSBarry Smith if (nnz > 0) {\ 274003131ecSBarry Smith switch (nnz & 0x3) {\ 275003131ecSBarry Smith case 3: sum -= *xv++ * r[*xi++];\ 276003131ecSBarry Smith case 2: sum -= *xv++ * r[*xi++];\ 277003131ecSBarry Smith case 1: sum -= *xv++ * r[*xi++];\ 278003131ecSBarry Smith nnz -= 4;}\ 279003131ecSBarry Smith while (nnz > 0) {\ 280003131ecSBarry Smith sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\ 281003131ecSBarry Smith xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\ 282003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 283003131ecSBarry Smith 284003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 285003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 286003131ecSBarry Smith PetscInt __i,__i1,__i2;\ 287003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 288003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 289003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 290003131ecSBarry Smith 291003131ecSBarry Smith #else 292003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\ 293003131ecSBarry Smith PetscInt __i;\ 294003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];} 295003131ecSBarry Smith #endif 296003131ecSBarry Smith 297003131ecSBarry Smith 298003131ecSBarry Smith 299003131ecSBarry Smith /* 300003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 301003131ecSBarry Smith 302003131ecSBarry Smith Input Parameters: 303003131ecSBarry Smith + nnz - the number of entries 304003131ecSBarry Smith . r - the array of vector values 305003131ecSBarry Smith . xv - the matrix values for the row 306003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 307003131ecSBarry Smith 308003131ecSBarry Smith Output Parameter: 309003131ecSBarry Smith . sum - the sum of results 310003131ecSBarry Smith 311003131ecSBarry Smith PETSc compile flags: 312b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 313e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 314003131ecSBarry Smith 315003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot() 316003131ecSBarry Smith 317003131ecSBarry Smith */ 318003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4 319003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 320003131ecSBarry Smith if (nnz > 0) {\ 321003131ecSBarry Smith switch (nnz & 0x3) {\ 322003131ecSBarry Smith case 3: sum += *xv++ * r[*xi++];\ 323003131ecSBarry Smith case 2: sum += *xv++ * r[*xi++];\ 324003131ecSBarry Smith case 1: sum += *xv++ * r[*xi++];\ 325003131ecSBarry Smith nnz -= 4;}\ 326003131ecSBarry Smith while (nnz > 0) {\ 327003131ecSBarry Smith sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\ 328003131ecSBarry Smith xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\ 329003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 330003131ecSBarry Smith 331003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 332003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 333003131ecSBarry Smith PetscInt __i,__i1,__i2;\ 334003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\ 335003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\ 336003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 337003131ecSBarry Smith 338003131ecSBarry Smith #else 339003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\ 340003131ecSBarry Smith PetscInt __i;\ 341003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];} 342003131ecSBarry Smith #endif 343003131ecSBarry Smith 3442d40f771SBarry Smith #endif 345