xref: /petsc/src/mat/impls/aij/seq/aij.h (revision e4e711185261551e8ef604f908c6bef2fd293644)
1b8a66259SBarry Smith 
22d40f771SBarry Smith #if !defined(__AIJ_H)
32d40f771SBarry Smith #define __AIJ_H
4e6b907acSBarry Smith 
5af0996ceSBarry Smith #include <petsc/private/matimpl.h>
6eca6b86aSHong Zhang #include <petscctable.h>
78f690400SShri Abhyankar 
8e6b907acSBarry Smith /*
9e6b907acSBarry Smith     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
10e6b907acSBarry Smith */
11421e10b8SBarry Smith #define SEQAIJHEADER(datatype)        \
12ace3abfcSBarry Smith   PetscBool roworiented;              /* if true, row-oriented input, default */ \
13e6b907acSBarry Smith   PetscInt  nonew;                    /* 1 don't add new nonzeros, -1 generate error on new */ \
1428b2fa4aSMatthew Knepley   PetscInt  nounused;                 /* -1 generate error on unused space */ \
15ace3abfcSBarry Smith   PetscBool singlemalloc;             /* if true a, i, and j have been obtained with one big malloc */ \
16e6b907acSBarry Smith   PetscInt  maxnz;                    /* allocated nonzeros */ \
17e6b907acSBarry Smith   PetscInt  *imax;                    /* maximum space allocated for each row */ \
18e6b907acSBarry Smith   PetscInt  *ilen;                    /* actual length of each row */ \
19846b4da1SFande Kong   PetscInt  *ipre;                    /* space preallocated for each row by user */ \
20ace3abfcSBarry Smith   PetscBool free_imax_ilen;  \
21e6b907acSBarry Smith   PetscInt  reallocs;                 /* number of mallocs done during MatSetValues() \
22e6b907acSBarry Smith                                         as more values are set than were prealloced */\
23e6b907acSBarry Smith   PetscInt          rmax;             /* max nonzeros in any row */ \
24ace3abfcSBarry Smith   PetscBool         keepnonzeropattern;   /* keeps matrix structure same in calls to MatZeroRows()*/ \
25ace3abfcSBarry Smith   PetscBool         ignorezeroentries; \
26ace3abfcSBarry Smith   PetscBool         free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
27ace3abfcSBarry Smith   PetscBool         free_a;           /* free the numerical values when matrix is destroy */ \
28e6b907acSBarry Smith   Mat_CompressedRow compressedrow;    /* use compressed row format */                      \
29e6b907acSBarry Smith   PetscInt          nz;               /* nonzeros */                                       \
30e6b907acSBarry Smith   PetscInt          *i;               /* pointer to beginning of each row */               \
31e6b907acSBarry Smith   PetscInt          *j;               /* column values: j + i[k] - 1 is start of row k */  \
32e6b907acSBarry Smith   PetscInt          *diag;            /* pointers to diagonal elements */                  \
337b083b7cSBarry Smith   PetscInt          nonzerorowcnt;    /* how many rows have nonzero entries */             \
34ace3abfcSBarry Smith   PetscBool         free_diag;         \
35421e10b8SBarry Smith   datatype          *a;               /* nonzero elements */                               \
36e6b907acSBarry Smith   PetscScalar       *solve_work;      /* work space used in MatSolve */                    \
374fd072dbSBarry Smith   IS                row, col, icol;   /* index sets, used for reorderings */ \
38ace3abfcSBarry Smith   PetscBool         pivotinblocks;    /* pivot inside factorization of each diagonal block */ \
3917df9f7cSHong Zhang   Mat               parent;           /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
4017df9f7cSHong Zhang                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */\
415c39f6d9SHong Zhang   Mat_SubSppt       *submatis1         /* used by MatCreateSubMatrices_MPIXAIJ_Local */
42e6b907acSBarry Smith 
4353565b12SHong Zhang typedef struct {
4453565b12SHong Zhang   MatTransposeColoring matcoloring;
4553565b12SHong Zhang   Mat                  Bt_den;       /* dense matrix of B^T */
4653565b12SHong Zhang   Mat                  ABt_den;      /* dense matrix of A*B^T */
4753565b12SHong Zhang   PetscBool            usecoloring;
4853565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
4953565b12SHong Zhang } Mat_MatMatTransMult;
5053565b12SHong Zhang 
516d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */
526d373c3eSHong Zhang   Mat          At;           /* transpose of the first matrix */
532cff0574SHong Zhang   Mat          mA;           /* maij matrix of A */
542cff0574SHong Zhang   Vec          bt,ct;        /* vectors to hold locally transposed arrays of B and C */
554222ddf1SHong Zhang   PetscBool    updateAt;     /* flg to avoid recomputing At in MatProductNumeric_AtB_SeqAIJ_SeqAIJ() */
562cff0574SHong Zhang   PetscErrorCode (*destroy)(Mat);
572cff0574SHong Zhang } Mat_MatTransMatMult;
582cff0574SHong Zhang 
5953565b12SHong Zhang typedef struct {
6053565b12SHong Zhang   PetscInt    *api,*apj;       /* symbolic structure of A*P */
6153565b12SHong Zhang   PetscScalar *apa;            /* temporary array for storing one row of A*P */
6253565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
633cdca5ebSHong Zhang } Mat_AP;
6453565b12SHong Zhang 
6553565b12SHong Zhang typedef struct {
6653565b12SHong Zhang   MatTransposeColoring matcoloring;
67257c235dSHong Zhang   Mat                  Rt;    /* sparse or dense matrix of R^T */
6853565b12SHong Zhang   Mat                  RARt;  /* dense matrix of R*A*R^T */
693b1b9624SHong Zhang   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
7053565b12SHong Zhang   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
7153565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
7253565b12SHong Zhang } Mat_RARt;
7353565b12SHong Zhang 
746d0b6147SHong Zhang typedef struct {
756d0b6147SHong Zhang   Mat BC;               /* temp matrix for storing B*C */
766d0b6147SHong Zhang   PetscErrorCode (*destroy)(Mat);
776d0b6147SHong Zhang } Mat_MatMatMatMult;
786d0b6147SHong Zhang 
792d40f771SBarry Smith /*
80ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
81e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
82dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
835768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
842d40f771SBarry Smith */
85d35516d3SLois Curfman McInnes 
86e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
87b8a66259SBarry Smith typedef struct {
884108e4d5SBarry Smith   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
89f0d39aaaSBarry Smith   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
90ace3abfcSBarry Smith   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
91f0d39aaaSBarry Smith 
92ace3abfcSBarry Smith   PetscBool        use;
93e6b907acSBarry Smith   PetscInt         node_count;                     /* number of inodes */
94e6b907acSBarry Smith   PetscInt         *size;                          /* size of each inode */
95e6b907acSBarry Smith   PetscInt         limit;                          /* inode limit */
96e6b907acSBarry Smith   PetscInt         max_limit;                      /* maximum supported inode limit */
97ace3abfcSBarry Smith   PetscBool        checked;                        /* if inodes have been checked for */
98a02bda8eSBarry Smith   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
994108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
100e6b907acSBarry Smith 
1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
1065a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
1095a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
110f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat,PetscScalar**);
111f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat,PetscScalar**);
112e6b907acSBarry Smith 
113e6b907acSBarry Smith typedef struct {
11454f21887SBarry Smith   SEQAIJHEADER(MatScalar);
1154108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
11654f21887SBarry Smith   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
11771f1c65dSBarry Smith 
11871f1c65dSBarry Smith   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
119ace3abfcSBarry Smith   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
120bbead8a2SBarry Smith   PetscScalar *ibdiag;                        /* inverses of block diagonals */
1212291e4fdSJed Brown   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
12261ecd0c6SBarry Smith   PetscBool   diagonaldense;                  /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
12371f1c65dSBarry Smith   PetscScalar fshift,omega;                   /* last used omega and fshift */
12471f1c65dSBarry Smith 
1253a7fca6bSBarry Smith   ISColoring  coloring;                       /* set with MatADSetColoring() used by MatADSetValues() */
1260b7e3e3dSHong Zhang 
1270b7e3e3dSHong Zhang   PetscScalar         *matmult_abdense;    /* used by MatMatMult() */
1283cdca5ebSHong Zhang   Mat_AP              *ap;                 /* used by MatPtAP() */
1296d0b6147SHong Zhang   Mat_MatMatMatMult   *matmatmatmult;      /* used by MatMatMatMult() */
130257c235dSHong Zhang   Mat_RARt            *rart;               /* used by MatRARt() */
13140192850SHong Zhang   Mat_MatMatTransMult *abt;                /* used by MatMatTransposeMult() */
1326d373c3eSHong Zhang   Mat_MatTransMatMult *atb;                /* used by MatTransposeMatMult() */
133ec8511deSBarry Smith } Mat_SeqAIJ;
134b8a66259SBarry Smith 
135e6b907acSBarry Smith /*
136e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
137e6b907acSBarry Smith */
138d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
139d0f46423SBarry Smith {
140e6b907acSBarry Smith   PetscErrorCode ierr;
141e6b907acSBarry Smith   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
142e6b907acSBarry Smith   if (A->singlemalloc) {
143e6b907acSBarry Smith     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
144e6b907acSBarry Smith   } else {
145c31cb41cSBarry Smith     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
146c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
147c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
148e6b907acSBarry Smith   }
149e6b907acSBarry Smith   return 0;
150e6b907acSBarry Smith }
151e6b907acSBarry Smith /*
152e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
153357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
154e6b907acSBarry Smith */
155fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
156fef13f97SBarry Smith   if (NROW >= RMAX) { \
157fef13f97SBarry Smith     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
158fef13f97SBarry Smith     /* there is no extra room in row, therefore enlarge */ \
159fef13f97SBarry Smith     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
160fef13f97SBarry Smith     datatype *new_a; \
161fef13f97SBarry Smith  \
1625aae435aSMatthew G. Knepley     if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \
163fef13f97SBarry Smith     /* malloc new storage space */ \
164dcca6d9dSJed Brown     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
165fef13f97SBarry Smith  \
166fef13f97SBarry Smith     /* copy over old data into new slots */ \
167fef13f97SBarry Smith     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
168fef13f97SBarry Smith     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
169580bdb30SBarry Smith     ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \
170fef13f97SBarry Smith     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
171580bdb30SBarry Smith     ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \
172580bdb30SBarry Smith     ierr = PetscArraycpy(new_a,AA,BS2*(AI[ROW]+NROW));CHKERRQ(ierr);    \
173580bdb30SBarry Smith     ierr = PetscArrayzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE);CHKERRQ(ierr); \
174580bdb30SBarry Smith     ierr = PetscArraycpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len);CHKERRQ(ierr);  \
175fef13f97SBarry Smith     /* free up old matrix storage */ \
176fef13f97SBarry Smith     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
177fef13f97SBarry Smith     AA                = new_a; \
178fef13f97SBarry Smith     Ain->a            = (MatScalar*) new_a;                   \
179fef13f97SBarry Smith     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
180fef13f97SBarry Smith     Ain->singlemalloc = PETSC_TRUE; \
181fef13f97SBarry Smith  \
182fef13f97SBarry Smith     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
183fef13f97SBarry Smith     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
184fef13f97SBarry Smith     Ain->maxnz += BS2*CHUNKSIZE; \
185fef13f97SBarry Smith     Ain->reallocs++; \
186fef13f97SBarry Smith   } \
18717454e89SShri Abhyankar 
188876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \
189720833daSHong Zhang   if (NROW >= RMAX) { \
190720833daSHong Zhang     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
191720833daSHong Zhang     /* there is no extra room in row, therefore enlarge */ \
192720833daSHong Zhang     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
193720833daSHong Zhang  \
194720833daSHong Zhang     if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \
195720833daSHong Zhang     /* malloc new storage space */ \
196720833daSHong Zhang     ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \
197720833daSHong Zhang     ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\
198720833daSHong Zhang  \
199720833daSHong Zhang     /* copy over old data into new slots */ \
200720833daSHong Zhang     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
201720833daSHong Zhang     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
202580bdb30SBarry Smith     ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \
203720833daSHong Zhang     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
204580bdb30SBarry Smith     ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \
205876c6284SHong Zhang  \
206720833daSHong Zhang     /* free up old matrix storage */ \
207720833daSHong Zhang     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
208876c6284SHong Zhang     Ain->a            = NULL;                   \
209720833daSHong Zhang     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
210720833daSHong Zhang     Ain->singlemalloc = PETSC_FALSE; \
211876c6284SHong Zhang     Ain->free_a       = PETSC_FALSE;             \
212720833daSHong Zhang  \
213876c6284SHong Zhang     RP          = AJ + AI[ROW];    \
214720833daSHong Zhang     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
215720833daSHong Zhang     Ain->maxnz += BS2*CHUNKSIZE; \
216720833daSHong Zhang     Ain->reallocs++; \
217720833daSHong Zhang   } \
218e6b907acSBarry Smith 
219cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
2205a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2225a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
2231df811f5SHong Zhang 
2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
23508480c60SBarry Smith 
2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
2375a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2385a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
2395a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2405a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
24108480c60SBarry Smith 
2425a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
2433a7fca6bSBarry Smith 
2445a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
2485008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
2492462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
2505a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2525a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2535a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2545a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
2555a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
2705a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
27179c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
27293dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
273f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
274a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
27552f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat,PetscViewer);
27652f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat,PetscViewer);
2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
2785a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
2797bab7c10SHong Zhang 
280df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE)
2814222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
282df97dc6dSFande Kong #endif
2834222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);
2844222ddf1SHong Zhang 
2854222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
2864222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
2874222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);
2884222ddf1SHong Zhang 
2894222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
2904222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,PetscReal,Mat);
2914222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat);
2924222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat);
2934222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat);
2944222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat);
2954222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat);
2964222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat);
2974222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat);
2984222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE)
2994222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
3004222ddf1SHong Zhang #endif
3014222ddf1SHong Zhang 
3025a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
303df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,Mat);
3044222ddf1SHong Zhang 
3054099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
3065a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
3072b8ad9a3SHong Zhang 
3084222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat);
3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
31153565b12SHong Zhang 
3124222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
3134222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat);
3144222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat);
3155a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
31655bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
31755bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
3185df89d91SHong Zhang 
3194222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
3205a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
3216d373c3eSHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat);
3223bf78175SHong Zhang 
3234222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
3243bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
3253bf78175SHong Zhang 
3264222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
3275a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
3285a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
3295a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
3305a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
3315df89d91SHong Zhang 
3324222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat);
3335a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
3347bab7c10SHong Zhang 
335679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat,PetscInt,PetscInt,PetscRandom);
3365a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3375a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
3385a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
339db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
34087c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
34187c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
3425a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
3435a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3445a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3455a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3465a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3476378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
3486378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
3495a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
3505a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
3515a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
3529af31e4aSHong Zhang 
3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
3545a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
355a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
356a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
357019b515eSShri Abhyankar 
3585a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
3599f5f6813SShri Abhyankar 
360f2fbf96bSVaclav Hapla #if defined(PETSC_HAVE_MATLAB_ENGINE)
361388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
362388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
363f2fbf96bSVaclav Hapla #endif
364cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
365cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
366388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
367388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
368388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
369388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
370cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
3714dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*);
3724a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
373388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
3745a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
3755a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
3775a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
3782f947c57SVictor Minden 
379b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
3809c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
3819c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
38270f19b1fSKris Buschelman 
38353dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
384f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
3850fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
386f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
38763a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]);
388feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
38953dd7562SDmitry Karpeev 
390a3bb6f32SFande Kong PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat,ISLocalToGlobalMapping*);
391*e4e71118SStefano Zampini PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],MatType,Mat);
392a3bb6f32SFande Kong 
393003131ecSBarry Smith /*
394003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
395003131ecSBarry Smith 
396003131ecSBarry Smith   Input Parameters:
397003131ecSBarry Smith +  nnz - the number of entries
398003131ecSBarry Smith .  r - the array of vector values
399003131ecSBarry Smith .  xv - the matrix values for the row
400003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
401003131ecSBarry Smith 
402003131ecSBarry Smith   Output Parameter:
403003131ecSBarry Smith .  sum - negative the sum of results
404003131ecSBarry Smith 
405003131ecSBarry Smith   PETSc compile flags:
4067b42bb93SJunchao Zhang +   PETSC_KERNEL_USE_UNROLL_4
4077b42bb93SJunchao Zhang -   PETSC_KERNEL_USE_UNROLL_2
4087b42bb93SJunchao Zhang 
4097b42bb93SJunchao Zhang   Developer Notes:
4107b42bb93SJunchao Zhang     The macro changes sum but not other parameters
411003131ecSBarry Smith 
412003131ecSBarry Smith .seealso: PetscSparseDensePlusDot()
413003131ecSBarry Smith 
414003131ecSBarry Smith */
415519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
416003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
417003131ecSBarry Smith     if (nnz > 0) { \
4187b42bb93SJunchao Zhang       PetscInt nnz2=nnz,rem=nnz&0x3; \
4197b42bb93SJunchao Zhang       switch (rem) { \
420003131ecSBarry Smith       case 3: sum -= *xv++ *r[*xi++]; \
421003131ecSBarry Smith       case 2: sum -= *xv++ *r[*xi++]; \
422003131ecSBarry Smith       case 1: sum -= *xv++ *r[*xi++]; \
4237b42bb93SJunchao Zhang         nnz2      -= rem;} \
4247b42bb93SJunchao Zhang       while (nnz2 > 0) { \
4257b42bb93SJunchao Zhang         sum -=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
4267b42bb93SJunchao Zhang                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
4277b42bb93SJunchao Zhang         xv += 4; xi += 4; nnz2 -= 4; \
4287b42bb93SJunchao Zhang       } \
4297b42bb93SJunchao Zhang       xv -= nnz; xi -= nnz; \
4307b42bb93SJunchao Zhang     } \
4317b42bb93SJunchao Zhang   }
432003131ecSBarry Smith 
433003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
434003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
435003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
436003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
437003131ecSBarry Smith                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
438003131ecSBarry Smith     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
439003131ecSBarry Smith 
440003131ecSBarry Smith #else
441003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
442003131ecSBarry Smith     PetscInt __i; \
443003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
444003131ecSBarry Smith #endif
445003131ecSBarry Smith 
446003131ecSBarry Smith 
447003131ecSBarry Smith 
448003131ecSBarry Smith /*
449003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
450003131ecSBarry Smith 
451003131ecSBarry Smith   Input Parameters:
452003131ecSBarry Smith +  nnz - the number of entries
453003131ecSBarry Smith .  r - the array of vector values
454003131ecSBarry Smith .  xv - the matrix values for the row
455003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
456003131ecSBarry Smith 
457003131ecSBarry Smith   Output Parameter:
458003131ecSBarry Smith .  sum - the sum of results
459003131ecSBarry Smith 
460003131ecSBarry Smith   PETSc compile flags:
4617b42bb93SJunchao Zhang +   PETSC_KERNEL_USE_UNROLL_4
4627b42bb93SJunchao Zhang -   PETSC_KERNEL_USE_UNROLL_2
4637b42bb93SJunchao Zhang 
4647b42bb93SJunchao Zhang   Developer Notes:
4657b42bb93SJunchao Zhang     The macro changes sum but not other parameters
466003131ecSBarry Smith 
467003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot()
468003131ecSBarry Smith 
469003131ecSBarry Smith */
470519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
471003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
472003131ecSBarry Smith     if (nnz > 0) { \
4737b42bb93SJunchao Zhang       PetscInt nnz2=nnz,rem=nnz&0x3; \
4747b42bb93SJunchao Zhang       switch (rem) { \
475003131ecSBarry Smith       case 3: sum += *xv++ *r[*xi++]; \
476003131ecSBarry Smith       case 2: sum += *xv++ *r[*xi++]; \
477003131ecSBarry Smith       case 1: sum += *xv++ *r[*xi++]; \
4787b42bb93SJunchao Zhang         nnz2      -= rem;} \
4797b42bb93SJunchao Zhang       while (nnz2 > 0) { \
480003131ecSBarry Smith         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
481003131ecSBarry Smith                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
4827b42bb93SJunchao Zhang         xv += 4; xi += 4; nnz2 -= 4; \
4837b42bb93SJunchao Zhang       } \
4847b42bb93SJunchao Zhang       xv -= nnz; xi -= nnz; \
4857b42bb93SJunchao Zhang     } \
4867b42bb93SJunchao Zhang   }
487003131ecSBarry Smith 
488003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
489003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
490003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
491003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
492003131ecSBarry Smith                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
493003131ecSBarry Smith     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
494003131ecSBarry Smith 
49599acd6aaSStefano 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)
49654e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))
497d68ff82bSRichard Tran Mills 
498003131ecSBarry Smith #else
499003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
500003131ecSBarry Smith     PetscInt __i; \
501003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
502003131ecSBarry Smith #endif
503003131ecSBarry Smith 
50499acd6aaSStefano 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)
50554e8760dSRichard Tran Mills   #include <immintrin.h>
50654e8760dSRichard Tran Mills   #if !defined(_MM_SCALE_8)
50754e8760dSRichard Tran Mills   #define _MM_SCALE_8    8
50854e8760dSRichard Tran Mills   #endif
50954e8760dSRichard Tran Mills 
51054e8760dSRichard Tran Mills PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
51154e8760dSRichard Tran Mills {
51254e8760dSRichard Tran Mills   __m512d  vec_x,vec_y,vec_vals;
51354e8760dSRichard Tran Mills   __m256i  vec_idx;
51454e8760dSRichard Tran Mills   __mmask8 mask;
51554e8760dSRichard Tran Mills   PetscInt j;
51654e8760dSRichard Tran Mills 
51754e8760dSRichard Tran Mills   vec_y = _mm512_setzero_pd();
51854e8760dSRichard Tran Mills   for (j=0; j<(n>>3); j++) {
519ef588d5cSRichard Tran Mills     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
520ef588d5cSRichard Tran Mills     vec_vals = _mm512_loadu_pd(aa);
52154e8760dSRichard Tran Mills     vec_x    = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8);
52254e8760dSRichard Tran Mills     vec_y    = _mm512_fmadd_pd(vec_x,vec_vals,vec_y);
52354e8760dSRichard Tran Mills     aj += 8; aa += 8;
52454e8760dSRichard Tran Mills   }
52554e8760dSRichard Tran Mills   /* masked load does not work on KNL, it requires avx512vl */
52654e8760dSRichard Tran Mills   if ((n&0x07)>2) {
52754e8760dSRichard Tran Mills     mask     = (__mmask8)(0xff >> (8-(n&0x07)));
528ef588d5cSRichard Tran Mills     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
529ef588d5cSRichard Tran Mills     vec_vals = _mm512_loadu_pd(aa);
53054e8760dSRichard Tran Mills     vec_x    = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8);
53154e8760dSRichard Tran Mills     vec_y    = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask);
53254e8760dSRichard Tran Mills   } else if ((n&0x07)==2) {
53354e8760dSRichard Tran Mills     *sum += aa[0]*x[aj[0]];
53454e8760dSRichard Tran Mills     *sum += aa[1]*x[aj[1]];
53554e8760dSRichard Tran Mills   } else if ((n&0x07)==1) {
53654e8760dSRichard Tran Mills     *sum += aa[0]*x[aj[0]];
53754e8760dSRichard Tran Mills   }
53854e8760dSRichard Tran Mills   if (n>2) *sum += _mm512_reduce_add_pd(vec_y);
53954e8760dSRichard Tran Mills /*
54054e8760dSRichard Tran Mills   for(j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
54154e8760dSRichard Tran Mills */
54254e8760dSRichard Tran Mills }
54354e8760dSRichard Tran Mills #endif
544b434eb95SMatthew G. Knepley 
545b434eb95SMatthew G. Knepley /*
546b434eb95SMatthew G. Knepley     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
547b434eb95SMatthew G. Knepley 
548b434eb95SMatthew G. Knepley   Input Parameters:
549b434eb95SMatthew G. Knepley +  nnz - the number of entries
550b434eb95SMatthew G. Knepley .  r - the array of vector values
551b434eb95SMatthew G. Knepley .  xv - the matrix values for the row
552b434eb95SMatthew G. Knepley -  xi - the column indices of the nonzeros in the row
553b434eb95SMatthew G. Knepley 
554b434eb95SMatthew G. Knepley   Output Parameter:
555b434eb95SMatthew G. Knepley .  max - the max of results
556b434eb95SMatthew G. Knepley 
557b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
558b434eb95SMatthew G. Knepley 
559b434eb95SMatthew G. Knepley */
560b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
561b434eb95SMatthew G. Knepley     PetscInt __i; \
5625e792707SMatthew G. Knepley     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
563b434eb95SMatthew G. Knepley 
5644b38b95cSHong Zhang /*
5654b38b95cSHong Zhang  Add column indices into table for counting the max nonzeros of merged rows
5664b38b95cSHong Zhang  */
5674b38b95cSHong Zhang #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) {       \
5684b38b95cSHong Zhang     PetscInt _j,_row,_nz,*_col;                     \
569ec07b8f8SHong Zhang     if (mat) { \
5704b38b95cSHong Zhang       for (_row=0; _row<nrows; _row++) {   \
5714b38b95cSHong Zhang         _nz = mat->i[_row+1] - mat->i[_row];    \
5724b38b95cSHong Zhang         for (_j=0; _j<_nz; _j++) {               \
5734b38b95cSHong Zhang           _col = _j + mat->j + mat->i[_row];       \
5744b38b95cSHong Zhang           PetscTableAdd(ta,*_col+1,1,INSERT_VALUES);                    \
5754b38b95cSHong Zhang         }                                                               \
5764b38b95cSHong Zhang       }                                                                 \
577ec07b8f8SHong Zhang     }    \
5784b38b95cSHong Zhang }
5794b38b95cSHong Zhang 
5800ca7d551SHong Zhang /*
5810ca7d551SHong Zhang  Add column indices into table for counting the nonzeros of merged rows
5820ca7d551SHong Zhang  */
5830ca7d551SHong Zhang #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) {    \
5840ca7d551SHong Zhang   PetscInt _j,_row,_nz,*_col,_i;                      \
5850ca7d551SHong Zhang     for (_i=0; _i<nrows; _i++) {\
5860ca7d551SHong Zhang       _row = rows[_i]; \
5870ca7d551SHong Zhang       _nz = mat->i[_row+1] - mat->i[_row]; \
5880ca7d551SHong Zhang       for (_j=0; _j<_nz; _j++) {                \
5890ca7d551SHong Zhang         _col = _j + mat->j + mat->i[_row];       \
5900ca7d551SHong Zhang         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
5910ca7d551SHong Zhang       }                                                                 \
5920ca7d551SHong Zhang     }                                                                   \
5930ca7d551SHong Zhang }
5940ca7d551SHong Zhang 
5952d40f771SBarry Smith #endif
596