xref: /petsc/src/mat/impls/aij/seq/aij.h (revision ef588d5c511ff21a987e266aaf5c648de461dbe2)
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 */
552cff0574SHong Zhang   PetscErrorCode (*destroy)(Mat);
562cff0574SHong Zhang } Mat_MatTransMatMult;
572cff0574SHong Zhang 
5853565b12SHong Zhang typedef struct {
5953565b12SHong Zhang   PetscInt    *api,*apj;       /* symbolic structure of A*P */
6053565b12SHong Zhang   PetscScalar *apa;            /* temporary array for storing one row of A*P */
6153565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
6253565b12SHong Zhang } Mat_PtAP;
6353565b12SHong Zhang 
6453565b12SHong Zhang typedef struct {
6553565b12SHong Zhang   MatTransposeColoring matcoloring;
66257c235dSHong Zhang   Mat                  Rt;    /* sparse or dense matrix of R^T */
6753565b12SHong Zhang   Mat                  RARt;  /* dense matrix of R*A*R^T */
683b1b9624SHong Zhang   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
6953565b12SHong Zhang   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
7053565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
7153565b12SHong Zhang } Mat_RARt;
7253565b12SHong Zhang 
736d0b6147SHong Zhang typedef struct {
746d0b6147SHong Zhang   Mat BC;               /* temp matrix for storing B*C */
756d0b6147SHong Zhang   PetscErrorCode (*destroy)(Mat);
766d0b6147SHong Zhang } Mat_MatMatMatMult;
776d0b6147SHong Zhang 
782d40f771SBarry Smith /*
79ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
80e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
81dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
825768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
832d40f771SBarry Smith */
84d35516d3SLois Curfman McInnes 
85e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
86b8a66259SBarry Smith typedef struct {
874108e4d5SBarry Smith   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
88f0d39aaaSBarry Smith   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
89ace3abfcSBarry Smith   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
90f0d39aaaSBarry Smith 
91ace3abfcSBarry Smith   PetscBool        use;
92e6b907acSBarry Smith   PetscInt         node_count;                     /* number of inodes */
93e6b907acSBarry Smith   PetscInt         *size;                          /* size of each inode */
94e6b907acSBarry Smith   PetscInt         limit;                          /* inode limit */
95e6b907acSBarry Smith   PetscInt         max_limit;                      /* maximum supported inode limit */
96ace3abfcSBarry Smith   PetscBool        checked;                        /* if inodes have been checked for */
97a02bda8eSBarry Smith   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
984108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
99e6b907acSBarry Smith 
1005a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
1065a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
109e6b907acSBarry Smith 
110e6b907acSBarry Smith typedef struct {
11154f21887SBarry Smith   SEQAIJHEADER(MatScalar);
1124108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
11354f21887SBarry Smith   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
11471f1c65dSBarry Smith 
11571f1c65dSBarry Smith   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
116ace3abfcSBarry Smith   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
117bbead8a2SBarry Smith   PetscScalar *ibdiag;                        /* inverses of block diagonals */
1182291e4fdSJed Brown   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
11971f1c65dSBarry Smith   PetscScalar fshift,omega;                   /* last used omega and fshift */
12071f1c65dSBarry Smith 
1213a7fca6bSBarry Smith   ISColoring  coloring;                       /* set with MatADSetColoring() used by MatADSetValues() */
1220b7e3e3dSHong Zhang 
1230b7e3e3dSHong Zhang   PetscScalar         *matmult_abdense;    /* used by MatMatMult() */
12453565b12SHong Zhang   Mat_PtAP            *ptap;               /* used by MatPtAP() */
1256d0b6147SHong Zhang   Mat_MatMatMatMult   *matmatmatmult;      /* used by MatMatMatMult() */
126257c235dSHong Zhang   Mat_RARt            *rart;               /* used by MatRARt() */
12740192850SHong Zhang   Mat_MatMatTransMult *abt;                /* used by MatMatTransposeMult() */
1286d373c3eSHong Zhang   Mat_MatTransMatMult *atb;                /* used by MatTransposeMatMult() */
129ec8511deSBarry Smith } Mat_SeqAIJ;
130b8a66259SBarry Smith 
131e6b907acSBarry Smith /*
132e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
133e6b907acSBarry Smith */
134d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
135d0f46423SBarry Smith {
136e6b907acSBarry Smith   PetscErrorCode ierr;
137e6b907acSBarry Smith   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
138e6b907acSBarry Smith   if (A->singlemalloc) {
139e6b907acSBarry Smith     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
140e6b907acSBarry Smith   } else {
141c31cb41cSBarry Smith     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
142c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
143c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
144e6b907acSBarry Smith   }
145e6b907acSBarry Smith   return 0;
146e6b907acSBarry Smith }
147e6b907acSBarry Smith /*
148e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
149357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
150e6b907acSBarry Smith */
151fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
152fef13f97SBarry Smith   if (NROW >= RMAX) { \
153fef13f97SBarry Smith     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
154fef13f97SBarry Smith     /* there is no extra room in row, therefore enlarge */ \
155fef13f97SBarry Smith     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
156fef13f97SBarry Smith     datatype *new_a; \
157fef13f97SBarry Smith  \
1585aae435aSMatthew 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); \
159fef13f97SBarry Smith     /* malloc new storage space */ \
160dcca6d9dSJed Brown     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
161fef13f97SBarry Smith  \
162fef13f97SBarry Smith     /* copy over old data into new slots */ \
163fef13f97SBarry Smith     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
164fef13f97SBarry Smith     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
165fef13f97SBarry Smith     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
166fef13f97SBarry Smith     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
167fef13f97SBarry Smith     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
168fef13f97SBarry Smith     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
169fef13f97SBarry Smith     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
170fef13f97SBarry Smith     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
171fef13f97SBarry Smith     /* free up old matrix storage */ \
172fef13f97SBarry Smith     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
173fef13f97SBarry Smith     AA                = new_a; \
174fef13f97SBarry Smith     Ain->a            = (MatScalar*) new_a;                   \
175fef13f97SBarry Smith     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
176fef13f97SBarry Smith     Ain->singlemalloc = PETSC_TRUE; \
177fef13f97SBarry Smith  \
178fef13f97SBarry Smith     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
179fef13f97SBarry Smith     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
180fef13f97SBarry Smith     Ain->maxnz += BS2*CHUNKSIZE; \
181fef13f97SBarry Smith     Ain->reallocs++; \
182fef13f97SBarry Smith   } \
18317454e89SShri Abhyankar 
184876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \
185720833daSHong Zhang   if (NROW >= RMAX) { \
186720833daSHong Zhang     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
187720833daSHong Zhang     /* there is no extra room in row, therefore enlarge */ \
188720833daSHong Zhang     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
189720833daSHong Zhang  \
190720833daSHong 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); \
191720833daSHong Zhang     /* malloc new storage space */ \
192720833daSHong Zhang     ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \
193720833daSHong Zhang     ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\
194720833daSHong Zhang  \
195720833daSHong Zhang     /* copy over old data into new slots */ \
196720833daSHong Zhang     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
197720833daSHong Zhang     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
198720833daSHong Zhang     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
199720833daSHong Zhang     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
200720833daSHong Zhang     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
201876c6284SHong Zhang  \
202720833daSHong Zhang     /* free up old matrix storage */ \
203720833daSHong Zhang     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
204876c6284SHong Zhang     Ain->a            = NULL;                   \
205720833daSHong Zhang     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
206720833daSHong Zhang     Ain->singlemalloc = PETSC_FALSE; \
207876c6284SHong Zhang     Ain->free_a       = PETSC_FALSE;             \
208720833daSHong Zhang  \
209876c6284SHong Zhang     RP          = AJ + AI[ROW];    \
210720833daSHong Zhang     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
211720833daSHong Zhang     Ain->maxnz += BS2*CHUNKSIZE; \
212720833daSHong Zhang     Ain->reallocs++; \
213720833daSHong Zhang   } \
214e6b907acSBarry Smith 
215cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
2165a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2175a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2185a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
2191df811f5SHong Zhang 
2205a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
2225a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
23108480c60SBarry Smith 
2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
23708480c60SBarry Smith 
2385a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
2393a7fca6bSBarry Smith 
2405a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
2425a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2435a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
2445008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
2452462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2485a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2495a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2505a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
2525a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
2535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
2545a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
2555a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
26779c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
26893dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
269f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
270a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
2715a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
2725a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
2737bab7c10SHong Zhang 
2745a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2755a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2764099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
2785a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
2795a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
2805a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
2815a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2824099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
2842b8ad9a3SHong Zhang 
2855a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2885a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
28953565b12SHong Zhang 
2905a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
29155bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
29255bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
2935a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
29455bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
29555bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
2965df89d91SHong Zhang 
2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2985a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2995a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
3006d373c3eSHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat);
3013bf78175SHong Zhang 
3023bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
3033bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
3043bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
3053bf78175SHong Zhang 
3065a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3075a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
3085a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
3115a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
3125df89d91SHong Zhang 
3135a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
3145a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
3155a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
3167bab7c10SHong Zhang 
3175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3185a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
3195a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
320db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
32187c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
32287c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
3235a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
3245a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3255a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3265a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3275a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3286378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
3296378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
3305a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
3315a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
3325a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
3339af31e4aSHong Zhang 
3345a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
3355a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
336a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
337a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
338019b515eSShri Abhyankar 
3395a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
3409f5f6813SShri Abhyankar 
341388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
342388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
343388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
344cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
345cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
346388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
347388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
348388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
349388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
350cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
3514a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
352388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
3545a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3555a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
3575a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
3582f947c57SVictor Minden 
359b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
3609c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
3619c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
36270f19b1fSKris Buschelman 
36353dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
364f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
3650fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
366f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
367feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
36853dd7562SDmitry Karpeev 
369003131ecSBarry Smith /*
370003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
371003131ecSBarry Smith 
372003131ecSBarry Smith   Input Parameters:
373003131ecSBarry Smith +  nnz - the number of entries
374003131ecSBarry Smith .  r - the array of vector values
375003131ecSBarry Smith .  xv - the matrix values for the row
376003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
377003131ecSBarry Smith 
378003131ecSBarry Smith   Output Parameter:
379003131ecSBarry Smith .  sum - negative the sum of results
380003131ecSBarry Smith 
381003131ecSBarry Smith   PETSc compile flags:
382b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
383e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
384003131ecSBarry Smith 
385003131ecSBarry Smith .seealso: PetscSparseDensePlusDot()
386003131ecSBarry Smith 
387003131ecSBarry Smith */
388519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
389003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
390003131ecSBarry Smith     if (nnz > 0) { \
391003131ecSBarry Smith       switch (nnz & 0x3) { \
392003131ecSBarry Smith       case 3: sum -= *xv++ *r[*xi++]; \
393003131ecSBarry Smith       case 2: sum -= *xv++ *r[*xi++]; \
394003131ecSBarry Smith       case 1: sum -= *xv++ *r[*xi++]; \
395003131ecSBarry Smith         nnz       -= 4;} \
396003131ecSBarry Smith       while (nnz > 0) { \
397003131ecSBarry Smith         sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \
398003131ecSBarry Smith                xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \
399003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
400003131ecSBarry Smith 
401003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
402003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
403003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
404003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
405003131ecSBarry Smith                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
406003131ecSBarry Smith     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
407003131ecSBarry Smith 
408003131ecSBarry Smith #else
409003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
410003131ecSBarry Smith     PetscInt __i; \
411003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
412003131ecSBarry Smith #endif
413003131ecSBarry Smith 
414003131ecSBarry Smith 
415003131ecSBarry Smith 
416003131ecSBarry Smith /*
417003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
418003131ecSBarry Smith 
419003131ecSBarry Smith   Input Parameters:
420003131ecSBarry Smith +  nnz - the number of entries
421003131ecSBarry Smith .  r - the array of vector values
422003131ecSBarry Smith .  xv - the matrix values for the row
423003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
424003131ecSBarry Smith 
425003131ecSBarry Smith   Output Parameter:
426003131ecSBarry Smith .  sum - the sum of results
427003131ecSBarry Smith 
428003131ecSBarry Smith   PETSc compile flags:
429b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
430e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
431003131ecSBarry Smith 
432003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot()
433003131ecSBarry Smith 
434003131ecSBarry Smith */
435519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
436003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
437003131ecSBarry Smith     if (nnz > 0) { \
438003131ecSBarry Smith       switch (nnz & 0x3) { \
439003131ecSBarry Smith       case 3: sum += *xv++ *r[*xi++]; \
440003131ecSBarry Smith       case 2: sum += *xv++ *r[*xi++]; \
441003131ecSBarry Smith       case 1: sum += *xv++ *r[*xi++]; \
442003131ecSBarry Smith         nnz       -= 4;} \
443003131ecSBarry Smith       while (nnz > 0) { \
444003131ecSBarry Smith         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
445003131ecSBarry Smith                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
446003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
447003131ecSBarry Smith 
448003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
449003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
450003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
451003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
452003131ecSBarry Smith                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
453003131ecSBarry Smith     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
454003131ecSBarry Smith 
455ef4c6e81SRichard Tran Mills #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)
45654e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))
457d68ff82bSRichard Tran Mills 
458003131ecSBarry Smith #else
459003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
460003131ecSBarry Smith     PetscInt __i; \
461003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
462003131ecSBarry Smith #endif
463003131ecSBarry Smith 
464ef4c6e81SRichard Tran Mills #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)
46554e8760dSRichard Tran Mills   #include <immintrin.h>
46654e8760dSRichard Tran Mills   #if !defined(_MM_SCALE_8)
46754e8760dSRichard Tran Mills   #define _MM_SCALE_8    8
46854e8760dSRichard Tran Mills   #endif
46954e8760dSRichard Tran Mills 
47054e8760dSRichard Tran Mills PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
47154e8760dSRichard Tran Mills {
47254e8760dSRichard Tran Mills   __m512d           vec_x,vec_y,vec_vals;
47354e8760dSRichard Tran Mills   __m256i           vec_idx;
47454e8760dSRichard Tran Mills   __mmask8          mask;
47554e8760dSRichard Tran Mills   PetscInt          j;
47654e8760dSRichard Tran Mills 
47754e8760dSRichard Tran Mills   vec_y = _mm512_setzero_pd();
47854e8760dSRichard Tran Mills   for (j=0; j<(n>>3); j++) {
479*ef588d5cSRichard Tran Mills     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
480*ef588d5cSRichard Tran Mills     vec_vals = _mm512_loadu_pd(aa);
48154e8760dSRichard Tran Mills     vec_x    = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8);
48254e8760dSRichard Tran Mills     vec_y    = _mm512_fmadd_pd(vec_x,vec_vals,vec_y);
48354e8760dSRichard Tran Mills     aj += 8; aa += 8;
48454e8760dSRichard Tran Mills   }
48554e8760dSRichard Tran Mills   /* masked load does not work on KNL, it requires avx512vl */
48654e8760dSRichard Tran Mills   if ((n&0x07)>2) {
48754e8760dSRichard Tran Mills     mask     = (__mmask8)(0xff >> (8-(n&0x07)));
488*ef588d5cSRichard Tran Mills     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
489*ef588d5cSRichard Tran Mills     vec_vals = _mm512_loadu_pd(aa);
49054e8760dSRichard Tran Mills     vec_x    = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8);
49154e8760dSRichard Tran Mills     vec_y    = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask);
49254e8760dSRichard Tran Mills   } else if ((n&0x07)==2) {
49354e8760dSRichard Tran Mills     *sum += aa[0]*x[aj[0]];
49454e8760dSRichard Tran Mills     *sum += aa[1]*x[aj[1]];
49554e8760dSRichard Tran Mills   } else if ((n&0x07)==1) {
49654e8760dSRichard Tran Mills     *sum += aa[0]*x[aj[0]];
49754e8760dSRichard Tran Mills   }
49854e8760dSRichard Tran Mills   if (n>2) *sum += _mm512_reduce_add_pd(vec_y);
49954e8760dSRichard Tran Mills /*
50054e8760dSRichard Tran Mills   for(j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
50154e8760dSRichard Tran Mills */
50254e8760dSRichard Tran Mills }
50354e8760dSRichard Tran Mills #endif
504b434eb95SMatthew G. Knepley 
505b434eb95SMatthew G. Knepley /*
506b434eb95SMatthew G. Knepley     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
507b434eb95SMatthew G. Knepley 
508b434eb95SMatthew G. Knepley   Input Parameters:
509b434eb95SMatthew G. Knepley +  nnz - the number of entries
510b434eb95SMatthew G. Knepley .  r - the array of vector values
511b434eb95SMatthew G. Knepley .  xv - the matrix values for the row
512b434eb95SMatthew G. Knepley -  xi - the column indices of the nonzeros in the row
513b434eb95SMatthew G. Knepley 
514b434eb95SMatthew G. Knepley   Output Parameter:
515b434eb95SMatthew G. Knepley .  max - the max of results
516b434eb95SMatthew G. Knepley 
517b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
518b434eb95SMatthew G. Knepley 
519b434eb95SMatthew G. Knepley */
520b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
521b434eb95SMatthew G. Knepley     PetscInt __i; \
5225e792707SMatthew G. Knepley     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
523b434eb95SMatthew G. Knepley 
5244b38b95cSHong Zhang /*
5254b38b95cSHong Zhang  Add column indices into table for counting the max nonzeros of merged rows
5264b38b95cSHong Zhang  */
5274b38b95cSHong Zhang #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) {       \
5284b38b95cSHong Zhang     PetscInt _j,_row,_nz,*_col;                     \
529ec07b8f8SHong Zhang     if (mat) { \
5304b38b95cSHong Zhang       for (_row=0; _row<nrows; _row++) {   \
5314b38b95cSHong Zhang         _nz = mat->i[_row+1] - mat->i[_row];    \
5324b38b95cSHong Zhang         for (_j=0; _j<_nz; _j++) {               \
5334b38b95cSHong Zhang           _col = _j + mat->j + mat->i[_row];       \
5344b38b95cSHong Zhang           PetscTableAdd(ta,*_col+1,1,INSERT_VALUES);                    \
5354b38b95cSHong Zhang         }                                                               \
5364b38b95cSHong Zhang       }                                                                 \
537ec07b8f8SHong Zhang     }    \
5384b38b95cSHong Zhang }
5394b38b95cSHong Zhang 
5400ca7d551SHong Zhang /*
5410ca7d551SHong Zhang  Add column indices into table for counting the nonzeros of merged rows
5420ca7d551SHong Zhang  */
5430ca7d551SHong Zhang #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) {    \
5440ca7d551SHong Zhang   PetscInt _j,_row,_nz,*_col,_i;                      \
5450ca7d551SHong Zhang     for (_i=0; _i<nrows; _i++) {\
5460ca7d551SHong Zhang       _row = rows[_i]; \
5470ca7d551SHong Zhang       _nz = mat->i[_row+1] - mat->i[_row]; \
5480ca7d551SHong Zhang       for (_j=0; _j<_nz; _j++) {                \
5490ca7d551SHong Zhang         _col = _j + mat->j + mat->i[_row];       \
5500ca7d551SHong Zhang         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
5510ca7d551SHong Zhang       }                                                                 \
5520ca7d551SHong Zhang     }                                                                   \
5530ca7d551SHong Zhang }
5540ca7d551SHong Zhang 
5552d40f771SBarry Smith #endif
556