xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 2d30e087755efd99e28fdfe792ffbeb2ee1ea928)
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   PetscErrorCode (*getcsrandmemtype)(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *);
17 } Mat_SeqAIJOps;
18 
19 /*
20     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
21 */
22 #define SEQAIJHEADER(datatype) \
23   PetscBool         roworiented;  /* if true, row-oriented input, default */ \
24   PetscInt          nonew;        /* 1 don't add new nonzeros, -1 generate error on new */ \
25   PetscInt          nounused;     /* -1 generate error on unused space */ \
26   PetscBool         singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \
27   PetscInt          maxnz;        /* allocated nonzeros */ \
28   PetscInt         *imax;         /* maximum space allocated for each row */ \
29   PetscInt         *ilen;         /* actual length of each row */ \
30   PetscInt         *ipre;         /* space preallocated for each row by user */ \
31   PetscBool         free_imax_ilen; \
32   PetscInt          reallocs;           /* number of mallocs done during MatSetValues() \
33                                         as more values are set than were prealloced */ \
34   PetscInt          rmax;               /* max nonzeros in any row */ \
35   PetscBool         keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \
36   PetscBool         ignorezeroentries; \
37   PetscBool         free_ij;       /* free the column indices j and row offsets i when the matrix is destroyed */ \
38   PetscBool         free_a;        /* free the numerical values when matrix is destroy */ \
39   Mat_CompressedRow compressedrow; /* use compressed row format */ \
40   PetscInt          nz;            /* nonzeros */ \
41   PetscInt         *i;             /* pointer to beginning of each row */ \
42   PetscInt         *j;             /* column values: j + i[k] - 1 is start of row k */ \
43   PetscInt         *diag;          /* pointers to diagonal elements */ \
44   PetscInt          nonzerorowcnt; /* how many rows have nonzero entries */ \
45   PetscBool         free_diag; \
46   datatype         *a;              /* nonzero elements */ \
47   PetscScalar      *solve_work;     /* work space used in MatSolve */ \
48   IS                row, col, icol; /* index sets, used for reorderings */ \
49   PetscBool         pivotinblocks;  /* pivot inside factorization of each diagonal block */ \
50   Mat               parent;         /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
51                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \
52   Mat_SubSppt      *submatis1;      /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \
53   Mat_SeqAIJOps     ops[1]          /* operations for SeqAIJ and its subclasses */
54 
55 typedef struct {
56   MatTransposeColoring matcoloring;
57   Mat                  Bt_den;  /* dense matrix of B^T */
58   Mat                  ABt_den; /* dense matrix of A*B^T */
59   PetscBool            usecoloring;
60 } Mat_MatMatTransMult;
61 
62 typedef struct { /* used by MatTransposeMatMult() */
63   Mat   At;      /* transpose of the first matrix */
64   Mat   mA;      /* maij matrix of A */
65   Vec   bt, ct;  /* vectors to hold locally transposed arrays of B and C */
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 
137   /* MatSetValuesCOO() related fields on host */
138   PetscCount  coo_n; /* Number of entries in MatSetPreallocationCOO() */
139   PetscCount  Atot;  /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */
140   PetscCount *jmap;  /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */
141   PetscCount *perm;  /* The permutation array in sorting (i,j) by row and then by col */
142 } Mat_SeqAIJ;
143 
144 /*
145   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
146 */
147 static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i) {
148   Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data;
149   if (A->singlemalloc) {
150     PetscCall(PetscFree3(*a, *j, *i));
151   } else {
152     if (A->free_a) PetscCall(PetscFree(*a));
153     if (A->free_ij) PetscCall(PetscFree(*j));
154     if (A->free_ij) PetscCall(PetscFree(*i));
155   }
156   return 0;
157 }
158 /*
159     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
160     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
161 */
162 #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \
163   if (NROW >= RMAX) { \
164     Mat_SeqAIJ *Ain       = (Mat_SeqAIJ *)Amat->data; \
165     /* there is no extra room in row, therefore enlarge */ \
166     PetscInt    CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
167     datatype   *new_a; \
168 \
169     PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \
170     /* malloc new storage space */ \
171     PetscCall(PetscMalloc3(BS2 *new_nz, &new_a, new_nz, &new_j, AM + 1, &new_i)); \
172 \
173     /* copy over old data into new slots */ \
174     for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
175     for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
176     PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
177     len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
178     PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
179     PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \
180     PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \
181     PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), AA + BS2 * (AI[ROW] + NROW), BS2 * len)); \
182     /* free up old matrix storage */ \
183     PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
184     AA     = new_a; \
185     Ain->a = (MatScalar *)new_a; \
186     AI = Ain->i = new_i; \
187     AJ = Ain->j       = new_j; \
188     Ain->singlemalloc = PETSC_TRUE; \
189 \
190     RP   = AJ + AI[ROW]; \
191     AP   = AA + BS2 * AI[ROW]; \
192     RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
193     Ain->maxnz += BS2 * CHUNKSIZE; \
194     Ain->reallocs++; \
195   }
196 
197 #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \
198   if (NROW >= RMAX) { \
199     Mat_SeqAIJ *Ain       = (Mat_SeqAIJ *)Amat->data; \
200     /* there is no extra room in row, therefore enlarge */ \
201     PetscInt    CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
202 \
203     PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \
204     /* malloc new storage space */ \
205     PetscCall(PetscMalloc1(new_nz, &new_j)); \
206     PetscCall(PetscMalloc1(AM + 1, &new_i)); \
207 \
208     /* copy over old data into new slots */ \
209     for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
210     for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
211     PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
212     len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
213     PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
214 \
215     /* free up old matrix storage */ \
216     PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
217     Ain->a = NULL; \
218     AI = Ain->i = new_i; \
219     AJ = Ain->j       = new_j; \
220     Ain->singlemalloc = PETSC_FALSE; \
221     Ain->free_a       = PETSC_FALSE; \
222 \
223     RP   = AJ + AI[ROW]; \
224     RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
225     Ain->maxnz += BS2 * CHUNKSIZE; \
226     Ain->reallocs++; \
227   }
228 
229 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *);
230 PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);
231 PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_SeqAIJ(Mat);
232 
233 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *);
234 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
235 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *);
236 
237 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
238 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
239 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
240 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
241 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
242 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
243 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *);
244 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure);
245 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *);
246 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
247 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **);
248 
249 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec);
250 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec);
251 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec);
252 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec);
253 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec);
254 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
255 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
256 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
257 
258 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool);
259 
260 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
261 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
262 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]);
263 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *);
264 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *);
265 
266 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **);
267 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *);
268 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
269 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
270 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
271 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *);
272 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *);
273 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec);
274 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec);
275 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat, Vec, Vec);
276 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec);
277 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat, Vec, Vec);
278 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec);
279 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat, Vec, Vec);
280 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
281 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec);
282 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec);
283 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec);
284 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
285 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
286 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat, Mat, Mat);
287 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat);
288 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *);
289 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring);
290 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring);
291 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt);
292 PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer);
293 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer);
294 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer);
295 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
296 
297 #if defined(PETSC_HAVE_HYPRE)
298 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
299 #endif
300 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);
301 
302 PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
303 PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
304 PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);
305 
306 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
307 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat);
308 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat);
309 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat);
310 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat);
311 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat);
312 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat);
313 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat);
314 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat);
315 #if defined(PETSC_HAVE_HYPRE)
316 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
317 #endif
318 
319 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
320 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat);
321 
322 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat);
323 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat);
324 
325 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat);
326 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
327 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat);
328 
329 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
330 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat);
331 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat);
332 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
333 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat);
334 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat);
335 
336 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
337 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
338 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *);
339 
340 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
341 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
342 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring);
343 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat);
344 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat);
345 
346 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat);
347 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat);
348 
349 PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom);
350 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
351 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
352 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
353 PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar);
354 PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec);
355 PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode);
356 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure);
357 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
358 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
359 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
360 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
361 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
362 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
363 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
364 PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
365 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer);
366 
367 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
368 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
369 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
370 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
371 
372 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *);
373 
374 #if defined(PETSC_HAVE_MATLAB_ENGINE)
375 PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *);
376 PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *);
377 #endif
378 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *);
379 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *);
380 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *);
381 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
382 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
383 #if defined(PETSC_HAVE_SCALAPACK)
384 PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
385 #endif
386 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *);
387 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *);
388 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *);
389 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *);
390 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *);
391 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS);
392 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *);
393 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
394 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType);
395 PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat);
396 
397 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *);
398 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
399 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
400 
401 PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat);
402 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *);
403 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
404 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
405 PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]);
406 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *);
407 
408 PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *);
409 PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat);
410 
411 /*
412     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
413 
414   Input Parameters:
415 +  nnz - the number of entries
416 .  r - the array of vector values
417 .  xv - the matrix values for the row
418 -  xi - the column indices of the nonzeros in the row
419 
420   Output Parameter:
421 .  sum - negative the sum of results
422 
423   PETSc compile flags:
424 +   PETSC_KERNEL_USE_UNROLL_4
425 -   PETSC_KERNEL_USE_UNROLL_2
426 
427   Developer Note:
428     The macro changes sum but not other parameters
429 
430 .seealso: `PetscSparseDensePlusDot()`
431 */
432 #if defined(PETSC_KERNEL_USE_UNROLL_4)
433 #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
434   { \
435     if (nnz > 0) { \
436       PetscInt nnz2 = nnz, rem = nnz & 0x3; \
437       switch (rem) { \
438       case 3: sum -= *xv++ * r[*xi++]; \
439       case 2: sum -= *xv++ * r[*xi++]; \
440       case 1: sum -= *xv++ * r[*xi++]; nnz2 -= rem; \
441       } \
442       while (nnz2 > 0) { \
443         sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
444         xv += 4; \
445         xi += 4; \
446         nnz2 -= 4; \
447       } \
448       xv -= nnz; \
449       xi -= nnz; \
450     } \
451   }
452 
453 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
454 #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
455   { \
456     PetscInt __i, __i1, __i2; \
457     for (__i = 0; __i < nnz - 1; __i += 2) { \
458       __i1 = xi[__i]; \
459       __i2 = xi[__i + 1]; \
460       sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
461     } \
462     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \
463   }
464 
465 #else
466 #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
467   { \
468     PetscInt __i; \
469     for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \
470   }
471 #endif
472 
473 /*
474     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
475 
476   Input Parameters:
477 +  nnz - the number of entries
478 .  r - the array of vector values
479 .  xv - the matrix values for the row
480 -  xi - the column indices of the nonzeros in the row
481 
482   Output Parameter:
483 .  sum - the sum of results
484 
485   PETSc compile flags:
486 +   PETSC_KERNEL_USE_UNROLL_4
487 -   PETSC_KERNEL_USE_UNROLL_2
488 
489   Developer Note:
490     The macro changes sum but not other parameters
491 
492 .seealso: `PetscSparseDenseMinusDot()`
493 */
494 #if defined(PETSC_KERNEL_USE_UNROLL_4)
495 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
496   { \
497     if (nnz > 0) { \
498       PetscInt nnz2 = nnz, rem = nnz & 0x3; \
499       switch (rem) { \
500       case 3: sum += *xv++ * r[*xi++]; \
501       case 2: sum += *xv++ * r[*xi++]; \
502       case 1: sum += *xv++ * r[*xi++]; nnz2 -= rem; \
503       } \
504       while (nnz2 > 0) { \
505         sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
506         xv += 4; \
507         xi += 4; \
508         nnz2 -= 4; \
509       } \
510       xv -= nnz; \
511       xi -= nnz; \
512     } \
513   }
514 
515 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
516 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
517   { \
518     PetscInt __i, __i1, __i2; \
519     for (__i = 0; __i < nnz - 1; __i += 2) { \
520       __i1 = xi[__i]; \
521       __i2 = xi[__i + 1]; \
522       sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
523     } \
524     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \
525   }
526 
527 #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)
528 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz))
529 
530 #else
531 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
532   { \
533     PetscInt __i; \
534     for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \
535   }
536 #endif
537 
538 #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)
539 #include <immintrin.h>
540 #if !defined(_MM_SCALE_8)
541 #define _MM_SCALE_8 8
542 #endif
543 
544 static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n) {
545   __m512d  vec_x, vec_y, vec_vals;
546   __m256i  vec_idx;
547   __mmask8 mask;
548   PetscInt j;
549 
550   vec_y = _mm512_setzero_pd();
551   for (j = 0; j < (n >> 3); j++) {
552     vec_idx  = _mm256_loadu_si256((__m256i const *)aj);
553     vec_vals = _mm512_loadu_pd(aa);
554     vec_x    = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
555     vec_y    = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
556     aj += 8;
557     aa += 8;
558   }
559   /* masked load does not work on KNL, it requires avx512vl */
560   if ((n & 0x07) > 2) {
561     mask     = (__mmask8)(0xff >> (8 - (n & 0x07)));
562     vec_idx  = _mm256_loadu_si256((__m256i const *)aj);
563     vec_vals = _mm512_loadu_pd(aa);
564     vec_x    = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8);
565     vec_y    = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask);
566   } else if ((n & 0x07) == 2) {
567     *sum += aa[0] * x[aj[0]];
568     *sum += aa[1] * x[aj[1]];
569   } else if ((n & 0x07) == 1) {
570     *sum += aa[0] * x[aj[0]];
571   }
572   if (n > 2) *sum += _mm512_reduce_add_pd(vec_y);
573   /*
574   for (j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
575 */
576 }
577 #endif
578 
579 /*
580     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
581 
582   Input Parameters:
583 +  nnz - the number of entries
584 .  r - the array of vector values
585 .  xv - the matrix values for the row
586 -  xi - the column indices of the nonzeros in the row
587 
588   Output Parameter:
589 .  max - the max of results
590 
591 .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()`
592 */
593 #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \
594   { \
595     PetscInt __i; \
596     for (__i = 0; __i < nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]])); \
597   }
598 
599 /*
600  Add column indices into table for counting the max nonzeros of merged rows
601  */
602 #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \
603   { \
604     PetscInt _j, _row, _nz, *_col; \
605     if (mat) { \
606       for (_row = 0; _row < nrows; _row++) { \
607         _nz = mat->i[_row + 1] - mat->i[_row]; \
608         for (_j = 0; _j < _nz; _j++) { \
609           _col = _j + mat->j + mat->i[_row]; \
610           PetscTableAdd(ta, *_col + 1, 1, INSERT_VALUES); \
611         } \
612       } \
613     } \
614   }
615 
616 /*
617  Add column indices into table for counting the nonzeros of merged rows
618  */
619 #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \
620   { \
621     PetscInt _j, _row, _nz, *_col, _i; \
622     for (_i = 0; _i < nrows; _i++) { \
623       _row = rows[_i]; \
624       _nz  = mat->i[_row + 1] - mat->i[_row]; \
625       for (_j = 0; _j < _nz; _j++) { \
626         _col = _j + mat->j + mat->i[_row]; \
627         PetscTableAdd(ta, *_col + 1, 1, INSERT_VALUES); \
628       } \
629     } \
630   }
631 
632 #endif
633