1b8a66259SBarry Smith 26524c165SJacob Faibussowitsch #ifndef __AIJ_H 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5af0996ceSBarry Smith #include <petsc/private/matimpl.h> 6eca6b86aSHong Zhang #include <petscctable.h> 78f690400SShri Abhyankar 8d67d9f35SJunchao Zhang /* Operations provided by MATSEQAIJ and its subclasses */ 9d67d9f35SJunchao Zhang typedef struct { 10d67d9f35SJunchao Zhang PetscErrorCode (*getarray)(Mat, PetscScalar **); 11d67d9f35SJunchao Zhang PetscErrorCode (*restorearray)(Mat, PetscScalar **); 12d67d9f35SJunchao Zhang PetscErrorCode (*getarrayread)(Mat, const PetscScalar **); 13d67d9f35SJunchao Zhang PetscErrorCode (*restorearrayread)(Mat, const PetscScalar **); 14d67d9f35SJunchao Zhang PetscErrorCode (*getarraywrite)(Mat, PetscScalar **); 15d67d9f35SJunchao Zhang PetscErrorCode (*restorearraywrite)(Mat, PetscScalar **); 167ee59b9bSJunchao Zhang PetscErrorCode (*getcsrandmemtype)(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *); 17d67d9f35SJunchao Zhang } Mat_SeqAIJOps; 18d67d9f35SJunchao Zhang 19e6b907acSBarry Smith /* 20e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 21e6b907acSBarry Smith */ 22421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 23ace3abfcSBarry Smith PetscBool roworiented; /* if true, row-oriented input, default */ \ 24e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 2528b2fa4aSMatthew Knepley PetscInt nounused; /* -1 generate error on unused space */ \ 26ace3abfcSBarry Smith PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \ 27e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */ \ 28e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */ \ 29e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */ \ 30846b4da1SFande Kong PetscInt *ipre; /* space preallocated for each row by user */ \ 31ace3abfcSBarry Smith PetscBool free_imax_ilen; \ 32e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 33e6b907acSBarry Smith as more values are set than were prealloced */ \ 34e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */ \ 35ace3abfcSBarry Smith PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \ 36ace3abfcSBarry Smith PetscBool ignorezeroentries; \ 37ace3abfcSBarry Smith PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 38ace3abfcSBarry Smith PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 39e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 40e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 41e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 42e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 43e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 447b083b7cSBarry Smith PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 45ace3abfcSBarry Smith PetscBool free_diag; \ 46421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 47e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 484fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 49ace3abfcSBarry Smith PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 5017df9f7cSHong Zhang Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \ 5117df9f7cSHong Zhang means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \ 52d67d9f35SJunchao Zhang Mat_SubSppt *submatis1; /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \ 53d67d9f35SJunchao Zhang Mat_SeqAIJOps ops[1] /* operations for SeqAIJ and its subclasses */ 54e6b907acSBarry Smith 5553565b12SHong Zhang typedef struct { 5653565b12SHong Zhang MatTransposeColoring matcoloring; 5753565b12SHong Zhang Mat Bt_den; /* dense matrix of B^T */ 5853565b12SHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 5953565b12SHong Zhang PetscBool usecoloring; 6053565b12SHong Zhang } Mat_MatMatTransMult; 6153565b12SHong Zhang 626d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */ 636d373c3eSHong Zhang Mat At; /* transpose of the first matrix */ 642cff0574SHong Zhang Mat mA; /* maij matrix of A */ 652cff0574SHong Zhang Vec bt, ct; /* vectors to hold locally transposed arrays of B and C */ 666718818eSStefano Zampini /* used by PtAP */ 676718818eSStefano Zampini void *data; 686718818eSStefano Zampini PetscErrorCode (*destroy)(void *); 692cff0574SHong Zhang } Mat_MatTransMatMult; 702cff0574SHong Zhang 7153565b12SHong Zhang typedef struct { 7253565b12SHong Zhang PetscInt *api, *apj; /* symbolic structure of A*P */ 7353565b12SHong Zhang PetscScalar *apa; /* temporary array for storing one row of A*P */ 743cdca5ebSHong Zhang } Mat_AP; 7553565b12SHong Zhang 7653565b12SHong Zhang typedef struct { 7753565b12SHong Zhang MatTransposeColoring matcoloring; 78257c235dSHong Zhang Mat Rt; /* sparse or dense matrix of R^T */ 7953565b12SHong Zhang Mat RARt; /* dense matrix of R*A*R^T */ 803b1b9624SHong Zhang Mat ARt; /* A*R^T used for the case -matrart_color_art */ 8153565b12SHong Zhang MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 826718818eSStefano Zampini /* free intermediate products needed for PtAP */ 836718818eSStefano Zampini void *data; 846718818eSStefano Zampini PetscErrorCode (*destroy)(void *); 8553565b12SHong Zhang } Mat_RARt; 8653565b12SHong Zhang 876d0b6147SHong Zhang typedef struct { 886d0b6147SHong Zhang Mat BC; /* temp matrix for storing B*C */ 896d0b6147SHong Zhang } Mat_MatMatMatMult; 906d0b6147SHong Zhang 912d40f771SBarry Smith /* 92ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 93e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 94dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 955768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 962d40f771SBarry Smith */ 97d35516d3SLois Curfman McInnes 98e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 99b8a66259SBarry Smith typedef struct { 1004108e4d5SBarry Smith MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 101f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 102ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 103f0d39aaaSBarry Smith 104ace3abfcSBarry Smith PetscBool use; 105e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 106e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 107e6b907acSBarry Smith PetscInt limit; /* inode limit */ 108e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 109ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 110a02bda8eSBarry Smith PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 1114108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 112e6b907acSBarry Smith 1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer); 1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType); 1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 1175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool); 1185a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *); 1195a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool); 1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat, Mat, const MatFactorInfo *); 1215a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *); 122f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **); 123f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **); 124e6b907acSBarry Smith 125e6b907acSBarry Smith typedef struct { 12654f21887SBarry Smith SEQAIJHEADER(MatScalar); 1274108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 12854f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 12971f1c65dSBarry Smith 13071f1c65dSBarry Smith PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 131ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 132bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 1332291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 13461ecd0c6SBarry Smith PetscBool diagonaldense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */ 13571f1c65dSBarry Smith PetscScalar fshift, omega; /* last used omega and fshift */ 136394ed5ebSJunchao Zhang 137394ed5ebSJunchao Zhang /* MatSetValuesCOO() related fields on host */ 138394ed5ebSJunchao Zhang PetscCount coo_n; /* Number of entries in MatSetPreallocationCOO() */ 139394ed5ebSJunchao Zhang PetscCount Atot; /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */ 140394ed5ebSJunchao Zhang PetscCount *jmap; /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */ 141394ed5ebSJunchao Zhang PetscCount *perm; /* The permutation array in sorting (i,j) by row and then by col */ 142ec8511deSBarry Smith } Mat_SeqAIJ; 143b8a66259SBarry Smith 144e6b907acSBarry Smith /* 145e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 146e6b907acSBarry Smith */ 147d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i) 148d71ae5a4SJacob Faibussowitsch { 149e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data; 150e6b907acSBarry Smith if (A->singlemalloc) { 1519566063dSJacob Faibussowitsch PetscCall(PetscFree3(*a, *j, *i)); 152e6b907acSBarry Smith } else { 1539566063dSJacob Faibussowitsch if (A->free_a) PetscCall(PetscFree(*a)); 1549566063dSJacob Faibussowitsch if (A->free_ij) PetscCall(PetscFree(*j)); 1559566063dSJacob Faibussowitsch if (A->free_ij) PetscCall(PetscFree(*i)); 156e6b907acSBarry Smith } 157e6b907acSBarry Smith return 0; 158e6b907acSBarry Smith } 159e6b907acSBarry Smith /* 160e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 161357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 162e6b907acSBarry Smith */ 163fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \ 164fef13f97SBarry Smith if (NROW >= RMAX) { \ 165fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 166fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 167f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 168fef13f97SBarry Smith datatype *new_a; \ 169fef13f97SBarry Smith \ 17008401ef6SPierre Jolivet PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \ 171fef13f97SBarry Smith /* malloc new storage space */ \ 1729566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(BS2 *new_nz, &new_a, new_nz, &new_j, AM + 1, &new_i)); \ 173fef13f97SBarry Smith \ 174fef13f97SBarry Smith /* copy over old data into new slots */ \ 175ad540459SPierre Jolivet for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 176ad540459SPierre Jolivet for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 1779566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 178fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 1799566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \ 1809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \ 1819566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \ 1829566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), AA + BS2 * (AI[ROW] + NROW), BS2 * len)); \ 183fef13f97SBarry Smith /* free up old matrix storage */ \ 1849566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 185fef13f97SBarry Smith AA = new_a; \ 186fef13f97SBarry Smith Ain->a = (MatScalar *)new_a; \ 1879371c9d4SSatish Balay AI = Ain->i = new_i; \ 1889371c9d4SSatish Balay AJ = Ain->j = new_j; \ 189fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 190fef13f97SBarry Smith \ 1919371c9d4SSatish Balay RP = AJ + AI[ROW]; \ 1929371c9d4SSatish Balay AP = AA + BS2 * AI[ROW]; \ 193fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 194fef13f97SBarry Smith Ain->maxnz += BS2 * CHUNKSIZE; \ 195fef13f97SBarry Smith Ain->reallocs++; \ 1969371c9d4SSatish Balay } 19717454e89SShri Abhyankar 198876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \ 199720833daSHong Zhang if (NROW >= RMAX) { \ 200720833daSHong Zhang Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 201720833daSHong Zhang /* there is no extra room in row, therefore enlarge */ \ 202f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 203720833daSHong Zhang \ 20408401ef6SPierre Jolivet PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \ 205720833daSHong Zhang /* malloc new storage space */ \ 2069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_nz, &new_j)); \ 2079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(AM + 1, &new_i)); \ 208720833daSHong Zhang \ 209720833daSHong Zhang /* copy over old data into new slots */ \ 210ad540459SPierre Jolivet for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 211ad540459SPierre Jolivet for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 2129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 213720833daSHong Zhang len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 2149566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \ 215876c6284SHong Zhang \ 216720833daSHong Zhang /* free up old matrix storage */ \ 2179566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 218876c6284SHong Zhang Ain->a = NULL; \ 2199371c9d4SSatish Balay AI = Ain->i = new_i; \ 2209371c9d4SSatish Balay AJ = Ain->j = new_j; \ 221720833daSHong Zhang Ain->singlemalloc = PETSC_FALSE; \ 222876c6284SHong Zhang Ain->free_a = PETSC_FALSE; \ 223720833daSHong Zhang \ 224876c6284SHong Zhang RP = AJ + AI[ROW]; \ 225720833daSHong Zhang RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 226720833daSHong Zhang Ain->maxnz += BS2 * CHUNKSIZE; \ 227720833daSHong Zhang Ain->reallocs++; \ 2289371c9d4SSatish Balay } 229e6b907acSBarry Smith 230cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *); 231e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]); 232cbc6b225SStefano Zampini PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_SeqAIJ(Mat); 233cbc6b225SStefano Zampini 2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *); 2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *); 2371df811f5SHong Zhang 2385a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *); 2395a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 2405a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *); 2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 2425a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 2435a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 2445a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *); 2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure); 2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *); 2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2485a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **); 24908480c60SBarry Smith 250b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec); 251b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec); 252b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec); 253b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec); 254b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec); 255b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 257b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 25808480c60SBarry Smith 2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool); 2603a7fca6bSBarry Smith 2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 2637fb60732SBarry Smith PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]); 2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *); 2655008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *); 2667fb60732SBarry Smith 2672462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **); 2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *); 2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 2705a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 2715a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 2725a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *); 2735a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *); 2745a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec); 2755a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec); 2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat, Vec, Vec); 2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec); 2785a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat, Vec, Vec); 2795a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec); 2805a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat, Vec, Vec); 2815a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec); 2825a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec); 2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec); 2845a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec); 2855a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec); 2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat, Mat, Mat); 2885a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat); 28979c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *); 29093dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring); 291f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring); 292a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt); 29352f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer); 29452f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer); 2955a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer); 2965a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 2977bab7c10SHong Zhang 298df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE) 2994222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat); 300df97dc6dSFande Kong #endif 3014222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat); 3024222ddf1SHong Zhang 3034222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat); 3044222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat); 3054222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat); 3064222ddf1SHong Zhang 3074222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3084222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat); 3094222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat); 3104222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat); 3114222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat); 3124222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat); 3134222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat); 3144222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat); 3154222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat); 3164222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE) 3174222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat); 3184222ddf1SHong Zhang #endif 3194222ddf1SHong Zhang 3205a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 321df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat); 3224222ddf1SHong Zhang 3234099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat); 3245a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat); 3252b8ad9a3SHong Zhang 3264222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat); 3275a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3285a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat); 32953565b12SHong Zhang 3304222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3314222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat); 3324222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat); 3335a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 33455bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat); 33555bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat); 3365df89d91SHong Zhang 3374222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3385a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3396718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *); 3403bf78175SHong Zhang 3414222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3425a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3435a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring); 3445a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat); 3455a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat); 3465df89d91SHong Zhang 3474222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat); 3485a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat); 3497bab7c10SHong Zhang 350679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom); 3515a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 3525a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 354db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar); 35587c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec); 35687c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode); 3575a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure); 3585a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3595a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3605a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3615a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3626378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 3636378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 3645a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 3655a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat); 3665a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer); 3679af31e4aSHong Zhang 3685a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 3695a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 370a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 371a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 372019b515eSShri Abhyankar 3735a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *); 3749f5f6813SShri Abhyankar 375d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 376388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *); 377388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *); 378f2fbf96bSVaclav Hapla #endif 379cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *); 380cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *); 381388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *); 382388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *); 383388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *); 384d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK) 385d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *); 386d24d4204SJose E. Roman #endif 387388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *); 388cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *); 3894dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *); 3904a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *); 391388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *); 3925a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS); 3935a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *); 3948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 3955a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType); 3963fa6b06aSMark Adams PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat); 3972f947c57SVictor Minden 398b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *); 3999c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 4009c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 40170f19b1fSKris Buschelman 40253dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat); 403f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *); 4040fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat); 405f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat); 40663a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]); 407feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *); 40853dd7562SDmitry Karpeev 409a3bb6f32SFande Kong PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *); 410e4e71118SStefano Zampini PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat); 411a3bb6f32SFande Kong 412003131ecSBarry Smith /* 413003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 414003131ecSBarry Smith 415003131ecSBarry Smith Input Parameters: 416003131ecSBarry Smith + nnz - the number of entries 417003131ecSBarry Smith . r - the array of vector values 418003131ecSBarry Smith . xv - the matrix values for the row 419003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 420003131ecSBarry Smith 421003131ecSBarry Smith Output Parameter: 422003131ecSBarry Smith . sum - negative the sum of results 423003131ecSBarry Smith 424003131ecSBarry Smith PETSc compile flags: 4257b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4267b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4277b42bb93SJunchao Zhang 42811a5261eSBarry Smith Developer Note: 4297b42bb93SJunchao Zhang The macro changes sum but not other parameters 430003131ecSBarry Smith 431db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()` 432003131ecSBarry Smith */ 433519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 4349371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 4359371c9d4SSatish Balay { \ 436003131ecSBarry Smith if (nnz > 0) { \ 4377b42bb93SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 4387b42bb93SJunchao Zhang switch (rem) { \ 439d71ae5a4SJacob Faibussowitsch case 3: \ 440d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 441d71ae5a4SJacob Faibussowitsch case 2: \ 442d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 443d71ae5a4SJacob Faibussowitsch case 1: \ 444d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 445d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 4467b42bb93SJunchao Zhang } \ 4479371c9d4SSatish Balay while (nnz2 > 0) { \ 4489371c9d4SSatish Balay sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 4499371c9d4SSatish Balay xv += 4; \ 4509371c9d4SSatish Balay xi += 4; \ 4519371c9d4SSatish Balay nnz2 -= 4; \ 4529371c9d4SSatish Balay } \ 4539371c9d4SSatish Balay xv -= nnz; \ 4549371c9d4SSatish Balay xi -= nnz; \ 4557b42bb93SJunchao Zhang } \ 4567b42bb93SJunchao Zhang } 457003131ecSBarry Smith 458003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 4599371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 4609371c9d4SSatish Balay { \ 461003131ecSBarry Smith PetscInt __i, __i1, __i2; \ 4629371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 4639371c9d4SSatish Balay __i1 = xi[__i]; \ 4649371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 4659371c9d4SSatish Balay sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 4669371c9d4SSatish Balay } \ 4679371c9d4SSatish Balay if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \ 4689371c9d4SSatish Balay } 469003131ecSBarry Smith 470003131ecSBarry Smith #else 4719371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 4729371c9d4SSatish Balay { \ 473003131ecSBarry Smith PetscInt __i; \ 4749371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \ 4759371c9d4SSatish Balay } 476003131ecSBarry Smith #endif 477003131ecSBarry Smith 478003131ecSBarry Smith /* 479003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 480003131ecSBarry Smith 481003131ecSBarry Smith Input Parameters: 482003131ecSBarry Smith + nnz - the number of entries 483003131ecSBarry Smith . r - the array of vector values 484003131ecSBarry Smith . xv - the matrix values for the row 485003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 486003131ecSBarry Smith 487003131ecSBarry Smith Output Parameter: 488003131ecSBarry Smith . sum - the sum of results 489003131ecSBarry Smith 490003131ecSBarry Smith PETSc compile flags: 4917b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4927b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4937b42bb93SJunchao Zhang 49411a5261eSBarry Smith Developer Note: 4957b42bb93SJunchao Zhang The macro changes sum but not other parameters 496003131ecSBarry Smith 497db781477SPatrick Sanan .seealso: `PetscSparseDenseMinusDot()` 498003131ecSBarry Smith */ 499519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 5009371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 5019371c9d4SSatish Balay { \ 502003131ecSBarry Smith if (nnz > 0) { \ 5037b42bb93SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 5047b42bb93SJunchao Zhang switch (rem) { \ 505d71ae5a4SJacob Faibussowitsch case 3: \ 506d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 507d71ae5a4SJacob Faibussowitsch case 2: \ 508d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 509d71ae5a4SJacob Faibussowitsch case 1: \ 510d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 511d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 5127b42bb93SJunchao Zhang } \ 5139371c9d4SSatish Balay while (nnz2 > 0) { \ 5149371c9d4SSatish Balay sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 5159371c9d4SSatish Balay xv += 4; \ 5169371c9d4SSatish Balay xi += 4; \ 5179371c9d4SSatish Balay nnz2 -= 4; \ 5189371c9d4SSatish Balay } \ 5199371c9d4SSatish Balay xv -= nnz; \ 5209371c9d4SSatish Balay xi -= nnz; \ 5217b42bb93SJunchao Zhang } \ 5227b42bb93SJunchao Zhang } 523003131ecSBarry Smith 524003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 5259371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 5269371c9d4SSatish Balay { \ 527003131ecSBarry Smith PetscInt __i, __i1, __i2; \ 5289371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 5299371c9d4SSatish Balay __i1 = xi[__i]; \ 5309371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 5319371c9d4SSatish Balay sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 5329371c9d4SSatish Balay } \ 5339371c9d4SSatish Balay if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \ 5349371c9d4SSatish Balay } 535003131ecSBarry Smith 53699acd6aaSStefano Zampini #elif defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND) 53754e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz)) 538d68ff82bSRichard Tran Mills 539003131ecSBarry Smith #else 5409371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 5419371c9d4SSatish Balay { \ 542003131ecSBarry Smith PetscInt __i; \ 5439371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \ 5449371c9d4SSatish Balay } 545003131ecSBarry Smith #endif 546003131ecSBarry Smith 54799acd6aaSStefano Zampini #if defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES) && !defined(PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND) 54854e8760dSRichard Tran Mills #include <immintrin.h> 54954e8760dSRichard Tran Mills #if !defined(_MM_SCALE_8) 55054e8760dSRichard Tran Mills #define _MM_SCALE_8 8 55154e8760dSRichard Tran Mills #endif 55254e8760dSRichard Tran Mills 553d71ae5a4SJacob Faibussowitsch static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n) 554d71ae5a4SJacob Faibussowitsch { 55554e8760dSRichard Tran Mills __m512d vec_x, vec_y, vec_vals; 55654e8760dSRichard Tran Mills __m256i vec_idx; 55754e8760dSRichard Tran Mills PetscInt j; 55854e8760dSRichard Tran Mills 55954e8760dSRichard Tran Mills vec_y = _mm512_setzero_pd(); 56054e8760dSRichard Tran Mills for (j = 0; j < (n >> 3); j++) { 561ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const *)aj); 562ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 56354e8760dSRichard Tran Mills vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8); 56454e8760dSRichard Tran Mills vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y); 5659371c9d4SSatish Balay aj += 8; 5669371c9d4SSatish Balay aa += 8; 56754e8760dSRichard Tran Mills } 568851da9e1SHong Zhang #if defined(__AVX512VL__) 569851da9e1SHong Zhang /* masked load requires avx512vl, which is not supported by KNL */ 570851da9e1SHong Zhang if (n & 0x07) { 571*47a5f4f7SHong Zhang __mmask8 mask; 57254e8760dSRichard Tran Mills mask = (__mmask8)(0xff >> (8 - (n & 0x07))); 573851da9e1SHong Zhang vec_idx = _mm256_mask_loadu_epi32(vec_idx, mask, aj); 574851da9e1SHong Zhang vec_vals = _mm512_mask_loadu_pd(vec_vals, mask, aa); 57554e8760dSRichard Tran Mills vec_x = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8); 57654e8760dSRichard Tran Mills vec_y = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask); 57754e8760dSRichard Tran Mills } 578851da9e1SHong Zhang *sum += _mm512_reduce_add_pd(vec_y); 579851da9e1SHong Zhang #else 580851da9e1SHong Zhang *sum += _mm512_reduce_add_pd(vec_y); 58154e8760dSRichard Tran Mills for (j = 0; j < (n & 0x07); j++) *sum += aa[j] * x[aj[j]]; 582851da9e1SHong Zhang #endif 58354e8760dSRichard Tran Mills } 58454e8760dSRichard Tran Mills #endif 585b434eb95SMatthew G. Knepley 586b434eb95SMatthew G. Knepley /* 587b434eb95SMatthew G. Knepley PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 588b434eb95SMatthew G. Knepley 589b434eb95SMatthew G. Knepley Input Parameters: 590b434eb95SMatthew G. Knepley + nnz - the number of entries 591b434eb95SMatthew G. Knepley . r - the array of vector values 592b434eb95SMatthew G. Knepley . xv - the matrix values for the row 593b434eb95SMatthew G. Knepley - xi - the column indices of the nonzeros in the row 594b434eb95SMatthew G. Knepley 595b434eb95SMatthew G. Knepley Output Parameter: 596b434eb95SMatthew G. Knepley . max - the max of results 597b434eb95SMatthew G. Knepley 598db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()` 599b434eb95SMatthew G. Knepley */ 6009371c9d4SSatish Balay #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \ 6019371c9d4SSatish Balay { \ 602b434eb95SMatthew G. Knepley PetscInt __i; \ 6039371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]])); \ 6049371c9d4SSatish Balay } 605b434eb95SMatthew G. Knepley 6064b38b95cSHong Zhang /* 6074b38b95cSHong Zhang Add column indices into table for counting the max nonzeros of merged rows 6084b38b95cSHong Zhang */ 6099371c9d4SSatish Balay #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \ 6109371c9d4SSatish Balay { \ 6114b38b95cSHong Zhang PetscInt _j, _row, _nz, *_col; \ 612ec07b8f8SHong Zhang if (mat) { \ 6134b38b95cSHong Zhang for (_row = 0; _row < nrows; _row++) { \ 6144b38b95cSHong Zhang _nz = mat->i[_row + 1] - mat->i[_row]; \ 6154b38b95cSHong Zhang for (_j = 0; _j < _nz; _j++) { \ 6164b38b95cSHong Zhang _col = _j + mat->j + mat->i[_row]; \ 6174b38b95cSHong Zhang PetscTableAdd(ta, *_col + 1, 1, INSERT_VALUES); \ 6184b38b95cSHong Zhang } \ 6194b38b95cSHong Zhang } \ 620ec07b8f8SHong Zhang } \ 6214b38b95cSHong Zhang } 6224b38b95cSHong Zhang 6230ca7d551SHong Zhang /* 6240ca7d551SHong Zhang Add column indices into table for counting the nonzeros of merged rows 6250ca7d551SHong Zhang */ 6269371c9d4SSatish Balay #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \ 6279371c9d4SSatish Balay { \ 6280ca7d551SHong Zhang PetscInt _j, _row, _nz, *_col, _i; \ 6290ca7d551SHong Zhang for (_i = 0; _i < nrows; _i++) { \ 6300ca7d551SHong Zhang _row = rows[_i]; \ 6310ca7d551SHong Zhang _nz = mat->i[_row + 1] - mat->i[_row]; \ 6320ca7d551SHong Zhang for (_j = 0; _j < _nz; _j++) { \ 6330ca7d551SHong Zhang _col = _j + mat->j + mat->i[_row]; \ 6340ca7d551SHong Zhang PetscTableAdd(ta, *_col + 1, 1, INSERT_VALUES); \ 6350ca7d551SHong Zhang } \ 6360ca7d551SHong Zhang } \ 6370ca7d551SHong Zhang } 6380ca7d551SHong Zhang 6392d40f771SBarry Smith #endif 640