1b8a66259SBarry Smith 26524c165SJacob Faibussowitsch #ifndef __AIJ_H 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5af0996ceSBarry Smith #include <petsc/private/matimpl.h> 6eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h> 78f690400SShri Abhyankar 8*5b17b1ebSStefano Zampini /* 9*5b17b1ebSStefano Zampini Used by MatCreateSubMatrices_MPIXAIJ_Local() 10*5b17b1ebSStefano Zampini */ 11*5b17b1ebSStefano Zampini typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */ 12*5b17b1ebSStefano Zampini PetscInt id; /* index of submats, only submats[0] is responsible for deleting some arrays below */ 13*5b17b1ebSStefano Zampini PetscInt nrqs, nrqr; 14*5b17b1ebSStefano Zampini PetscInt **rbuf1, **rbuf2, **rbuf3, **sbuf1, **sbuf2; 15*5b17b1ebSStefano Zampini PetscInt **ptr; 16*5b17b1ebSStefano Zampini PetscInt *tmp; 17*5b17b1ebSStefano Zampini PetscInt *ctr; 18*5b17b1ebSStefano Zampini PetscInt *pa; /* proc array */ 19*5b17b1ebSStefano Zampini PetscInt *req_size, *req_source1, *req_source2; 20*5b17b1ebSStefano Zampini PetscBool allcolumns, allrows; 21*5b17b1ebSStefano Zampini PetscBool singleis; 22*5b17b1ebSStefano Zampini PetscInt *row2proc; /* row to proc map */ 23*5b17b1ebSStefano Zampini PetscInt nstages; 24*5b17b1ebSStefano Zampini #if defined(PETSC_USE_CTABLE) 25*5b17b1ebSStefano Zampini PetscHMapI cmap, rmap; 26*5b17b1ebSStefano Zampini PetscInt *cmap_loc, *rmap_loc; 27*5b17b1ebSStefano Zampini #else 28*5b17b1ebSStefano Zampini PetscInt *cmap, *rmap; 29*5b17b1ebSStefano Zampini #endif 30*5b17b1ebSStefano Zampini PetscErrorCode (*destroy)(Mat); 31*5b17b1ebSStefano Zampini } Mat_SubSppt; 32*5b17b1ebSStefano Zampini 33d67d9f35SJunchao Zhang /* Operations provided by MATSEQAIJ and its subclasses */ 34d67d9f35SJunchao Zhang typedef struct { 35d67d9f35SJunchao Zhang PetscErrorCode (*getarray)(Mat, PetscScalar **); 36d67d9f35SJunchao Zhang PetscErrorCode (*restorearray)(Mat, PetscScalar **); 37d67d9f35SJunchao Zhang PetscErrorCode (*getarrayread)(Mat, const PetscScalar **); 38d67d9f35SJunchao Zhang PetscErrorCode (*restorearrayread)(Mat, const PetscScalar **); 39d67d9f35SJunchao Zhang PetscErrorCode (*getarraywrite)(Mat, PetscScalar **); 40d67d9f35SJunchao Zhang PetscErrorCode (*restorearraywrite)(Mat, PetscScalar **); 417ee59b9bSJunchao Zhang PetscErrorCode (*getcsrandmemtype)(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *); 42d67d9f35SJunchao Zhang } Mat_SeqAIJOps; 43d67d9f35SJunchao Zhang 44e6b907acSBarry Smith /* 45e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 46e6b907acSBarry Smith */ 47421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 48ace3abfcSBarry Smith PetscBool roworiented; /* if true, row-oriented input, default */ \ 49e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 5028b2fa4aSMatthew Knepley PetscInt nounused; /* -1 generate error on unused space */ \ 51ace3abfcSBarry Smith PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \ 52e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */ \ 53e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */ \ 54e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */ \ 55846b4da1SFande Kong PetscInt *ipre; /* space preallocated for each row by user */ \ 56ace3abfcSBarry Smith PetscBool free_imax_ilen; \ 57e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 58e6b907acSBarry Smith as more values are set than were prealloced */ \ 59e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */ \ 60ace3abfcSBarry Smith PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \ 61ace3abfcSBarry Smith PetscBool ignorezeroentries; \ 62ace3abfcSBarry Smith PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 63ace3abfcSBarry Smith PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 64e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 65e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 66e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 67e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 68e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 697b083b7cSBarry Smith PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 70ace3abfcSBarry Smith PetscBool free_diag; \ 71421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 72e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 734fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 74ace3abfcSBarry Smith PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 7517df9f7cSHong Zhang Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \ 7617df9f7cSHong Zhang means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \ 77d67d9f35SJunchao Zhang Mat_SubSppt *submatis1; /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \ 78d67d9f35SJunchao Zhang Mat_SeqAIJOps ops[1] /* operations for SeqAIJ and its subclasses */ 79e6b907acSBarry Smith 8053565b12SHong Zhang typedef struct { 8153565b12SHong Zhang MatTransposeColoring matcoloring; 8253565b12SHong Zhang Mat Bt_den; /* dense matrix of B^T */ 8353565b12SHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 8453565b12SHong Zhang PetscBool usecoloring; 8553565b12SHong Zhang } Mat_MatMatTransMult; 8653565b12SHong Zhang 876d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */ 886d373c3eSHong Zhang Mat At; /* transpose of the first matrix */ 892cff0574SHong Zhang Mat mA; /* maij matrix of A */ 902cff0574SHong Zhang Vec bt, ct; /* vectors to hold locally transposed arrays of B and C */ 916718818eSStefano Zampini /* used by PtAP */ 926718818eSStefano Zampini void *data; 936718818eSStefano Zampini PetscErrorCode (*destroy)(void *); 942cff0574SHong Zhang } Mat_MatTransMatMult; 952cff0574SHong Zhang 9653565b12SHong Zhang typedef struct { 9753565b12SHong Zhang PetscInt *api, *apj; /* symbolic structure of A*P */ 9853565b12SHong Zhang PetscScalar *apa; /* temporary array for storing one row of A*P */ 993cdca5ebSHong Zhang } Mat_AP; 10053565b12SHong Zhang 10153565b12SHong Zhang typedef struct { 10253565b12SHong Zhang MatTransposeColoring matcoloring; 103257c235dSHong Zhang Mat Rt; /* sparse or dense matrix of R^T */ 10453565b12SHong Zhang Mat RARt; /* dense matrix of R*A*R^T */ 1053b1b9624SHong Zhang Mat ARt; /* A*R^T used for the case -matrart_color_art */ 10653565b12SHong Zhang MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 1076718818eSStefano Zampini /* free intermediate products needed for PtAP */ 1086718818eSStefano Zampini void *data; 1096718818eSStefano Zampini PetscErrorCode (*destroy)(void *); 11053565b12SHong Zhang } Mat_RARt; 11153565b12SHong Zhang 1126d0b6147SHong Zhang typedef struct { 1136d0b6147SHong Zhang Mat BC; /* temp matrix for storing B*C */ 1146d0b6147SHong Zhang } Mat_MatMatMatMult; 1156d0b6147SHong Zhang 1162d40f771SBarry Smith /* 117ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 118e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 119dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 1205768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 1212d40f771SBarry Smith */ 122d35516d3SLois Curfman McInnes 123e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 124b8a66259SBarry Smith typedef struct { 1254108e4d5SBarry Smith MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 126f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 127ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 128f0d39aaaSBarry Smith 129ace3abfcSBarry Smith PetscBool use; 130e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 131e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 132e6b907acSBarry Smith PetscInt limit; /* inode limit */ 133e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 134ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 135a02bda8eSBarry Smith PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 1364108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 137e6b907acSBarry Smith 1385a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer); 1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType); 1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 1425a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool); 1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *); 1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool); 1455a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat, Mat, const MatFactorInfo *); 1465a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *); 147f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **); 148f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **); 149e6b907acSBarry Smith 150e6b907acSBarry Smith typedef struct { 15154f21887SBarry Smith SEQAIJHEADER(MatScalar); 1524108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 15354f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 15471f1c65dSBarry Smith 15571f1c65dSBarry Smith PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 156ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 157bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 1582291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 15961ecd0c6SBarry Smith PetscBool diagonaldense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */ 16071f1c65dSBarry Smith PetscScalar fshift, omega; /* last used omega and fshift */ 161394ed5ebSJunchao Zhang 162394ed5ebSJunchao Zhang /* MatSetValuesCOO() related fields on host */ 163394ed5ebSJunchao Zhang PetscCount coo_n; /* Number of entries in MatSetPreallocationCOO() */ 164394ed5ebSJunchao Zhang PetscCount Atot; /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */ 165394ed5ebSJunchao Zhang PetscCount *jmap; /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */ 166394ed5ebSJunchao Zhang PetscCount *perm; /* The permutation array in sorting (i,j) by row and then by col */ 167ec8511deSBarry Smith } Mat_SeqAIJ; 168b8a66259SBarry Smith 169e6b907acSBarry Smith /* 170e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 171e6b907acSBarry Smith */ 172d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i) 173d71ae5a4SJacob Faibussowitsch { 174e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data; 175e6b907acSBarry Smith if (A->singlemalloc) { 1769566063dSJacob Faibussowitsch PetscCall(PetscFree3(*a, *j, *i)); 177e6b907acSBarry Smith } else { 1789566063dSJacob Faibussowitsch if (A->free_a) PetscCall(PetscFree(*a)); 1799566063dSJacob Faibussowitsch if (A->free_ij) PetscCall(PetscFree(*j)); 1809566063dSJacob Faibussowitsch if (A->free_ij) PetscCall(PetscFree(*i)); 181e6b907acSBarry Smith } 182e6b907acSBarry Smith return 0; 183e6b907acSBarry Smith } 184e6b907acSBarry Smith /* 185e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 186357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 187e6b907acSBarry Smith */ 188fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \ 189fef13f97SBarry Smith if (NROW >= RMAX) { \ 190fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 191fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 192f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 193fef13f97SBarry Smith datatype *new_a; \ 194fef13f97SBarry Smith \ 19508401ef6SPierre 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); \ 196fef13f97SBarry Smith /* malloc new storage space */ \ 1979566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(BS2 *new_nz, &new_a, new_nz, &new_j, AM + 1, &new_i)); \ 198fef13f97SBarry Smith \ 199fef13f97SBarry Smith /* copy over old data into new slots */ \ 200ad540459SPierre Jolivet for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 201ad540459SPierre Jolivet for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 2029566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 203fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 2049566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \ 2059566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \ 2069566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \ 2079566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), AA + BS2 * (AI[ROW] + NROW), BS2 * len)); \ 208fef13f97SBarry Smith /* free up old matrix storage */ \ 2099566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 210fef13f97SBarry Smith AA = new_a; \ 211fef13f97SBarry Smith Ain->a = (MatScalar *)new_a; \ 2129371c9d4SSatish Balay AI = Ain->i = new_i; \ 2139371c9d4SSatish Balay AJ = Ain->j = new_j; \ 214fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 215fef13f97SBarry Smith \ 2169371c9d4SSatish Balay RP = AJ + AI[ROW]; \ 2179371c9d4SSatish Balay AP = AA + BS2 * AI[ROW]; \ 218fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 219fef13f97SBarry Smith Ain->maxnz += BS2 * CHUNKSIZE; \ 220fef13f97SBarry Smith Ain->reallocs++; \ 2219371c9d4SSatish Balay } 22217454e89SShri Abhyankar 223876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \ 224720833daSHong Zhang if (NROW >= RMAX) { \ 225720833daSHong Zhang Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 226720833daSHong Zhang /* there is no extra room in row, therefore enlarge */ \ 227f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 228720833daSHong Zhang \ 22908401ef6SPierre 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); \ 230720833daSHong Zhang /* malloc new storage space */ \ 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_nz, &new_j)); \ 2329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(AM + 1, &new_i)); \ 233720833daSHong Zhang \ 234720833daSHong Zhang /* copy over old data into new slots */ \ 235ad540459SPierre Jolivet for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 236ad540459SPierre Jolivet for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 2379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 238720833daSHong Zhang len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 2399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \ 240876c6284SHong Zhang \ 241720833daSHong Zhang /* free up old matrix storage */ \ 2429566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 243876c6284SHong Zhang Ain->a = NULL; \ 2449371c9d4SSatish Balay AI = Ain->i = new_i; \ 2459371c9d4SSatish Balay AJ = Ain->j = new_j; \ 246720833daSHong Zhang Ain->singlemalloc = PETSC_FALSE; \ 247876c6284SHong Zhang Ain->free_a = PETSC_FALSE; \ 248720833daSHong Zhang \ 249876c6284SHong Zhang RP = AJ + AI[ROW]; \ 250720833daSHong Zhang RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 251720833daSHong Zhang Ain->maxnz += BS2 * CHUNKSIZE; \ 252720833daSHong Zhang Ain->reallocs++; \ 2539371c9d4SSatish Balay } 254e6b907acSBarry Smith 255cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *); 256e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]); 257cbc6b225SStefano Zampini PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_SeqAIJ(Mat); 258cbc6b225SStefano Zampini 2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *); 2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *); 2621df811f5SHong Zhang 2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *); 2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *); 2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *); 2705a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure); 2715a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *); 2725a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2735a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **); 27408480c60SBarry Smith 275b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec); 276b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec); 277b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec); 278b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec); 279b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec); 280b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 2815a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 282b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 28308480c60SBarry Smith 2845a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool); 2853a7fca6bSBarry Smith 2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 2887fb60732SBarry Smith PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]); 2895a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *); 2905008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *); 2917fb60732SBarry Smith 2922462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **); 2935a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *); 2945a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 2955a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *); 2985a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *); 2995a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec); 3005a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec); 3015a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat, Vec, Vec); 3025a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec); 3035a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat, Vec, Vec); 3045a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec); 3055a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat, Vec, Vec); 3065a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec); 3075a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec); 3085a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec); 3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec); 3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec); 3115a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 3125a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat, Mat, Mat); 3135a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat); 31479c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *); 31593dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring); 316f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring); 317a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt); 31852f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer); 31952f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer); 3205a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer); 3215a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 3227bab7c10SHong Zhang 323df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE) 3244222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat); 325df97dc6dSFande Kong #endif 3264222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat); 3274222ddf1SHong Zhang 3284222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat); 3294222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat); 3304222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat); 3314222ddf1SHong Zhang 3324222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3334222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat); 3344222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat); 3354222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat); 3364222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat); 3374222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat); 3384222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat); 3394222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat); 3404222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat); 3414222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE) 3424222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat); 3434222ddf1SHong Zhang #endif 3444222ddf1SHong Zhang 3455a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 346df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat); 3474222ddf1SHong Zhang 3484099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat); 3495a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat); 3502b8ad9a3SHong Zhang 3514222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat); 3525a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat); 35453565b12SHong Zhang 3554222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3564222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat); 3574222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat); 3585a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 35955bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat); 36055bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat); 3615df89d91SHong Zhang 3624222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3635a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3646718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *); 3653bf78175SHong Zhang 3664222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3675a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3685a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring); 3695a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat); 3705a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat); 3715df89d91SHong Zhang 3724222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat); 3735a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat); 3747bab7c10SHong Zhang 375679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom); 3765a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 3775a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 3785a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 379db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar); 38087c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec); 38187c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode); 3825a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure); 3835a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3845a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3855a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3865a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 3876378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 3886378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 3895a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 3905a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat); 3915a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer); 3929af31e4aSHong Zhang 3935a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 3945a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 395a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 396a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 397019b515eSShri Abhyankar 3985a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *); 3999f5f6813SShri Abhyankar 400d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 401388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *); 402388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *); 403f2fbf96bSVaclav Hapla #endif 404cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *); 405cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *); 406388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *); 407388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *); 408388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *); 409d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK) 410d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *); 411d24d4204SJose E. Roman #endif 412388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *); 413cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *); 4144dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *); 4154a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *); 416388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *); 4175a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS); 4185a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *); 4198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 4205a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType); 4213fa6b06aSMark Adams PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat); 4222f947c57SVictor Minden 423b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *); 4249c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 4259c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 42670f19b1fSKris Buschelman 42753dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat); 428dec0b466SHong Zhang PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A); 429f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *); 4300fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat); 431f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat); 43263a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]); 433feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *); 43453dd7562SDmitry Karpeev 435a3bb6f32SFande Kong PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *); 436e4e71118SStefano Zampini PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat); 437a3bb6f32SFande Kong 438003131ecSBarry Smith /* 439003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 440003131ecSBarry Smith 441003131ecSBarry Smith Input Parameters: 442003131ecSBarry Smith + nnz - the number of entries 443003131ecSBarry Smith . r - the array of vector values 444003131ecSBarry Smith . xv - the matrix values for the row 445003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 446003131ecSBarry Smith 447003131ecSBarry Smith Output Parameter: 448003131ecSBarry Smith . sum - negative the sum of results 449003131ecSBarry Smith 450003131ecSBarry Smith PETSc compile flags: 4517b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4527b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4537b42bb93SJunchao Zhang 45411a5261eSBarry Smith Developer Note: 4557b42bb93SJunchao Zhang The macro changes sum but not other parameters 456003131ecSBarry Smith 457db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()` 458003131ecSBarry Smith */ 459519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 4609371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 4619371c9d4SSatish Balay { \ 462003131ecSBarry Smith if (nnz > 0) { \ 4637b42bb93SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 4647b42bb93SJunchao Zhang switch (rem) { \ 465d71ae5a4SJacob Faibussowitsch case 3: \ 466d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 467d71ae5a4SJacob Faibussowitsch case 2: \ 468d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 469d71ae5a4SJacob Faibussowitsch case 1: \ 470d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 471d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 4727b42bb93SJunchao Zhang } \ 4739371c9d4SSatish Balay while (nnz2 > 0) { \ 4749371c9d4SSatish Balay sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 4759371c9d4SSatish Balay xv += 4; \ 4769371c9d4SSatish Balay xi += 4; \ 4779371c9d4SSatish Balay nnz2 -= 4; \ 4789371c9d4SSatish Balay } \ 4799371c9d4SSatish Balay xv -= nnz; \ 4809371c9d4SSatish Balay xi -= nnz; \ 4817b42bb93SJunchao Zhang } \ 4827b42bb93SJunchao Zhang } 483003131ecSBarry Smith 484003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 4859371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 4869371c9d4SSatish Balay { \ 487003131ecSBarry Smith PetscInt __i, __i1, __i2; \ 4889371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 4899371c9d4SSatish Balay __i1 = xi[__i]; \ 4909371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 4919371c9d4SSatish Balay sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 4929371c9d4SSatish Balay } \ 4939371c9d4SSatish Balay if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \ 4949371c9d4SSatish Balay } 495003131ecSBarry Smith 496003131ecSBarry Smith #else 4979371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 4989371c9d4SSatish Balay { \ 499003131ecSBarry Smith PetscInt __i; \ 5009371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \ 5019371c9d4SSatish Balay } 502003131ecSBarry Smith #endif 503003131ecSBarry Smith 504003131ecSBarry Smith /* 505003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 506003131ecSBarry Smith 507003131ecSBarry Smith Input Parameters: 508003131ecSBarry Smith + nnz - the number of entries 509003131ecSBarry Smith . r - the array of vector values 510003131ecSBarry Smith . xv - the matrix values for the row 511003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 512003131ecSBarry Smith 513003131ecSBarry Smith Output Parameter: 514003131ecSBarry Smith . sum - the sum of results 515003131ecSBarry Smith 516003131ecSBarry Smith PETSc compile flags: 5177b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 5187b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 5197b42bb93SJunchao Zhang 52011a5261eSBarry Smith Developer Note: 5217b42bb93SJunchao Zhang The macro changes sum but not other parameters 522003131ecSBarry Smith 523db781477SPatrick Sanan .seealso: `PetscSparseDenseMinusDot()` 524003131ecSBarry Smith */ 525519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 5269371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 5279371c9d4SSatish Balay { \ 528003131ecSBarry Smith if (nnz > 0) { \ 5297b42bb93SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 5307b42bb93SJunchao Zhang switch (rem) { \ 531d71ae5a4SJacob Faibussowitsch case 3: \ 532d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 533d71ae5a4SJacob Faibussowitsch case 2: \ 534d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 535d71ae5a4SJacob Faibussowitsch case 1: \ 536d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 537d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 5387b42bb93SJunchao Zhang } \ 5399371c9d4SSatish Balay while (nnz2 > 0) { \ 5409371c9d4SSatish Balay sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 5419371c9d4SSatish Balay xv += 4; \ 5429371c9d4SSatish Balay xi += 4; \ 5439371c9d4SSatish Balay nnz2 -= 4; \ 5449371c9d4SSatish Balay } \ 5459371c9d4SSatish Balay xv -= nnz; \ 5469371c9d4SSatish Balay xi -= nnz; \ 5477b42bb93SJunchao Zhang } \ 5487b42bb93SJunchao Zhang } 549003131ecSBarry Smith 550003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 5519371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 5529371c9d4SSatish Balay { \ 553003131ecSBarry Smith PetscInt __i, __i1, __i2; \ 5549371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 5559371c9d4SSatish Balay __i1 = xi[__i]; \ 5569371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 5579371c9d4SSatish Balay sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 5589371c9d4SSatish Balay } \ 5599371c9d4SSatish Balay if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \ 5609371c9d4SSatish Balay } 561003131ecSBarry Smith 56299acd6aaSStefano 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) 56354e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz)) 564d68ff82bSRichard Tran Mills 565003131ecSBarry Smith #else 5669371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 5679371c9d4SSatish Balay { \ 568003131ecSBarry Smith PetscInt __i; \ 5699371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \ 5709371c9d4SSatish Balay } 571003131ecSBarry Smith #endif 572003131ecSBarry Smith 57399acd6aaSStefano 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) 57454e8760dSRichard Tran Mills #include <immintrin.h> 57554e8760dSRichard Tran Mills #if !defined(_MM_SCALE_8) 57654e8760dSRichard Tran Mills #define _MM_SCALE_8 8 57754e8760dSRichard Tran Mills #endif 57854e8760dSRichard Tran Mills 579d71ae5a4SJacob Faibussowitsch static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n) 580d71ae5a4SJacob Faibussowitsch { 58154e8760dSRichard Tran Mills __m512d vec_x, vec_y, vec_vals; 58254e8760dSRichard Tran Mills __m256i vec_idx; 58354e8760dSRichard Tran Mills __mmask8 mask; 58454e8760dSRichard Tran Mills PetscInt j; 58554e8760dSRichard Tran Mills 58654e8760dSRichard Tran Mills vec_y = _mm512_setzero_pd(); 58754e8760dSRichard Tran Mills for (j = 0; j < (n >> 3); j++) { 588ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const *)aj); 589ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 59054e8760dSRichard Tran Mills vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8); 59154e8760dSRichard Tran Mills vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y); 5929371c9d4SSatish Balay aj += 8; 5939371c9d4SSatish Balay aa += 8; 59454e8760dSRichard Tran Mills } 59554e8760dSRichard Tran Mills /* masked load does not work on KNL, it requires avx512vl */ 59654e8760dSRichard Tran Mills if ((n & 0x07) > 2) { 59754e8760dSRichard Tran Mills mask = (__mmask8)(0xff >> (8 - (n & 0x07))); 598ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const *)aj); 599ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 60054e8760dSRichard Tran Mills vec_x = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8); 60154e8760dSRichard Tran Mills vec_y = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask); 60254e8760dSRichard Tran Mills } else if ((n & 0x07) == 2) { 60354e8760dSRichard Tran Mills *sum += aa[0] * x[aj[0]]; 60454e8760dSRichard Tran Mills *sum += aa[1] * x[aj[1]]; 60554e8760dSRichard Tran Mills } else if ((n & 0x07) == 1) { 60654e8760dSRichard Tran Mills *sum += aa[0] * x[aj[0]]; 60754e8760dSRichard Tran Mills } 60854e8760dSRichard Tran Mills if (n > 2) *sum += _mm512_reduce_add_pd(vec_y); 60954e8760dSRichard Tran Mills /* 61054e8760dSRichard Tran Mills for (j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]]; 61154e8760dSRichard Tran Mills */ 61254e8760dSRichard Tran Mills } 61354e8760dSRichard Tran Mills #endif 614b434eb95SMatthew G. Knepley 615b434eb95SMatthew G. Knepley /* 616b434eb95SMatthew G. Knepley PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 617b434eb95SMatthew G. Knepley 618b434eb95SMatthew G. Knepley Input Parameters: 619b434eb95SMatthew G. Knepley + nnz - the number of entries 620b434eb95SMatthew G. Knepley . r - the array of vector values 621b434eb95SMatthew G. Knepley . xv - the matrix values for the row 622b434eb95SMatthew G. Knepley - xi - the column indices of the nonzeros in the row 623b434eb95SMatthew G. Knepley 624b434eb95SMatthew G. Knepley Output Parameter: 625b434eb95SMatthew G. Knepley . max - the max of results 626b434eb95SMatthew G. Knepley 627db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()` 628b434eb95SMatthew G. Knepley */ 6299371c9d4SSatish Balay #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \ 630eec179cfSJacob Faibussowitsch do { \ 631eec179cfSJacob Faibussowitsch for (PetscInt __i = 0; __i < (nnz); __i++) { max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); } \ 632eec179cfSJacob Faibussowitsch } while (0) 633b434eb95SMatthew G. Knepley 6344b38b95cSHong Zhang /* 6354b38b95cSHong Zhang Add column indices into table for counting the max nonzeros of merged rows 6364b38b95cSHong Zhang */ 6379371c9d4SSatish Balay #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \ 638eec179cfSJacob Faibussowitsch do { \ 639eec179cfSJacob Faibussowitsch if ((mat)) { \ 640eec179cfSJacob Faibussowitsch for (PetscInt _row = 0; _row < (nrows); _row++) { \ 641eec179cfSJacob Faibussowitsch const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \ 642eec179cfSJacob Faibussowitsch for (PetscInt _j = 0; _j < _nz; _j++) { \ 643eec179cfSJacob Faibussowitsch PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \ 644c76ffc5fSJacob Faibussowitsch PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \ 6454b38b95cSHong Zhang } \ 6464b38b95cSHong Zhang } \ 647ec07b8f8SHong Zhang } \ 648eec179cfSJacob Faibussowitsch } while (0) 6494b38b95cSHong Zhang 6500ca7d551SHong Zhang /* 6510ca7d551SHong Zhang Add column indices into table for counting the nonzeros of merged rows 6520ca7d551SHong Zhang */ 6539371c9d4SSatish Balay #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \ 654eec179cfSJacob Faibussowitsch do { \ 655eec179cfSJacob Faibussowitsch for (PetscInt _i = 0; _i < (nrows); _i++) { \ 656eec179cfSJacob Faibussowitsch const PetscInt _row = (rows)[_i]; \ 657eec179cfSJacob Faibussowitsch const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \ 658eec179cfSJacob Faibussowitsch for (PetscInt _j = 0; _j < _nz; _j++) { \ 659eec179cfSJacob Faibussowitsch PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \ 660eec179cfSJacob Faibussowitsch PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \ 6610ca7d551SHong Zhang } \ 6620ca7d551SHong Zhang } \ 663eec179cfSJacob Faibussowitsch } while (0) 6640ca7d551SHong Zhang 6652d40f771SBarry Smith #endif 666