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