1b8a66259SBarry Smith 22d40f771SBarry Smith #if !defined(__AIJ_H) 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5e6b907acSBarry Smith #include "src/mat/matimpl.h" 6e6b907acSBarry Smith /* 7e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 8e6b907acSBarry Smith */ 9e6b907acSBarry Smith #define SEQAIJHEADER \ 10e6b907acSBarry Smith PetscTruth sorted; /* if true, rows are sorted by increasing columns */\ 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 */\ 13e6b907acSBarry Smith PetscTruth singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\ 14e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */\ 15e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */\ 16e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */\ 17e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 18e6b907acSBarry Smith as more values are set than were prealloced */\ 19e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */\ 20e6b907acSBarry Smith PetscTruth keepzeroedrows; /* keeps matrix structure same in calls to MatZeroRows()*/\ 21e6b907acSBarry Smith PetscTruth ignorezeroentries; \ 22e6b907acSBarry Smith PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 23e6b907acSBarry Smith Mat XtoY; /* used by MatAXPY() */\ 24e6b907acSBarry Smith PetscTruth free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 25e6b907acSBarry Smith PetscTruth free_a; /* free the numerical values when matrix is destroy */ \ 26e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 27e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 28e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 29e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 30e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 31e6b907acSBarry Smith PetscScalar *a; /* nonzero elements */ \ 32e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 33e6b907acSBarry Smith IS row, col, icol /* index sets, used for reorderings */ 34e6b907acSBarry Smith 352d40f771SBarry Smith /* 36ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 37e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 38dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 395768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 402d40f771SBarry Smith */ 41d35516d3SLois Curfman McInnes 42e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 43b8a66259SBarry Smith typedef struct { 44e6b907acSBarry Smith PetscTruth use; 45e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 46e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 47e6b907acSBarry Smith PetscInt limit; /* inode limit */ 48e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 49e6b907acSBarry Smith PetscTruth checked; /* if inodes have been checked for */ 50e6b907acSBarry Smith } Mat_Inode; 51e6b907acSBarry Smith 52e6b907acSBarry Smith EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer); 53e6b907acSBarry Smith EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType); 54e6b907acSBarry Smith EXTERN PetscErrorCode MatDestroy_Inode(Mat); 55e6b907acSBarry Smith EXTERN PetscErrorCode MatCreate_Inode(Mat); 56e6b907acSBarry Smith EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption); 57e6b907acSBarry Smith EXTERN PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *B); 58e6b907acSBarry Smith EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact); 59e6b907acSBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact); 60e6b907acSBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact); 61e6b907acSBarry Smith 62e6b907acSBarry Smith 63e6b907acSBarry Smith typedef struct { 64e6b907acSBarry Smith SEQAIJHEADER; 65e6b907acSBarry Smith Mat_Inode inode; 66ea709b57SSatish Balay PetscScalar *saved_values; /* location for stashing nonzero values of matrix */ 67ea709b57SSatish Balay PetscScalar *idiag,*ssor; /* inverse of diagonal entries; space for eisen */ 683a7fca6bSBarry Smith ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 69ec8511deSBarry Smith } Mat_SeqAIJ; 70b8a66259SBarry Smith 71e6b907acSBarry Smith /* 72e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 73e6b907acSBarry Smith */ 74e6b907acSBarry Smith #undef __FUNCT__ 75e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ" 76e6b907acSBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,PetscScalar **a,PetscInt **j,PetscInt **i) { 77e6b907acSBarry Smith PetscErrorCode ierr; 78e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 79e6b907acSBarry Smith if (A->singlemalloc) { 80e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 81e6b907acSBarry Smith } else { 82e6b907acSBarry Smith if (A->free_a && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 83e6b907acSBarry Smith if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);} 84e6b907acSBarry Smith if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);} 85e6b907acSBarry Smith } 86e6b907acSBarry Smith *a = 0; *j = 0; *i = 0; 87e6b907acSBarry Smith return 0; 88e6b907acSBarry Smith } 89e6b907acSBarry Smith 90e6b907acSBarry Smith /* 91e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 92e6b907acSBarry Smith */ 93e6b907acSBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW) \ 94e6b907acSBarry Smith if (NROW >= RMAX) { \ 95e6b907acSBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 96e6b907acSBarry Smith /* there is no extra room in row, therefore enlarge */ \ 97e6b907acSBarry Smith PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 98e6b907acSBarry Smith PetscScalar *new_a; \ 99e6b907acSBarry Smith \ 100e6b907acSBarry Smith if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 101e6b907acSBarry Smith /* malloc new storage space */ \ 102e6b907acSBarry Smith ierr = PetscMalloc3(BS2*new_nz,PetscScalar,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 103e6b907acSBarry Smith \ 104e6b907acSBarry Smith /* copy over old data into new slots */ \ 105e6b907acSBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 106e6b907acSBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 107e6b907acSBarry Smith ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 108e6b907acSBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 109e6b907acSBarry Smith ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 110e6b907acSBarry Smith ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(PetscScalar));CHKERRQ(ierr); \ 111e6b907acSBarry Smith ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);\ 112e6b907acSBarry Smith ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(PetscScalar));CHKERRQ(ierr); \ 113e6b907acSBarry Smith /* free up old matrix storage */ \ 114e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 115e6b907acSBarry Smith AA = Ain->a = new_a; AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 116e6b907acSBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 117e6b907acSBarry Smith \ 118e6b907acSBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 119e6b907acSBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 120e6b907acSBarry Smith Ain->maxnz += CHUNKSIZE; \ 121e6b907acSBarry Smith Ain->reallocs++; \ 122e6b907acSBarry Smith } \ 123e6b907acSBarry Smith 124f349cc81SSatish Balay EXTERN_C_BEGIN 125ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*); 126f349cc81SSatish Balay EXTERN_C_END 127dfbe8321SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat *); 128dfbe8321SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat *); 129dfbe8321SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat*); 130af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*); 131dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 13209f38230SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*); 133dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 13408480c60SBarry Smith 135dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 136dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 137dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 138dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 13938baddfdSBarry Smith EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 14008480c60SBarry Smith 141dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 142dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 14338baddfdSBarry Smith EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 1443a7fca6bSBarry Smith 14538baddfdSBarry Smith EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 1460390132cSHong Zhang EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 14738baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 14838baddfdSBarry Smith EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 149dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*); 150af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*); 151*137fb511SHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,MatFactorInfo*,Mat*); 152dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*); 153dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 154*137fb511SHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 155dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 156dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 157dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 1582bebee5dSHong Zhang EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 159dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 160dfbe8321SBarry Smith EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1617529efd4SKris Buschelman EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*); 162f69a0ea3SMatthew Knepley EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, MatType,Mat*); 163dfbe8321SBarry Smith EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 164dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 165dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 166dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 1677ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 1687ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 169dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 170dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 171bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 172bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 173bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 17438baddfdSBarry Smith EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 17538baddfdSBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 17638baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 177f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 1789af31e4aSHong Zhang 17997304618SKris Buschelman EXTERN_C_BEGIN 180f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType,MatReuse,Mat*); 181f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 182ba03775dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,MatType,MatReuse,Mat*); 183be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 18497304618SKris Buschelman EXTERN_C_END 18570f19b1fSKris Buschelman 1862d40f771SBarry Smith #endif 187