1b8a66259SBarry Smith 22d40f771SBarry Smith #if !defined(__AIJ_H) 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5*7c4f633dSBarry Smith #include "private/matimpl.h" 6e6b907acSBarry Smith /* 7e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 8e6b907acSBarry Smith */ 9421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 10e6b907acSBarry Smith PetscTruth roworiented; /* if true, row-oriented input, default */\ 11e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */\ 12e6b907acSBarry Smith PetscTruth singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\ 13e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */\ 14e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */\ 15e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */\ 16e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 17e6b907acSBarry Smith as more values are set than were prealloced */\ 18e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */\ 19e6b907acSBarry Smith PetscTruth keepzeroedrows; /* keeps matrix structure same in calls to MatZeroRows()*/\ 20e6b907acSBarry Smith PetscTruth ignorezeroentries; \ 21e6b907acSBarry Smith PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\ 22e6b907acSBarry Smith Mat XtoY; /* used by MatAXPY() */\ 23e6b907acSBarry Smith PetscTruth free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 24e6b907acSBarry Smith PetscTruth free_a; /* free the numerical values when matrix is destroy */ \ 25e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 26e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 27e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 28e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 29e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 30421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 31e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 32e6b907acSBarry Smith IS row, col, icol /* index sets, used for reorderings */ 33e6b907acSBarry Smith 342d40f771SBarry Smith /* 35ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 36e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 37dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 385768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 392d40f771SBarry Smith */ 40d35516d3SLois Curfman McInnes 41e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 42b8a66259SBarry Smith typedef struct { 43a77337e4SBarry Smith MatScalar *bdiag,*ibdiag; /* diagonal blocks of matrix used for MatRelax_Inode() */ 44f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 4571f1c65dSBarry Smith PetscTruth ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 46f0d39aaaSBarry Smith 47e6b907acSBarry Smith PetscTruth use; 48e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 49e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 50e6b907acSBarry Smith PetscInt limit; /* inode limit */ 51e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 52e6b907acSBarry Smith PetscTruth checked; /* if inodes have been checked for */ 53e6b907acSBarry Smith } Mat_Inode; 54e6b907acSBarry Smith 55e6b907acSBarry Smith EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer); 56e6b907acSBarry Smith EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType); 57e6b907acSBarry Smith EXTERN PetscErrorCode MatDestroy_Inode(Mat); 58e6b907acSBarry Smith EXTERN PetscErrorCode MatCreate_Inode(Mat); 594e0d8c25SBarry Smith EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth); 60719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicate_Inode(Mat,MatDuplicateOption,Mat*); 610481f469SBarry Smith EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat,IS,IS,const MatFactorInfo*,Mat*); 620481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 630481f469SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*); 64e6b907acSBarry Smith 65e6b907acSBarry Smith 66e6b907acSBarry Smith typedef struct { 6754f21887SBarry Smith SEQAIJHEADER(MatScalar); 68e6b907acSBarry Smith Mat_Inode inode; 6954f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 7071f1c65dSBarry Smith 7171f1c65dSBarry Smith PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 7271f1c65dSBarry Smith PetscTruth idiagvalid; /* current idiag[] and mdiag[] are valid */ 7371f1c65dSBarry Smith PetscScalar fshift,omega; /* last used omega and fshift */ 7471f1c65dSBarry Smith 753a7fca6bSBarry Smith ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 76ec8511deSBarry Smith } Mat_SeqAIJ; 77b8a66259SBarry Smith 78e6b907acSBarry Smith /* 79e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 80e6b907acSBarry Smith */ 81e6b907acSBarry Smith #undef __FUNCT__ 82e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ" 83d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 84d0f46423SBarry Smith { 85e6b907acSBarry Smith PetscErrorCode ierr; 86e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 87e6b907acSBarry Smith if (A->singlemalloc) { 88e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 89e6b907acSBarry Smith } else { 90e6b907acSBarry Smith if (A->free_a && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 91e6b907acSBarry Smith if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);} 92e6b907acSBarry Smith if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);} 93e6b907acSBarry Smith } 94e6b907acSBarry Smith *a = 0; *j = 0; *i = 0; 95e6b907acSBarry Smith return 0; 96e6b907acSBarry Smith } 97e6b907acSBarry Smith 98e6b907acSBarry Smith /* 99e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 100e6b907acSBarry Smith */ 101421e10b8SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 102e6b907acSBarry Smith if (NROW >= RMAX) {\ 103e6b907acSBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\ 104e6b907acSBarry Smith /* there is no extra room in row, therefore enlarge */ \ 105e6b907acSBarry Smith PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 106421e10b8SBarry Smith datatype *new_a; \ 107e6b907acSBarry Smith \ 108e6b907acSBarry Smith if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \ 109e6b907acSBarry Smith /* malloc new storage space */ \ 110421e10b8SBarry Smith ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\ 111e6b907acSBarry Smith \ 112e6b907acSBarry Smith /* copy over old data into new slots */ \ 113e6b907acSBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 114e6b907acSBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 115e6b907acSBarry Smith ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 116e6b907acSBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 117e6b907acSBarry Smith ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 118421e10b8SBarry Smith ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 119421e10b8SBarry Smith ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\ 120421e10b8SBarry Smith ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 121e6b907acSBarry Smith /* free up old matrix storage */ \ 122e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\ 123421e10b8SBarry Smith AA = new_a; \ 124421e10b8SBarry Smith Ain->a = (MatScalar*) new_a; \ 125421e10b8SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 126e6b907acSBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 127e6b907acSBarry Smith \ 128e6b907acSBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 129e6b907acSBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 130e6b907acSBarry Smith Ain->maxnz += CHUNKSIZE; \ 131e6b907acSBarry Smith Ain->reallocs++; \ 132e6b907acSBarry Smith } \ 133e6b907acSBarry Smith 134f349cc81SSatish Balay EXTERN_C_BEGIN 135ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*); 136f349cc81SSatish Balay EXTERN_C_END 1370481f469SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 1380481f469SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 1390481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 1400481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 141dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 14209f38230SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*); 143dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 14408480c60SBarry Smith 145dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 146dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 147dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 148dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 14938baddfdSBarry Smith EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 15008480c60SBarry Smith 151dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring); 152dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*); 15338baddfdSBarry Smith EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*); 1543a7fca6bSBarry Smith 15538baddfdSBarry Smith EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 1560390132cSHong Zhang EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 15738baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 15838baddfdSBarry Smith EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**); 1590481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 1600481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 1610481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 1620481f469SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 163dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 1648d8c7c43SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 165137fb511SHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 166dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 167dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 168dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 1692bebee5dSHong Zhang EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 170dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg); 171dfbe8321SBarry Smith EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1720481f469SBarry Smith EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*,Mat*); 173a313700dSBarry Smith EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*); 174dfbe8321SBarry Smith EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 175dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 176dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 177dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 1787ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*); 1797ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat); 180dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 181dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 182bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 183bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 184bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 18538baddfdSBarry Smith EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 186a77337e4SBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 187a77337e4SBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 188f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 1898f7157efSSatish Balay EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 1908f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 1918f7157efSSatish Balay EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 1928f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 193b24902e0SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 194b24902e0SBarry Smith EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 1959af31e4aSHong Zhang 19697304618SKris Buschelman EXTERN_C_BEGIN 1978cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*); 1988cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 1998cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*); 200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 20197304618SKris Buschelman EXTERN_C_END 20270f19b1fSKris Buschelman 2032d40f771SBarry Smith #endif 204