xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 5fedec97950c19de564efaecd0f125b1a6cb2b20)
1 
2 #if !defined(__AIJ_H)
3 #define __AIJ_H
4 
5 #include <petsc/private/matimpl.h>
6 #include <petscctable.h>
7 
8 /* Operations provided by MATSEQAIJ and its subclasses */
9 typedef struct {
10   PetscErrorCode (*getarray)(Mat,PetscScalar **);
11   PetscErrorCode (*restorearray)(Mat,PetscScalar **);
12   PetscErrorCode (*getarrayread)(Mat,const PetscScalar **);
13   PetscErrorCode (*restorearrayread)(Mat,const PetscScalar **);
14   PetscErrorCode (*getarraywrite)(Mat,PetscScalar **);
15   PetscErrorCode (*restorearraywrite)(Mat,PetscScalar **);
16 } Mat_SeqAIJOps;
17 
18 /*
19     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
20 */
21 #define SEQAIJHEADER(datatype)        \
22   PetscBool roworiented;              /* if true, row-oriented input, default */ \
23   PetscInt  nonew;                    /* 1 don't add new nonzeros, -1 generate error on new */ \
24   PetscInt  nounused;                 /* -1 generate error on unused space */ \
25   PetscBool singlemalloc;             /* if true a, i, and j have been obtained with one big malloc */ \
26   PetscInt  maxnz;                    /* allocated nonzeros */ \
27   PetscInt  *imax;                    /* maximum space allocated for each row */ \
28   PetscInt  *ilen;                    /* actual length of each row */ \
29   PetscInt  *ipre;                    /* space preallocated for each row by user */ \
30   PetscBool free_imax_ilen;  \
31   PetscInt  reallocs;                 /* number of mallocs done during MatSetValues() \
32                                         as more values are set than were prealloced */\
33   PetscInt          rmax;             /* max nonzeros in any row */ \
34   PetscBool         keepnonzeropattern;   /* keeps matrix structure same in calls to MatZeroRows()*/ \
35   PetscBool         ignorezeroentries; \
36   PetscBool         free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
37   PetscBool         free_a;           /* free the numerical values when matrix is destroy */ \
38   Mat_CompressedRow compressedrow;    /* use compressed row format */                      \
39   PetscInt          nz;               /* nonzeros */                                       \
40   PetscInt          *i;               /* pointer to beginning of each row */               \
41   PetscInt          *j;               /* column values: j + i[k] - 1 is start of row k */  \
42   PetscInt          *diag;            /* pointers to diagonal elements */                  \
43   PetscInt          nonzerorowcnt;    /* how many rows have nonzero entries */             \
44   PetscBool         free_diag;         \
45   datatype          *a;               /* nonzero elements */                               \
46   PetscScalar       *solve_work;      /* work space used in MatSolve */                    \
47   IS                row, col, icol;   /* index sets, used for reorderings */ \
48   PetscBool         pivotinblocks;    /* pivot inside factorization of each diagonal block */ \
49   Mat               parent;           /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
50                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */\
51   Mat_SubSppt       *submatis1;       /* used by MatCreateSubMatrices_MPIXAIJ_Local */    \
52   Mat_SeqAIJOps     ops[1]            /* operations for SeqAIJ and its subclasses */
53 
54 typedef struct {
55   MatTransposeColoring matcoloring;
56   Mat                  Bt_den;       /* dense matrix of B^T */
57   Mat                  ABt_den;      /* dense matrix of A*B^T */
58   PetscBool            usecoloring;
59 } Mat_MatMatTransMult;
60 
61 typedef struct { /* used by MatTransposeMatMult() */
62   Mat          At;           /* transpose of the first matrix */
63   Mat          mA;           /* maij matrix of A */
64   Vec          bt,ct;        /* vectors to hold locally transposed arrays of B and C */
65   PetscBool    updateAt;     /* flg to avoid recomputing At in MatProductNumeric_AtB_SeqAIJ_SeqAIJ() */
66   /* used by PtAP */
67   void           *data;
68   PetscErrorCode (*destroy)(void*);
69 } Mat_MatTransMatMult;
70 
71 typedef struct {
72   PetscInt    *api,*apj;       /* symbolic structure of A*P */
73   PetscScalar *apa;            /* temporary array for storing one row of A*P */
74 } Mat_AP;
75 
76 typedef struct {
77   MatTransposeColoring matcoloring;
78   Mat                  Rt;    /* sparse or dense matrix of R^T */
79   Mat                  RARt;  /* dense matrix of R*A*R^T */
80   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
81   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
82   /* free intermediate products needed for PtAP */
83   void                 *data;
84   PetscErrorCode       (*destroy)(void*);
85 } Mat_RARt;
86 
87 typedef struct {
88   Mat BC;               /* temp matrix for storing B*C */
89 } Mat_MatMatMatMult;
90 
91 /*
92   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
93   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
94   j[i[k]+p] is the pth column in row k.  Note that the diagonal
95   matrix elements are stored with the rest of the nonzeros (not separately).
96 */
97 
98 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
99 typedef struct {
100   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
101   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
102   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
103 
104   PetscBool        use;
105   PetscInt         node_count;                     /* number of inodes */
106   PetscInt         *size;                          /* size of each inode */
107   PetscInt         limit;                          /* inode limit */
108   PetscInt         max_limit;                      /* maximum supported inode limit */
109   PetscBool        checked;                        /* if inodes have been checked for */
110   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
111 } Mat_SeqAIJ_Inode;
112 
113 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
114 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
115 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
116 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
117 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
118 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
119 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
120 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
121 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
122 PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat,PetscScalar**);
123 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat,PetscScalar**);
124 
125 typedef struct {
126   SEQAIJHEADER(MatScalar);
127   Mat_SeqAIJ_Inode inode;
128   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
129 
130   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
131   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
132   PetscScalar *ibdiag;                        /* inverses of block diagonals */
133   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
134   PetscBool   diagonaldense;                  /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
135   PetscScalar fshift,omega;                   /* last used omega and fshift */
136 } Mat_SeqAIJ;
137 
138 /*
139   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
140 */
141 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
142 {
143   PetscErrorCode ierr;
144   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
145   if (A->singlemalloc) {
146     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
147   } else {
148     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
149     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
150     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
151   }
152   return 0;
153 }
154 /*
155     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
156     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
157 */
158 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
159   if (NROW >= RMAX) { \
160     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
161     /* there is no extra room in row, therefore enlarge */ \
162     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=NULL,*new_j=NULL; \
163     datatype *new_a; \
164  \
165     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); \
166     /* malloc new storage space */ \
167     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
168  \
169     /* copy over old data into new slots */ \
170     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
171     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
172     ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \
173     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
174     ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \
175     ierr = PetscArraycpy(new_a,AA,BS2*(AI[ROW]+NROW));CHKERRQ(ierr);    \
176     ierr = PetscArrayzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE);CHKERRQ(ierr); \
177     ierr = PetscArraycpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len);CHKERRQ(ierr);  \
178     /* free up old matrix storage */ \
179     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
180     AA                = new_a; \
181     Ain->a            = (MatScalar*) new_a;                   \
182     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
183     Ain->singlemalloc = PETSC_TRUE; \
184  \
185     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
186     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
187     Ain->maxnz += BS2*CHUNKSIZE; \
188     Ain->reallocs++; \
189   } \
190 
191 #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \
192   if (NROW >= RMAX) { \
193     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
194     /* there is no extra room in row, therefore enlarge */ \
195     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=NULL,*new_j=NULL; \
196  \
197     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); \
198     /* malloc new storage space */ \
199     ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \
200     ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\
201  \
202     /* copy over old data into new slots */ \
203     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
204     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
205     ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \
206     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
207     ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \
208  \
209     /* free up old matrix storage */ \
210     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
211     Ain->a            = NULL;                   \
212     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
213     Ain->singlemalloc = PETSC_FALSE; \
214     Ain->free_a       = PETSC_FALSE;             \
215  \
216     RP          = AJ + AI[ROW];    \
217     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
218     Ain->maxnz += BS2*CHUNKSIZE; \
219     Ain->reallocs++; \
220   } \
221 
222 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
223 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
224 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
225 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
226 
227 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
228 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
229 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
230 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
231 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
232 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
233 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
234 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
235 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
236 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
237 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
238 
239 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat,Vec,Vec);
240 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat,Vec,Vec);
241 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat,Vec,Vec,Vec);
242 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat,Vec,Vec,Vec);
243 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat,Vec,Vec);
244 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
245 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
246 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
247 
248 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
249 
250 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
251 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
252 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
253 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
254 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
255 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
256 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
257 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
258 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
259 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
260 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
261 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
262 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
263 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
264 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
265 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
266 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
267 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
268 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
269 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
270 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
271 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
272 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
273 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
274 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
275 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
276 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
277 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
278 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
279 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
280 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
281 PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat,PetscViewer);
282 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat,PetscViewer);
283 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
284 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
285 
286 #if defined(PETSC_HAVE_HYPRE)
287 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
288 #endif
289 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);
290 
291 PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
292 PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
293 PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);
294 
295 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
296 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,PetscReal,Mat);
297 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat);
298 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat);
299 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat);
300 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat);
301 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat);
302 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat);
303 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat);
304 #if defined(PETSC_HAVE_HYPRE)
305 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
306 #endif
307 
308 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
309 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,Mat);
310 
311 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
312 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
313 
314 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat);
315 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
316 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
317 
318 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
319 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat);
320 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat);
321 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
322 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
323 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
324 
325 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
326 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
327 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void*);
328 
329 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
330 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
331 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
332 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
333 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
334 
335 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat);
336 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
337 
338 PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat,PetscInt,PetscInt,PetscRandom);
339 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
340 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
341 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
342 PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
343 PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
344 PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
345 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
346 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
347 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
348 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
349 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
350 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
351 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
352 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
353 PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
354 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
355 
356 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
357 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
358 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
359 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
360 
361 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
362 
363 #if defined(PETSC_HAVE_MATLAB_ENGINE)
364 PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
365 PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
366 #endif
367 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
368 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
369 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
370 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
371 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
372 #if defined(PETSC_HAVE_SCALAPACK)
373 PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
374 #endif
375 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
376 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
377 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*);
378 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
379 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
380 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
381 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
382 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
383 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
384 PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat);
385 
386 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
387 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
388 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
389 
390 PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
391 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
392 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
393 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
394 PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]);
395 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
396 
397 PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat,ISLocalToGlobalMapping*);
398 PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],MatType,Mat);
399 
400 /*
401     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
402 
403   Input Parameters:
404 +  nnz - the number of entries
405 .  r - the array of vector values
406 .  xv - the matrix values for the row
407 -  xi - the column indices of the nonzeros in the row
408 
409   Output Parameter:
410 .  sum - negative the sum of results
411 
412   PETSc compile flags:
413 +   PETSC_KERNEL_USE_UNROLL_4
414 -   PETSC_KERNEL_USE_UNROLL_2
415 
416   Developer Notes:
417     The macro changes sum but not other parameters
418 
419 .seealso: PetscSparseDensePlusDot()
420 
421 */
422 #if defined(PETSC_KERNEL_USE_UNROLL_4)
423 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
424     if (nnz > 0) { \
425       PetscInt nnz2=nnz,rem=nnz&0x3; \
426       switch (rem) { \
427       case 3: sum -= *xv++ *r[*xi++]; \
428       case 2: sum -= *xv++ *r[*xi++]; \
429       case 1: sum -= *xv++ *r[*xi++]; \
430         nnz2      -= rem;} \
431       while (nnz2 > 0) { \
432         sum -=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
433                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
434         xv += 4; xi += 4; nnz2 -= 4; \
435       } \
436       xv -= nnz; xi -= nnz; \
437     } \
438   }
439 
440 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
441 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
442     PetscInt __i,__i1,__i2; \
443     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
444                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
445     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
446 
447 #else
448 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
449     PetscInt __i; \
450     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
451 #endif
452 
453 /*
454     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
455 
456   Input Parameters:
457 +  nnz - the number of entries
458 .  r - the array of vector values
459 .  xv - the matrix values for the row
460 -  xi - the column indices of the nonzeros in the row
461 
462   Output Parameter:
463 .  sum - the sum of results
464 
465   PETSc compile flags:
466 +   PETSC_KERNEL_USE_UNROLL_4
467 -   PETSC_KERNEL_USE_UNROLL_2
468 
469   Developer Notes:
470     The macro changes sum but not other parameters
471 
472 .seealso: PetscSparseDenseMinusDot()
473 
474 */
475 #if defined(PETSC_KERNEL_USE_UNROLL_4)
476 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
477     if (nnz > 0) { \
478       PetscInt nnz2=nnz,rem=nnz&0x3; \
479       switch (rem) { \
480       case 3: sum += *xv++ *r[*xi++]; \
481       case 2: sum += *xv++ *r[*xi++]; \
482       case 1: sum += *xv++ *r[*xi++]; \
483         nnz2      -= rem;} \
484       while (nnz2 > 0) { \
485         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
486                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
487         xv += 4; xi += 4; nnz2 -= 4; \
488       } \
489       xv -= nnz; xi -= nnz; \
490     } \
491   }
492 
493 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
494 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
495     PetscInt __i,__i1,__i2; \
496     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
497                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
498     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
499 
500 #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)
501 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))
502 
503 #else
504 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
505     PetscInt __i; \
506     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
507 #endif
508 
509 #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)
510   #include <immintrin.h>
511   #if !defined(_MM_SCALE_8)
512   #define _MM_SCALE_8    8
513   #endif
514 
515 PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
516 {
517   __m512d  vec_x,vec_y,vec_vals;
518   __m256i  vec_idx;
519   __mmask8 mask;
520   PetscInt j;
521 
522   vec_y = _mm512_setzero_pd();
523   for (j=0; j<(n>>3); j++) {
524     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
525     vec_vals = _mm512_loadu_pd(aa);
526     vec_x    = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8);
527     vec_y    = _mm512_fmadd_pd(vec_x,vec_vals,vec_y);
528     aj += 8; aa += 8;
529   }
530   /* masked load does not work on KNL, it requires avx512vl */
531   if ((n&0x07)>2) {
532     mask     = (__mmask8)(0xff >> (8-(n&0x07)));
533     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
534     vec_vals = _mm512_loadu_pd(aa);
535     vec_x    = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8);
536     vec_y    = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask);
537   } else if ((n&0x07)==2) {
538     *sum += aa[0]*x[aj[0]];
539     *sum += aa[1]*x[aj[1]];
540   } else if ((n&0x07)==1) {
541     *sum += aa[0]*x[aj[0]];
542   }
543   if (n>2) *sum += _mm512_reduce_add_pd(vec_y);
544 /*
545   for (j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
546 */
547 }
548 #endif
549 
550 /*
551     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
552 
553   Input Parameters:
554 +  nnz - the number of entries
555 .  r - the array of vector values
556 .  xv - the matrix values for the row
557 -  xi - the column indices of the nonzeros in the row
558 
559   Output Parameter:
560 .  max - the max of results
561 
562 .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
563 
564 */
565 #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
566     PetscInt __i; \
567     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
568 
569 /*
570  Add column indices into table for counting the max nonzeros of merged rows
571  */
572 #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) {       \
573     PetscInt _j,_row,_nz,*_col;                     \
574     if (mat) { \
575       for (_row=0; _row<nrows; _row++) {   \
576         _nz = mat->i[_row+1] - mat->i[_row];    \
577         for (_j=0; _j<_nz; _j++) {               \
578           _col = _j + mat->j + mat->i[_row];       \
579           PetscTableAdd(ta,*_col+1,1,INSERT_VALUES);                    \
580         }                                                               \
581       }                                                                 \
582     }    \
583 }
584 
585 /*
586  Add column indices into table for counting the nonzeros of merged rows
587  */
588 #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) {    \
589   PetscInt _j,_row,_nz,*_col,_i;                      \
590     for (_i=0; _i<nrows; _i++) {\
591       _row = rows[_i]; \
592       _nz = mat->i[_row+1] - mat->i[_row]; \
593       for (_j=0; _j<_nz; _j++) {                \
594         _col = _j + mat->j + mat->i[_row];       \
595         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
596       }                                                                 \
597     }                                                                   \
598 }
599 
600 #endif
601