1b8a66259SBarry Smith 22d40f771SBarry Smith #if !defined(__AIJ_H) 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5b45d2f2cSJed Brown #include <petsc-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 */ \ 337b083b7cSBarry Smith PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 34ace3abfcSBarry Smith PetscBool free_diag; \ 35421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 36e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 374fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 38ace3abfcSBarry Smith PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 394fd072dbSBarry Smith Mat parent /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); 404fd072dbSBarry Smith means that this shares some data structures with the parent including diag, ilen, imax, i, j */ 41e6b907acSBarry Smith 4253565b12SHong Zhang typedef struct { 4353565b12SHong Zhang MatTransposeColoring matcoloring; 4453565b12SHong Zhang Mat Bt_den; /* dense matrix of B^T */ 4553565b12SHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 4653565b12SHong Zhang PetscBool usecoloring; 4753565b12SHong Zhang PetscErrorCode (*destroy)(Mat); 4853565b12SHong Zhang } Mat_MatMatTransMult; 4953565b12SHong Zhang 502cff0574SHong Zhang typedef struct { /* for MatTransposeMatMult_SeqAIJ_SeqDense() */ 512cff0574SHong Zhang Mat mA; /* maij matrix of A */ 522cff0574SHong Zhang Vec bt,ct; /* vectors to hold locally transposed arrays of B and C */ 532cff0574SHong Zhang PetscErrorCode (*destroy)(Mat); 542cff0574SHong Zhang } Mat_MatTransMatMult; 552cff0574SHong Zhang 5653565b12SHong Zhang typedef struct { 5753565b12SHong Zhang PetscInt *api,*apj; /* symbolic structure of A*P */ 5853565b12SHong Zhang PetscScalar *apa; /* temporary array for storing one row of A*P */ 5953565b12SHong Zhang PetscErrorCode (*destroy)(Mat); 6053565b12SHong Zhang } Mat_PtAP; 6153565b12SHong Zhang 6253565b12SHong Zhang typedef struct { 6353565b12SHong Zhang MatTransposeColoring matcoloring; 64257c235dSHong Zhang Mat Rt; /* sparse or dense matrix of R^T */ 6553565b12SHong Zhang Mat RARt; /* dense matrix of R*A*R^T */ 663b1b9624SHong Zhang Mat ARt; /* A*R^T used for the case -matrart_color_art */ 6753565b12SHong Zhang MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 6853565b12SHong Zhang PetscErrorCode (*destroy)(Mat); 6953565b12SHong Zhang } Mat_RARt; 7053565b12SHong Zhang 716d0b6147SHong Zhang typedef struct { 726d0b6147SHong Zhang Mat BC; /* temp matrix for storing B*C */ 736d0b6147SHong Zhang PetscErrorCode (*destroy)(Mat); 746d0b6147SHong Zhang } Mat_MatMatMatMult; 756d0b6147SHong Zhang 762d40f771SBarry Smith /* 77ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 78e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 79dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 805768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 812d40f771SBarry Smith */ 82d35516d3SLois Curfman McInnes 83e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 84b8a66259SBarry Smith typedef struct { 854108e4d5SBarry Smith MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 86f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 87ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 88f0d39aaaSBarry Smith 89ace3abfcSBarry Smith PetscBool use; 90e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 91e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 92e6b907acSBarry Smith PetscInt limit; /* inode limit */ 93e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 94ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 954108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 96e6b907acSBarry Smith 975a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 985a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 995a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 1005a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool); 1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool); 1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 106e6b907acSBarry Smith 107e6b907acSBarry Smith typedef struct { 10854f21887SBarry Smith SEQAIJHEADER(MatScalar); 1094108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 11054f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 11171f1c65dSBarry Smith 11271f1c65dSBarry Smith PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 113ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 114bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 1152291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 11671f1c65dSBarry Smith PetscScalar fshift,omega; /* last used omega and fshift */ 11771f1c65dSBarry Smith 1183a7fca6bSBarry Smith ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 1190b7e3e3dSHong Zhang 1200b7e3e3dSHong Zhang PetscScalar *matmult_abdense; /* used by MatMatMult() */ 12153565b12SHong Zhang Mat_PtAP *ptap; /* used by MatPtAP() */ 1226d0b6147SHong Zhang Mat_MatMatMatMult *matmatmatmult; /* used by MatMatMatMult() */ 123257c235dSHong Zhang Mat_RARt *rart; /* used by MatRARt() */ 12440192850SHong Zhang Mat_MatMatTransMult *abt; /* used by MatMatTransposeMult() */ 125ec8511deSBarry Smith } Mat_SeqAIJ; 126b8a66259SBarry Smith 127e6b907acSBarry Smith /* 128e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 129e6b907acSBarry Smith */ 130e6b907acSBarry Smith #undef __FUNCT__ 131e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ" 132d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 133d0f46423SBarry Smith { 134e6b907acSBarry Smith PetscErrorCode ierr; 135e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 136e6b907acSBarry Smith if (A->singlemalloc) { 137e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 138e6b907acSBarry Smith } else { 139c31cb41cSBarry Smith if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 140c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 141c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 142e6b907acSBarry Smith } 143e6b907acSBarry Smith return 0; 144e6b907acSBarry Smith } 145e6b907acSBarry Smith /* 146e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 147357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 148e6b907acSBarry Smith */ 149fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 150fef13f97SBarry Smith if (NROW >= RMAX) { \ 151fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 152fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 153fef13f97SBarry Smith PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 154fef13f97SBarry Smith datatype *new_a; \ 155fef13f97SBarry Smith \ 1565aae435aSMatthew G. Knepley if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \ 157fef13f97SBarry Smith /* malloc new storage space */ \ 158dcca6d9dSJed Brown ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \ 159fef13f97SBarry Smith \ 160fef13f97SBarry Smith /* copy over old data into new slots */ \ 161fef13f97SBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 162fef13f97SBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 163fef13f97SBarry Smith ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 164fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 165fef13f97SBarry Smith ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 166fef13f97SBarry Smith ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 167fef13f97SBarry Smith ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \ 168fef13f97SBarry Smith ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 169fef13f97SBarry Smith /* free up old matrix storage */ \ 170fef13f97SBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 171fef13f97SBarry Smith AA = new_a; \ 172fef13f97SBarry Smith Ain->a = (MatScalar*) new_a; \ 173fef13f97SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 174fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 175fef13f97SBarry Smith \ 176fef13f97SBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 177fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 178fef13f97SBarry Smith Ain->maxnz += BS2*CHUNKSIZE; \ 179fef13f97SBarry Smith Ain->reallocs++; \ 180fef13f97SBarry Smith } \ 18117454e89SShri Abhyankar 182e6b907acSBarry Smith 1838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 1845a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 1855a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 1865a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 1871df811f5SHong Zhang 1885a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 1895a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 1905a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 1915a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 1925a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 1935a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 1945a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 1955a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 1965a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*); 1975a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 1985a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**); 19908480c60SBarry Smith 2005a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 2015a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 2025a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 2035a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 2045a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 20508480c60SBarry Smith 2065a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool); 2075a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 2085a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 2093a7fca6bSBarry Smith 2105a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 2115a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 2125a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 2135a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*); 2145008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*); 2155a576424SJed Brown PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 2165a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 2175a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 2185a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 2195a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 2205a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 2225a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 23779c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*); 23893dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring); 239f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring); 240a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt); 241d1e9a80fSBarry Smith PETSC_INTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*); 2425a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 2435a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 2447bab7c10SHong Zhang 2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*); 2485a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*); 2495a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*); 2505a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*); 2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 2525a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat); 2532b8ad9a3SHong Zhang 2545a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2554a1b09b7SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat,Mat,PetscReal,Mat*); 2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*); 2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat); 25953565b12SHong Zhang 2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 26155bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*); 26255bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*); 2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 26455bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat); 26555bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat); 2665df89d91SHong Zhang 2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 2703bf78175SHong Zhang 2713bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*); 2723bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*); 2733bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat); 2743bf78175SHong Zhang 2755a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 2785a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 2795a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 2805a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 2815df89d91SHong Zhang 2825a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*); 2845a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat); 2857bab7c10SHong Zhang 2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 2885a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 2895a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 2905a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 2915a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 2925a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 2935a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 2946378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 2956378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat); 2985a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 2999af31e4aSHong Zhang 3005a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 3015a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 3025a576424SJed Brown PETSC_INTERN PetscErrorCode Mat_CheckInode(Mat,PetscBool); 303abb87a52SBarry Smith PETSC_INTERN PetscErrorCode Mat_CheckInode_FactorLU(Mat); 304019b515eSShri Abhyankar 3055a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 3069f5f6813SShri Abhyankar 3078cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*); 3088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*); 3098cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*); 3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 3115a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3125a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 3145a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 3155a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 3162f947c57SVictor Minden 317b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 318de4b53ebSHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_Redundant(Mat_Redundant **); 319*9c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 320*9c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 32170f19b1fSKris Buschelman 322003131ecSBarry Smith /* 323003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 324003131ecSBarry Smith 325003131ecSBarry Smith Input Parameters: 326003131ecSBarry Smith + nnz - the number of entries 327003131ecSBarry Smith . r - the array of vector values 328003131ecSBarry Smith . xv - the matrix values for the row 329003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 330003131ecSBarry Smith 331003131ecSBarry Smith Output Parameter: 332003131ecSBarry Smith . sum - negative the sum of results 333003131ecSBarry Smith 334003131ecSBarry Smith PETSc compile flags: 335b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 336e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 337003131ecSBarry Smith 338003131ecSBarry Smith .seealso: PetscSparseDensePlusDot() 339003131ecSBarry Smith 340003131ecSBarry Smith */ 341519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 342003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 343003131ecSBarry Smith if (nnz > 0) { \ 344003131ecSBarry Smith switch (nnz & 0x3) { \ 345003131ecSBarry Smith case 3: sum -= *xv++ *r[*xi++]; \ 346003131ecSBarry Smith case 2: sum -= *xv++ *r[*xi++]; \ 347003131ecSBarry Smith case 1: sum -= *xv++ *r[*xi++]; \ 348003131ecSBarry Smith nnz -= 4;} \ 349003131ecSBarry Smith while (nnz > 0) { \ 350003131ecSBarry Smith sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \ 351003131ecSBarry Smith xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \ 352003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 353003131ecSBarry Smith 354003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 355003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 356003131ecSBarry Smith PetscInt __i,__i1,__i2; \ 357003131ecSBarry Smith for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 358003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 359003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 360003131ecSBarry Smith 361003131ecSBarry Smith #else 362003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 363003131ecSBarry Smith PetscInt __i; \ 364003131ecSBarry Smith for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];} 365003131ecSBarry Smith #endif 366003131ecSBarry Smith 367003131ecSBarry Smith 368003131ecSBarry Smith 369003131ecSBarry Smith /* 370003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 371003131ecSBarry Smith 372003131ecSBarry Smith Input Parameters: 373003131ecSBarry Smith + nnz - the number of entries 374003131ecSBarry Smith . r - the array of vector values 375003131ecSBarry Smith . xv - the matrix values for the row 376003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 377003131ecSBarry Smith 378003131ecSBarry Smith Output Parameter: 379003131ecSBarry Smith . sum - the sum of results 380003131ecSBarry Smith 381003131ecSBarry Smith PETSc compile flags: 382b2843cf1SBarry Smith + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 383e0e43300SBarry Smith - PETSC_KERNEL_USE_UNROLL_2 - 384003131ecSBarry Smith 385003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot() 386003131ecSBarry Smith 387003131ecSBarry Smith */ 388519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 389003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 390003131ecSBarry Smith if (nnz > 0) { \ 391003131ecSBarry Smith switch (nnz & 0x3) { \ 392003131ecSBarry Smith case 3: sum += *xv++ *r[*xi++]; \ 393003131ecSBarry Smith case 2: sum += *xv++ *r[*xi++]; \ 394003131ecSBarry Smith case 1: sum += *xv++ *r[*xi++]; \ 395003131ecSBarry Smith nnz -= 4;} \ 396003131ecSBarry Smith while (nnz > 0) { \ 397003131ecSBarry Smith sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 398003131ecSBarry Smith xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 399003131ecSBarry Smith xv += 4; xi += 4; nnz -= 4; }}} 400003131ecSBarry Smith 401003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 402003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 403003131ecSBarry Smith PetscInt __i,__i1,__i2; \ 404003131ecSBarry Smith for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 405003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 406003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 407003131ecSBarry Smith 408003131ecSBarry Smith #else 409003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 410003131ecSBarry Smith PetscInt __i; \ 411003131ecSBarry Smith for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];} 412003131ecSBarry Smith #endif 413003131ecSBarry Smith 414b434eb95SMatthew G. Knepley 415b434eb95SMatthew G. Knepley /* 416b434eb95SMatthew G. Knepley PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 417b434eb95SMatthew G. Knepley 418b434eb95SMatthew G. Knepley Input Parameters: 419b434eb95SMatthew G. Knepley + nnz - the number of entries 420b434eb95SMatthew G. Knepley . r - the array of vector values 421b434eb95SMatthew G. Knepley . xv - the matrix values for the row 422b434eb95SMatthew G. Knepley - xi - the column indices of the nonzeros in the row 423b434eb95SMatthew G. Knepley 424b434eb95SMatthew G. Knepley Output Parameter: 425b434eb95SMatthew G. Knepley . max - the max of results 426b434eb95SMatthew G. Knepley 427b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot() 428b434eb95SMatthew G. Knepley 429b434eb95SMatthew G. Knepley */ 430b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \ 431b434eb95SMatthew G. Knepley PetscInt __i; \ 4325e792707SMatthew G. Knepley for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));} 433b434eb95SMatthew G. Knepley 4342d40f771SBarry Smith #endif 435