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