xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1*3ba16761SJacob Faibussowitsch #ifndef PETSC_MATAIJ_IMPL_H
2*3ba16761SJacob Faibussowitsch #define PETSC_MATAIJ_IMPL_H
3e6b907acSBarry Smith 
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h>
68f690400SShri Abhyankar 
75b17b1ebSStefano Zampini /*
85b17b1ebSStefano Zampini  Used by MatCreateSubMatrices_MPIXAIJ_Local()
95b17b1ebSStefano Zampini */
105b17b1ebSStefano Zampini typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
115b17b1ebSStefano Zampini   PetscInt   id; /* index of submats, only submats[0] is responsible for deleting some arrays below */
125b17b1ebSStefano Zampini   PetscInt   nrqs, nrqr;
135b17b1ebSStefano Zampini   PetscInt **rbuf1, **rbuf2, **rbuf3, **sbuf1, **sbuf2;
145b17b1ebSStefano Zampini   PetscInt **ptr;
155b17b1ebSStefano Zampini   PetscInt  *tmp;
165b17b1ebSStefano Zampini   PetscInt  *ctr;
175b17b1ebSStefano Zampini   PetscInt  *pa; /* proc array */
185b17b1ebSStefano Zampini   PetscInt  *req_size, *req_source1, *req_source2;
195b17b1ebSStefano Zampini   PetscBool  allcolumns, allrows;
205b17b1ebSStefano Zampini   PetscBool  singleis;
215b17b1ebSStefano Zampini   PetscInt  *row2proc; /* row to proc map */
225b17b1ebSStefano Zampini   PetscInt   nstages;
235b17b1ebSStefano Zampini #if defined(PETSC_USE_CTABLE)
245b17b1ebSStefano Zampini   PetscHMapI cmap, rmap;
255b17b1ebSStefano Zampini   PetscInt  *cmap_loc, *rmap_loc;
265b17b1ebSStefano Zampini #else
275b17b1ebSStefano Zampini   PetscInt *cmap, *rmap;
285b17b1ebSStefano Zampini #endif
295b17b1ebSStefano Zampini   PetscErrorCode (*destroy)(Mat);
305b17b1ebSStefano Zampini } Mat_SubSppt;
315b17b1ebSStefano Zampini 
32d67d9f35SJunchao Zhang /* Operations provided by MATSEQAIJ and its subclasses */
33d67d9f35SJunchao Zhang typedef struct {
34d67d9f35SJunchao Zhang   PetscErrorCode (*getarray)(Mat, PetscScalar **);
35d67d9f35SJunchao Zhang   PetscErrorCode (*restorearray)(Mat, PetscScalar **);
36d67d9f35SJunchao Zhang   PetscErrorCode (*getarrayread)(Mat, const PetscScalar **);
37d67d9f35SJunchao Zhang   PetscErrorCode (*restorearrayread)(Mat, const PetscScalar **);
38d67d9f35SJunchao Zhang   PetscErrorCode (*getarraywrite)(Mat, PetscScalar **);
39d67d9f35SJunchao Zhang   PetscErrorCode (*restorearraywrite)(Mat, PetscScalar **);
407ee59b9bSJunchao Zhang   PetscErrorCode (*getcsrandmemtype)(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *);
41d67d9f35SJunchao Zhang } Mat_SeqAIJOps;
42d67d9f35SJunchao Zhang 
43e6b907acSBarry Smith /*
44e6b907acSBarry Smith     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
45e6b907acSBarry Smith */
46421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \
47ace3abfcSBarry Smith   PetscBool         roworiented;  /* if true, row-oriented input, default */ \
48e6b907acSBarry Smith   PetscInt          nonew;        /* 1 don't add new nonzeros, -1 generate error on new */ \
4928b2fa4aSMatthew Knepley   PetscInt          nounused;     /* -1 generate error on unused space */ \
50ace3abfcSBarry Smith   PetscBool         singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \
51e6b907acSBarry Smith   PetscInt          maxnz;        /* allocated nonzeros */ \
52e6b907acSBarry Smith   PetscInt         *imax;         /* maximum space allocated for each row */ \
53e6b907acSBarry Smith   PetscInt         *ilen;         /* actual length of each row */ \
54846b4da1SFande Kong   PetscInt         *ipre;         /* space preallocated for each row by user */ \
55ace3abfcSBarry Smith   PetscBool         free_imax_ilen; \
56e6b907acSBarry Smith   PetscInt          reallocs;           /* number of mallocs done during MatSetValues() \
57e6b907acSBarry Smith                                         as more values are set than were prealloced */ \
58e6b907acSBarry Smith   PetscInt          rmax;               /* max nonzeros in any row */ \
59ace3abfcSBarry Smith   PetscBool         keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \
60ace3abfcSBarry Smith   PetscBool         ignorezeroentries; \
61ace3abfcSBarry Smith   PetscBool         free_ij;       /* free the column indices j and row offsets i when the matrix is destroyed */ \
62ace3abfcSBarry Smith   PetscBool         free_a;        /* free the numerical values when matrix is destroy */ \
63e6b907acSBarry Smith   Mat_CompressedRow compressedrow; /* use compressed row format */ \
64e6b907acSBarry Smith   PetscInt          nz;            /* nonzeros */ \
65e6b907acSBarry Smith   PetscInt         *i;             /* pointer to beginning of each row */ \
66e6b907acSBarry Smith   PetscInt         *j;             /* column values: j + i[k] - 1 is start of row k */ \
67e6b907acSBarry Smith   PetscInt         *diag;          /* pointers to diagonal elements */ \
687b083b7cSBarry Smith   PetscInt          nonzerorowcnt; /* how many rows have nonzero entries */ \
69ace3abfcSBarry Smith   PetscBool         free_diag; \
70421e10b8SBarry Smith   datatype         *a;              /* nonzero elements */ \
71e6b907acSBarry Smith   PetscScalar      *solve_work;     /* work space used in MatSolve */ \
724fd072dbSBarry Smith   IS                row, col, icol; /* index sets, used for reorderings */ \
73ace3abfcSBarry Smith   PetscBool         pivotinblocks;  /* pivot inside factorization of each diagonal block */ \
7417df9f7cSHong Zhang   Mat               parent;         /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
7517df9f7cSHong Zhang                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \
76d67d9f35SJunchao Zhang   Mat_SubSppt      *submatis1;      /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \
77d67d9f35SJunchao Zhang   Mat_SeqAIJOps     ops[1]          /* operations for SeqAIJ and its subclasses */
78e6b907acSBarry Smith 
7953565b12SHong Zhang typedef struct {
8053565b12SHong Zhang   MatTransposeColoring matcoloring;
8153565b12SHong Zhang   Mat                  Bt_den;  /* dense matrix of B^T */
8253565b12SHong Zhang   Mat                  ABt_den; /* dense matrix of A*B^T */
8353565b12SHong Zhang   PetscBool            usecoloring;
8453565b12SHong Zhang } Mat_MatMatTransMult;
8553565b12SHong Zhang 
866d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */
876d373c3eSHong Zhang   Mat At;        /* transpose of the first matrix */
882cff0574SHong Zhang   Mat mA;        /* maij matrix of A */
892cff0574SHong Zhang   Vec bt, ct;    /* vectors to hold locally transposed arrays of B and C */
906718818eSStefano Zampini   /* used by PtAP */
916718818eSStefano Zampini   void *data;
926718818eSStefano Zampini   PetscErrorCode (*destroy)(void *);
932cff0574SHong Zhang } Mat_MatTransMatMult;
942cff0574SHong Zhang 
9553565b12SHong Zhang typedef struct {
9653565b12SHong Zhang   PetscInt    *api, *apj; /* symbolic structure of A*P */
9753565b12SHong Zhang   PetscScalar *apa;       /* temporary array for storing one row of A*P */
983cdca5ebSHong Zhang } Mat_AP;
9953565b12SHong Zhang 
10053565b12SHong Zhang typedef struct {
10153565b12SHong Zhang   MatTransposeColoring matcoloring;
102257c235dSHong Zhang   Mat                  Rt;   /* sparse or dense matrix of R^T */
10353565b12SHong Zhang   Mat                  RARt; /* dense matrix of R*A*R^T */
1043b1b9624SHong Zhang   Mat                  ARt;  /* A*R^T used for the case -matrart_color_art */
10553565b12SHong Zhang   MatScalar           *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
1066718818eSStefano Zampini   /* free intermediate products needed for PtAP */
1076718818eSStefano Zampini   void *data;
1086718818eSStefano Zampini   PetscErrorCode (*destroy)(void *);
10953565b12SHong Zhang } Mat_RARt;
11053565b12SHong Zhang 
1116d0b6147SHong Zhang typedef struct {
1126d0b6147SHong Zhang   Mat BC; /* temp matrix for storing B*C */
1136d0b6147SHong Zhang } Mat_MatMatMatMult;
1146d0b6147SHong Zhang 
1152d40f771SBarry Smith /*
116ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
117e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
118dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
1195768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
1202d40f771SBarry Smith */
121d35516d3SLois Curfman McInnes 
122e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
123b8a66259SBarry Smith typedef struct {
1244108e4d5SBarry Smith   MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
125f0d39aaaSBarry Smith   PetscInt   bdiagsize;                  /* length of bdiag and ibdiag */
126ace3abfcSBarry Smith   PetscBool  ibdiagvalid;                /* do ibdiag[] and bdiag[] contain the most recent values */
127f0d39aaaSBarry Smith 
128ace3abfcSBarry Smith   PetscBool        use;
129e6b907acSBarry Smith   PetscInt         node_count;       /* number of inodes */
130e6b907acSBarry Smith   PetscInt        *size;             /* size of each inode */
131e6b907acSBarry Smith   PetscInt         limit;            /* inode limit */
132e6b907acSBarry Smith   PetscInt         max_limit;        /* maximum supported inode limit */
133ace3abfcSBarry Smith   PetscBool        checked;          /* if inodes have been checked for */
134a02bda8eSBarry Smith   PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */
1354108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
136e6b907acSBarry Smith 
1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer);
1385a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType);
1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool);
1425a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *);
1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool);
1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat, Mat, const MatFactorInfo *);
1455a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *);
146f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **);
147f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **);
148e6b907acSBarry Smith 
149e6b907acSBarry Smith typedef struct {
15054f21887SBarry Smith   SEQAIJHEADER(MatScalar);
1514108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
15254f21887SBarry Smith   MatScalar       *saved_values; /* location for stashing nonzero values of matrix */
15371f1c65dSBarry Smith 
15471f1c65dSBarry Smith   PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
155ace3abfcSBarry Smith   PetscBool    idiagvalid;                /* current idiag[] and mdiag[] are valid */
156bbead8a2SBarry Smith   PetscScalar *ibdiag;                    /* inverses of block diagonals */
1572291e4fdSJed Brown   PetscBool    ibdiagvalid;               /* inverses of block diagonals are valid. */
15861ecd0c6SBarry Smith   PetscBool    diagonaldense;             /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
15971f1c65dSBarry Smith   PetscScalar  fshift, omega;             /* last used omega and fshift */
160394ed5ebSJunchao Zhang 
161394ed5ebSJunchao Zhang   /* MatSetValuesCOO() related fields on host */
162394ed5ebSJunchao Zhang   PetscCount  coo_n; /* Number of entries in MatSetPreallocationCOO() */
163394ed5ebSJunchao Zhang   PetscCount  Atot;  /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */
164394ed5ebSJunchao Zhang   PetscCount *jmap;  /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */
165394ed5ebSJunchao Zhang   PetscCount *perm;  /* The permutation array in sorting (i,j) by row and then by col */
166ec8511deSBarry Smith } Mat_SeqAIJ;
167b8a66259SBarry Smith 
168e6b907acSBarry Smith /*
169e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
170e6b907acSBarry Smith */
171d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i)
172d71ae5a4SJacob Faibussowitsch {
173e6b907acSBarry Smith   Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data;
174*3ba16761SJacob Faibussowitsch 
175*3ba16761SJacob Faibussowitsch   PetscFunctionBegin;
176e6b907acSBarry Smith   if (A->singlemalloc) {
1779566063dSJacob Faibussowitsch     PetscCall(PetscFree3(*a, *j, *i));
178e6b907acSBarry Smith   } else {
1799566063dSJacob Faibussowitsch     if (A->free_a) PetscCall(PetscFree(*a));
1809566063dSJacob Faibussowitsch     if (A->free_ij) PetscCall(PetscFree(*j));
1819566063dSJacob Faibussowitsch     if (A->free_ij) PetscCall(PetscFree(*i));
182e6b907acSBarry Smith   }
183*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184e6b907acSBarry Smith }
185e6b907acSBarry Smith /*
186e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
187357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
188e6b907acSBarry Smith */
189fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \
190fef13f97SBarry Smith   if (NROW >= RMAX) { \
191fef13f97SBarry Smith     Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
192fef13f97SBarry Smith     /* there is no extra room in row, therefore enlarge */ \
193f4259b30SLisandro Dalcin     PetscInt  CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
194fef13f97SBarry Smith     datatype *new_a; \
195fef13f97SBarry Smith \
19608401ef6SPierre 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); \
197fef13f97SBarry Smith     /* malloc new storage space */ \
1989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(BS2 *new_nz, &new_a, new_nz, &new_j, AM + 1, &new_i)); \
199fef13f97SBarry Smith \
200fef13f97SBarry Smith     /* copy over old data into new slots */ \
201ad540459SPierre Jolivet     for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
202ad540459SPierre Jolivet     for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
2039566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
204fef13f97SBarry Smith     len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
2059566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
2069566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \
2079566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \
2089566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), AA + BS2 * (AI[ROW] + NROW), BS2 * len)); \
209fef13f97SBarry Smith     /* free up old matrix storage */ \
2109566063dSJacob Faibussowitsch     PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
211fef13f97SBarry Smith     AA     = new_a; \
212fef13f97SBarry Smith     Ain->a = (MatScalar *)new_a; \
2139371c9d4SSatish Balay     AI = Ain->i = new_i; \
2149371c9d4SSatish Balay     AJ = Ain->j       = new_j; \
215fef13f97SBarry Smith     Ain->singlemalloc = PETSC_TRUE; \
216fef13f97SBarry Smith \
2179371c9d4SSatish Balay     RP   = AJ + AI[ROW]; \
2189371c9d4SSatish Balay     AP   = AA + BS2 * AI[ROW]; \
219fef13f97SBarry Smith     RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
220fef13f97SBarry Smith     Ain->maxnz += BS2 * CHUNKSIZE; \
221fef13f97SBarry Smith     Ain->reallocs++; \
2229371c9d4SSatish Balay   }
22317454e89SShri Abhyankar 
224876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \
225720833daSHong Zhang   if (NROW >= RMAX) { \
226720833daSHong Zhang     Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
227720833daSHong Zhang     /* there is no extra room in row, therefore enlarge */ \
228f4259b30SLisandro Dalcin     PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
229720833daSHong Zhang \
23008401ef6SPierre 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); \
231720833daSHong Zhang     /* malloc new storage space */ \
2329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(new_nz, &new_j)); \
2339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(AM + 1, &new_i)); \
234720833daSHong Zhang \
235720833daSHong Zhang     /* copy over old data into new slots */ \
236ad540459SPierre Jolivet     for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
237ad540459SPierre Jolivet     for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
2389566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
239720833daSHong Zhang     len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
2409566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
241876c6284SHong Zhang \
242720833daSHong Zhang     /* free up old matrix storage */ \
2439566063dSJacob Faibussowitsch     PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
244876c6284SHong Zhang     Ain->a = NULL; \
2459371c9d4SSatish Balay     AI = Ain->i = new_i; \
2469371c9d4SSatish Balay     AJ = Ain->j       = new_j; \
247720833daSHong Zhang     Ain->singlemalloc = PETSC_FALSE; \
248876c6284SHong Zhang     Ain->free_a       = PETSC_FALSE; \
249720833daSHong Zhang \
250876c6284SHong Zhang     RP   = AJ + AI[ROW]; \
251720833daSHong Zhang     RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
252720833daSHong Zhang     Ain->maxnz += BS2 * CHUNKSIZE; \
253720833daSHong Zhang     Ain->reallocs++; \
2549371c9d4SSatish Balay   }
255e6b907acSBarry Smith 
256cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *);
257e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);
258cbc6b225SStefano Zampini PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_SeqAIJ(Mat);
259cbc6b225SStefano Zampini 
2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *);
2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *);
2631df811f5SHong Zhang 
2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
2705a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *);
2715a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure);
2725a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *);
2735a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2745a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **);
27508480c60SBarry Smith 
276b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec);
277b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec);
278b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec);
279b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec);
280b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec);
281b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
2825a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
283b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
28408480c60SBarry Smith 
2855a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool);
2863a7fca6bSBarry Smith 
2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
2885a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
2897fb60732SBarry Smith PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]);
2905a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *);
2915008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *);
2927fb60732SBarry Smith 
2932462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **);
2945a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *);
2955a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
2985a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *);
2995a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *);
3005a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec);
3015a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec);
3025a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat, Vec, Vec);
3035a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec);
3045a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat, Vec, Vec);
3055a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec);
3065a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat, Vec, Vec);
3075a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
3085a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec);
3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec);
3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec);
3115a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
3125a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
3135a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat, Mat, Mat);
3145a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat);
31579c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *);
31693dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring);
317f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring);
318a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt);
31952f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer);
32052f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer);
3215a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer);
3225a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
3237bab7c10SHong Zhang 
324df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE)
3254222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
326df97dc6dSFande Kong #endif
3274222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);
3284222ddf1SHong Zhang 
3294222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
3304222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
3314222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);
3324222ddf1SHong Zhang 
3334222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
3344222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat);
3354222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat);
3364222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat);
3374222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat);
3384222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat);
3394222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat);
3404222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat);
3414222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat);
3424222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE)
3434222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
3444222ddf1SHong Zhang #endif
3454222ddf1SHong Zhang 
3465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
347df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat);
3484222ddf1SHong Zhang 
3494099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat);
3505a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat);
3512b8ad9a3SHong Zhang 
3524222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat);
3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
3545a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat);
35553565b12SHong Zhang 
3564222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
3574222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat);
3584222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat);
3595a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
36055bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat);
36155bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat);
3625df89d91SHong Zhang 
3634222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
3645a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
3656718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *);
3663bf78175SHong Zhang 
3674222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
3685a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
3695a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring);
3705a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat);
3715a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat);
3725df89d91SHong Zhang 
3734222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat);
3745a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat);
3757bab7c10SHong Zhang 
376679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom);
3775a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
3785a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
3795a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
380db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar);
38187c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec);
38287c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode);
3835a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure);
3845a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
3855a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
3865a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
3875a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
3886378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
3896378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
3905a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
3915a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
3925a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer);
3939af31e4aSHong Zhang 
3945a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
3955a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
396a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
397a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
398019b515eSShri Abhyankar 
3995a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *);
4009f5f6813SShri Abhyankar 
401d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
402388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *);
403388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *);
404f2fbf96bSVaclav Hapla #endif
405cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *);
406cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *);
407388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *);
408388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
409388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
410d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
411d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
412d24d4204SJose E. Roman #endif
413388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *);
414cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *);
4154dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *);
4164a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *);
417388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *);
4185a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS);
4195a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *);
4208cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
4215a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType);
4223fa6b06aSMark Adams PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat);
4232f947c57SVictor Minden 
424b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *);
4259c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
4269c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
42770f19b1fSKris Buschelman 
42853dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat);
429dec0b466SHong Zhang PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A);
430f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *);
4310fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
432f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
43363a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]);
434feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *);
43553dd7562SDmitry Karpeev 
436a3bb6f32SFande Kong PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *);
437e4e71118SStefano Zampini PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat);
438a3bb6f32SFande Kong 
439003131ecSBarry Smith /*
440003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
441003131ecSBarry Smith 
442003131ecSBarry Smith   Input Parameters:
443003131ecSBarry Smith +  nnz - the number of entries
444003131ecSBarry Smith .  r - the array of vector values
445003131ecSBarry Smith .  xv - the matrix values for the row
446003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
447003131ecSBarry Smith 
448003131ecSBarry Smith   Output Parameter:
449003131ecSBarry Smith .  sum - negative the sum of results
450003131ecSBarry Smith 
451003131ecSBarry Smith   PETSc compile flags:
4527b42bb93SJunchao Zhang +   PETSC_KERNEL_USE_UNROLL_4
4537b42bb93SJunchao Zhang -   PETSC_KERNEL_USE_UNROLL_2
4547b42bb93SJunchao Zhang 
45511a5261eSBarry Smith   Developer Note:
4567b42bb93SJunchao Zhang     The macro changes sum but not other parameters
457003131ecSBarry Smith 
458db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()`
459003131ecSBarry Smith */
460519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
4619371c9d4SSatish Balay   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
4629371c9d4SSatish Balay     { \
463003131ecSBarry Smith       if (nnz > 0) { \
4647b42bb93SJunchao Zhang         PetscInt nnz2 = nnz, rem = nnz & 0x3; \
4657b42bb93SJunchao Zhang         switch (rem) { \
466d71ae5a4SJacob Faibussowitsch         case 3: \
467d71ae5a4SJacob Faibussowitsch           sum -= *xv++ * r[*xi++]; \
468d71ae5a4SJacob Faibussowitsch         case 2: \
469d71ae5a4SJacob Faibussowitsch           sum -= *xv++ * r[*xi++]; \
470d71ae5a4SJacob Faibussowitsch         case 1: \
471d71ae5a4SJacob Faibussowitsch           sum -= *xv++ * r[*xi++]; \
472d71ae5a4SJacob Faibussowitsch           nnz2 -= rem; \
4737b42bb93SJunchao Zhang         } \
4749371c9d4SSatish Balay         while (nnz2 > 0) { \
4759371c9d4SSatish Balay           sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
4769371c9d4SSatish Balay           xv += 4; \
4779371c9d4SSatish Balay           xi += 4; \
4789371c9d4SSatish Balay           nnz2 -= 4; \
4799371c9d4SSatish Balay         } \
4809371c9d4SSatish Balay         xv -= nnz; \
4819371c9d4SSatish Balay         xi -= nnz; \
4827b42bb93SJunchao Zhang       } \
4837b42bb93SJunchao Zhang     }
484003131ecSBarry Smith 
485003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
4869371c9d4SSatish Balay   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
4879371c9d4SSatish Balay     { \
488003131ecSBarry Smith       PetscInt __i, __i1, __i2; \
4899371c9d4SSatish Balay       for (__i = 0; __i < nnz - 1; __i += 2) { \
4909371c9d4SSatish Balay         __i1 = xi[__i]; \
4919371c9d4SSatish Balay         __i2 = xi[__i + 1]; \
4929371c9d4SSatish Balay         sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
4939371c9d4SSatish Balay       } \
4949371c9d4SSatish Balay       if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \
4959371c9d4SSatish Balay     }
496003131ecSBarry Smith 
497003131ecSBarry Smith #else
4989371c9d4SSatish Balay   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
4999371c9d4SSatish Balay     { \
500003131ecSBarry Smith       PetscInt __i; \
5019371c9d4SSatish Balay       for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \
5029371c9d4SSatish Balay     }
503003131ecSBarry Smith #endif
504003131ecSBarry Smith 
505003131ecSBarry Smith /*
506003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
507003131ecSBarry Smith 
508003131ecSBarry Smith   Input Parameters:
509003131ecSBarry Smith +  nnz - the number of entries
510003131ecSBarry Smith .  r - the array of vector values
511003131ecSBarry Smith .  xv - the matrix values for the row
512003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
513003131ecSBarry Smith 
514003131ecSBarry Smith   Output Parameter:
515003131ecSBarry Smith .  sum - the sum of results
516003131ecSBarry Smith 
517003131ecSBarry Smith   PETSc compile flags:
5187b42bb93SJunchao Zhang +   PETSC_KERNEL_USE_UNROLL_4
5197b42bb93SJunchao Zhang -   PETSC_KERNEL_USE_UNROLL_2
5207b42bb93SJunchao Zhang 
52111a5261eSBarry Smith   Developer Note:
5227b42bb93SJunchao Zhang     The macro changes sum but not other parameters
523003131ecSBarry Smith 
524db781477SPatrick Sanan .seealso: `PetscSparseDenseMinusDot()`
525003131ecSBarry Smith */
526519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
5279371c9d4SSatish Balay   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
5289371c9d4SSatish Balay     { \
529003131ecSBarry Smith       if (nnz > 0) { \
5307b42bb93SJunchao Zhang         PetscInt nnz2 = nnz, rem = nnz & 0x3; \
5317b42bb93SJunchao Zhang         switch (rem) { \
532d71ae5a4SJacob Faibussowitsch         case 3: \
533d71ae5a4SJacob Faibussowitsch           sum += *xv++ * r[*xi++]; \
534d71ae5a4SJacob Faibussowitsch         case 2: \
535d71ae5a4SJacob Faibussowitsch           sum += *xv++ * r[*xi++]; \
536d71ae5a4SJacob Faibussowitsch         case 1: \
537d71ae5a4SJacob Faibussowitsch           sum += *xv++ * r[*xi++]; \
538d71ae5a4SJacob Faibussowitsch           nnz2 -= rem; \
5397b42bb93SJunchao Zhang         } \
5409371c9d4SSatish Balay         while (nnz2 > 0) { \
5419371c9d4SSatish Balay           sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
5429371c9d4SSatish Balay           xv += 4; \
5439371c9d4SSatish Balay           xi += 4; \
5449371c9d4SSatish Balay           nnz2 -= 4; \
5459371c9d4SSatish Balay         } \
5469371c9d4SSatish Balay         xv -= nnz; \
5479371c9d4SSatish Balay         xi -= nnz; \
5487b42bb93SJunchao Zhang       } \
5497b42bb93SJunchao Zhang     }
550003131ecSBarry Smith 
551003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
5529371c9d4SSatish Balay   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
5539371c9d4SSatish Balay     { \
554003131ecSBarry Smith       PetscInt __i, __i1, __i2; \
5559371c9d4SSatish Balay       for (__i = 0; __i < nnz - 1; __i += 2) { \
5569371c9d4SSatish Balay         __i1 = xi[__i]; \
5579371c9d4SSatish Balay         __i2 = xi[__i + 1]; \
5589371c9d4SSatish Balay         sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
5599371c9d4SSatish Balay       } \
5609371c9d4SSatish Balay       if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \
5619371c9d4SSatish Balay     }
562003131ecSBarry Smith 
56399acd6aaSStefano 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)
56454e8760dSRichard Tran Mills   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz))
565d68ff82bSRichard Tran Mills 
566003131ecSBarry Smith #else
5679371c9d4SSatish Balay   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
5689371c9d4SSatish Balay     { \
569003131ecSBarry Smith       PetscInt __i; \
5709371c9d4SSatish Balay       for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \
5719371c9d4SSatish Balay     }
572003131ecSBarry Smith #endif
573003131ecSBarry Smith 
57499acd6aaSStefano 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)
57554e8760dSRichard Tran Mills   #include <immintrin.h>
57654e8760dSRichard Tran Mills   #if !defined(_MM_SCALE_8)
57754e8760dSRichard Tran Mills     #define _MM_SCALE_8 8
57854e8760dSRichard Tran Mills   #endif
57954e8760dSRichard Tran Mills 
580d71ae5a4SJacob Faibussowitsch static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n)
581d71ae5a4SJacob Faibussowitsch {
58254e8760dSRichard Tran Mills   __m512d  vec_x, vec_y, vec_vals;
58354e8760dSRichard Tran Mills   __m256i  vec_idx;
58454e8760dSRichard Tran Mills   __mmask8 mask;
58554e8760dSRichard Tran Mills   PetscInt j;
58654e8760dSRichard Tran Mills 
58754e8760dSRichard Tran Mills   vec_y = _mm512_setzero_pd();
58854e8760dSRichard Tran Mills   for (j = 0; j < (n >> 3); j++) {
589ef588d5cSRichard Tran Mills     vec_idx  = _mm256_loadu_si256((__m256i const *)aj);
590ef588d5cSRichard Tran Mills     vec_vals = _mm512_loadu_pd(aa);
59154e8760dSRichard Tran Mills     vec_x    = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
59254e8760dSRichard Tran Mills     vec_y    = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
5939371c9d4SSatish Balay     aj += 8;
5949371c9d4SSatish Balay     aa += 8;
59554e8760dSRichard Tran Mills   }
596851da9e1SHong Zhang   #if defined(__AVX512VL__)
597851da9e1SHong Zhang   /* masked load requires avx512vl, which is not supported by KNL */
598851da9e1SHong Zhang   if (n & 0x07) {
59954e8760dSRichard Tran Mills     mask     = (__mmask8)(0xff >> (8 - (n & 0x07)));
600851da9e1SHong Zhang     vec_idx  = _mm256_mask_loadu_epi32(vec_idx, mask, aj);
601851da9e1SHong Zhang     vec_vals = _mm512_mask_loadu_pd(vec_vals, mask, aa);
60254e8760dSRichard Tran Mills     vec_x    = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8);
60354e8760dSRichard Tran Mills     vec_y    = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask);
60454e8760dSRichard Tran Mills   }
605851da9e1SHong Zhang   *sum += _mm512_reduce_add_pd(vec_y);
606851da9e1SHong Zhang   #else
607851da9e1SHong Zhang   *sum += _mm512_reduce_add_pd(vec_y);
60854e8760dSRichard Tran Mills   for (j = 0; j < (n & 0x07); j++) *sum += aa[j] * x[aj[j]];
609851da9e1SHong Zhang   #endif
61054e8760dSRichard Tran Mills }
61154e8760dSRichard Tran Mills #endif
612b434eb95SMatthew G. Knepley 
613b434eb95SMatthew G. Knepley /*
614b434eb95SMatthew G. Knepley     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
615b434eb95SMatthew G. Knepley 
616b434eb95SMatthew G. Knepley   Input Parameters:
617b434eb95SMatthew G. Knepley +  nnz - the number of entries
618b434eb95SMatthew G. Knepley .  r - the array of vector values
619b434eb95SMatthew G. Knepley .  xv - the matrix values for the row
620b434eb95SMatthew G. Knepley -  xi - the column indices of the nonzeros in the row
621b434eb95SMatthew G. Knepley 
622b434eb95SMatthew G. Knepley   Output Parameter:
623b434eb95SMatthew G. Knepley .  max - the max of results
624b434eb95SMatthew G. Knepley 
625db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()`
626b434eb95SMatthew G. Knepley */
6279371c9d4SSatish Balay #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \
628eec179cfSJacob Faibussowitsch   do { \
629eec179cfSJacob Faibussowitsch     for (PetscInt __i = 0; __i < (nnz); __i++) { max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); } \
630eec179cfSJacob Faibussowitsch   } while (0)
631b434eb95SMatthew G. Knepley 
6324b38b95cSHong Zhang /*
6334b38b95cSHong Zhang  Add column indices into table for counting the max nonzeros of merged rows
6344b38b95cSHong Zhang  */
6359371c9d4SSatish Balay #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \
636eec179cfSJacob Faibussowitsch   do { \
637eec179cfSJacob Faibussowitsch     if ((mat)) { \
638eec179cfSJacob Faibussowitsch       for (PetscInt _row = 0; _row < (nrows); _row++) { \
639eec179cfSJacob Faibussowitsch         const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \
640eec179cfSJacob Faibussowitsch         for (PetscInt _j = 0; _j < _nz; _j++) { \
641eec179cfSJacob Faibussowitsch           PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
642c76ffc5fSJacob Faibussowitsch           PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \
6434b38b95cSHong Zhang         } \
6444b38b95cSHong Zhang       } \
645ec07b8f8SHong Zhang     } \
646eec179cfSJacob Faibussowitsch   } while (0)
6474b38b95cSHong Zhang 
6480ca7d551SHong Zhang /*
6490ca7d551SHong Zhang  Add column indices into table for counting the nonzeros of merged rows
6500ca7d551SHong Zhang  */
6519371c9d4SSatish Balay #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \
652eec179cfSJacob Faibussowitsch   do { \
653eec179cfSJacob Faibussowitsch     for (PetscInt _i = 0; _i < (nrows); _i++) { \
654eec179cfSJacob Faibussowitsch       const PetscInt _row = (rows)[_i]; \
655eec179cfSJacob Faibussowitsch       const PetscInt _nz  = (mat)->i[_row + 1] - (mat)->i[_row]; \
656eec179cfSJacob Faibussowitsch       for (PetscInt _j = 0; _j < _nz; _j++) { \
657eec179cfSJacob Faibussowitsch         PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
658eec179cfSJacob Faibussowitsch         PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \
6590ca7d551SHong Zhang       } \
6600ca7d551SHong Zhang     } \
661eec179cfSJacob Faibussowitsch   } while (0)
6620ca7d551SHong Zhang 
6632d40f771SBarry Smith #endif
664