xref: /petsc/src/mat/impls/aij/seq/aij.h (revision cedd07cade5cbdfdad435c8172b7ec8972d9cd8d)
1 
2 #ifndef __AIJ_H
3 #define __AIJ_H
4 
5 #include <petsc/private/matimpl.h>
6 #include <petsc/private/hashmapi.h>
7 
8 /*
9  Used by MatCreateSubMatrices_MPIXAIJ_Local()
10 */
11 typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
12   PetscInt   id; /* index of submats, only submats[0] is responsible for deleting some arrays below */
13   PetscInt   nrqs, nrqr;
14   PetscInt **rbuf1, **rbuf2, **rbuf3, **sbuf1, **sbuf2;
15   PetscInt **ptr;
16   PetscInt  *tmp;
17   PetscInt  *ctr;
18   PetscInt  *pa; /* proc array */
19   PetscInt  *req_size, *req_source1, *req_source2;
20   PetscBool  allcolumns, allrows;
21   PetscBool  singleis;
22   PetscInt  *row2proc; /* row to proc 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   PetscBool         singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \
52   PetscInt          maxnz;        /* allocated nonzeros */ \
53   PetscInt         *imax;         /* maximum space allocated for each row */ \
54   PetscInt         *ilen;         /* actual length of each row */ \
55   PetscInt         *ipre;         /* space preallocated for each row by user */ \
56   PetscBool         free_imax_ilen; \
57   PetscInt          reallocs;           /* number of mallocs done during MatSetValues() \
58                                         as more values are set than were prealloced */ \
59   PetscInt          rmax;               /* max nonzeros in any row */ \
60   PetscBool         keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \
61   PetscBool         ignorezeroentries; \
62   PetscBool         free_ij;       /* free the column indices j and row offsets i when the matrix is destroyed */ \
63   PetscBool         free_a;        /* free the numerical values when matrix is destroy */ \
64   Mat_CompressedRow compressedrow; /* use compressed row format */ \
65   PetscInt          nz;            /* nonzeros */ \
66   PetscInt         *i;             /* pointer to beginning of each row */ \
67   PetscInt         *j;             /* column values: j + i[k] - 1 is start of row k */ \
68   PetscInt         *diag;          /* pointers to diagonal elements */ \
69   PetscInt          nonzerorowcnt; /* how many rows have nonzero entries */ \
70   PetscBool         free_diag; \
71   datatype         *a;              /* nonzero elements */ \
72   PetscScalar      *solve_work;     /* work space used in MatSolve */ \
73   IS                row, col, icol; /* index sets, used for reorderings */ \
74   PetscBool         pivotinblocks;  /* pivot inside factorization of each diagonal block */ \
75   Mat               parent;         /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
76                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \
77   Mat_SubSppt      *submatis1;      /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \
78   Mat_SeqAIJOps     ops[1]          /* operations for SeqAIJ and its subclasses */
79 
80 typedef struct {
81   MatTransposeColoring matcoloring;
82   Mat                  Bt_den;  /* dense matrix of B^T */
83   Mat                  ABt_den; /* dense matrix of A*B^T */
84   PetscBool            usecoloring;
85 } Mat_MatMatTransMult;
86 
87 typedef struct { /* used by MatTransposeMatMult() */
88   Mat At;        /* transpose of the first matrix */
89   Mat mA;        /* maij matrix of A */
90   Vec bt, ct;    /* vectors to hold locally transposed arrays of B and C */
91   /* used by PtAP */
92   void *data;
93   PetscErrorCode (*destroy)(void *);
94 } Mat_MatTransMatMult;
95 
96 typedef struct {
97   PetscInt    *api, *apj; /* symbolic structure of A*P */
98   PetscScalar *apa;       /* temporary array for storing one row of A*P */
99 } Mat_AP;
100 
101 typedef struct {
102   MatTransposeColoring matcoloring;
103   Mat                  Rt;   /* sparse or dense matrix of R^T */
104   Mat                  RARt; /* dense matrix of R*A*R^T */
105   Mat                  ARt;  /* A*R^T used for the case -matrart_color_art */
106   MatScalar           *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
107   /* free intermediate products needed for PtAP */
108   void *data;
109   PetscErrorCode (*destroy)(void *);
110 } Mat_RARt;
111 
112 typedef struct {
113   Mat BC; /* temp matrix for storing B*C */
114 } Mat_MatMatMatMult;
115 
116 /*
117   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
118   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
119   j[i[k]+p] is the pth column in row k.  Note that the diagonal
120   matrix elements are stored with the rest of the nonzeros (not separately).
121 */
122 
123 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
124 typedef struct {
125   MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
126   PetscInt   bdiagsize;                  /* length of bdiag and ibdiag */
127   PetscBool  ibdiagvalid;                /* do ibdiag[] and bdiag[] contain the most recent values */
128 
129   PetscBool        use;
130   PetscInt         node_count;       /* number of inodes */
131   PetscInt        *size;             /* size of each inode */
132   PetscInt         limit;            /* inode limit */
133   PetscInt         max_limit;        /* maximum supported inode limit */
134   PetscBool        checked;          /* if inodes have been checked for */
135   PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */
136 } Mat_SeqAIJ_Inode;
137 
138 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer);
139 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType);
140 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
141 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
142 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool);
143 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *);
144 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool);
145 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat, Mat, const MatFactorInfo *);
146 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *);
147 PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **);
148 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **);
149 
150 typedef struct {
151   SEQAIJHEADER(MatScalar);
152   Mat_SeqAIJ_Inode inode;
153   MatScalar       *saved_values; /* location for stashing nonzero values of matrix */
154 
155   PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
156   PetscBool    idiagvalid;                /* current idiag[] and mdiag[] are valid */
157   PetscScalar *ibdiag;                    /* inverses of block diagonals */
158   PetscBool    ibdiagvalid;               /* inverses of block diagonals are valid. */
159   PetscBool    diagonaldense;             /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
160   PetscScalar  fshift, omega;             /* last used omega and fshift */
161 
162   /* MatSetValuesCOO() related fields on host */
163   PetscCount  coo_n; /* Number of entries in MatSetPreallocationCOO() */
164   PetscCount  Atot;  /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */
165   PetscCount *jmap;  /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */
166   PetscCount *perm;  /* The permutation array in sorting (i,j) by row and then by col */
167 } Mat_SeqAIJ;
168 
169 /*
170   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
171 */
172 static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i)
173 {
174   Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data;
175   if (A->singlemalloc) {
176     PetscCall(PetscFree3(*a, *j, *i));
177   } else {
178     if (A->free_a) PetscCall(PetscFree(*a));
179     if (A->free_ij) PetscCall(PetscFree(*j));
180     if (A->free_ij) PetscCall(PetscFree(*i));
181   }
182   return 0;
183 }
184 /*
185     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
186     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
187 */
188 #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \
189   if (NROW >= RMAX) { \
190     Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
191     /* there is no extra room in row, therefore enlarge */ \
192     PetscInt  CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
193     datatype *new_a; \
194 \
195     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); \
196     /* malloc new storage space */ \
197     PetscCall(PetscMalloc3(BS2 *new_nz, &new_a, new_nz, &new_j, AM + 1, &new_i)); \
198 \
199     /* copy over old data into new slots */ \
200     for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
201     for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
202     PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
203     len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
204     PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
205     PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \
206     PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \
207     PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), AA + BS2 * (AI[ROW] + NROW), BS2 * len)); \
208     /* free up old matrix storage */ \
209     PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
210     AA     = new_a; \
211     Ain->a = (MatScalar *)new_a; \
212     AI = Ain->i = new_i; \
213     AJ = Ain->j       = new_j; \
214     Ain->singlemalloc = PETSC_TRUE; \
215 \
216     RP   = AJ + AI[ROW]; \
217     AP   = AA + BS2 * AI[ROW]; \
218     RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
219     Ain->maxnz += BS2 * CHUNKSIZE; \
220     Ain->reallocs++; \
221   }
222 
223 #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \
224   if (NROW >= RMAX) { \
225     Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
226     /* there is no extra room in row, therefore enlarge */ \
227     PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
228 \
229     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); \
230     /* malloc new storage space */ \
231     PetscCall(PetscMalloc1(new_nz, &new_j)); \
232     PetscCall(PetscMalloc1(AM + 1, &new_i)); \
233 \
234     /* copy over old data into new slots */ \
235     for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
236     for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
237     PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
238     len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
239     PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
240 \
241     /* free up old matrix storage */ \
242     PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
243     Ain->a = NULL; \
244     AI = Ain->i = new_i; \
245     AJ = Ain->j       = new_j; \
246     Ain->singlemalloc = PETSC_FALSE; \
247     Ain->free_a       = PETSC_FALSE; \
248 \
249     RP   = AJ + AI[ROW]; \
250     RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
251     Ain->maxnz += BS2 * CHUNKSIZE; \
252     Ain->reallocs++; \
253   }
254 
255 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *);
256 PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);
257 PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_SeqAIJ(Mat);
258 
259 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *);
260 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
261 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *);
262 
263 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
264 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
265 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, const MatFactorInfo *);
266 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
267 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
268 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
269 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *);
270 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure);
271 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *);
272 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
273 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **);
274 
275 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec);
276 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec);
277 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec);
278 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec);
279 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec);
280 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
281 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
282 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
283 
284 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool);
285 
286 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
287 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
288 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]);
289 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *);
290 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *);
291 
292 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **);
293 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat, Mat, IS, IS, const MatFactorInfo *);
294 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
295 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
296 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
297 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *);
298 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *);
299 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec);
300 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec);
301 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat, Vec, Vec);
302 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec);
303 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat, Vec, Vec);
304 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec);
305 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat, Vec, Vec);
306 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
307 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec);
308 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec);
309 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec);
310 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
311 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
312 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat, Mat, Mat);
313 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat);
314 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *);
315 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring);
316 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring);
317 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt);
318 PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer);
319 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer);
320 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer);
321 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
322 
323 #if defined(PETSC_HAVE_HYPRE)
324 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
325 #endif
326 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);
327 
328 PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
329 PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
330 PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);
331 
332 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
333 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat);
334 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat);
335 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat);
336 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat);
337 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat);
338 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat);
339 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat);
340 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat);
341 #if defined(PETSC_HAVE_HYPRE)
342 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
343 #endif
344 
345 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
346 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat);
347 
348 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat);
349 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat);
350 
351 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat);
352 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
353 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat);
354 
355 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
356 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat);
357 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat);
358 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
359 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat);
360 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat);
361 
362 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
363 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
364 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *);
365 
366 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
367 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
368 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring);
369 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat);
370 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat);
371 
372 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat);
373 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat);
374 
375 PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom);
376 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
377 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
378 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
379 PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar);
380 PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec);
381 PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode);
382 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure);
383 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
384 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
385 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
386 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
387 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
388 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
389 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
390 PETSC_INTERN PetscErrorCode MatSetUp_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 A);
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     { \
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     }
483 
484 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
485   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
486     { \
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     }
495 
496 #else
497   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
498     { \
499       PetscInt __i; \
500       for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \
501     }
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     { \
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     }
549 
550 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
551   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
552     { \
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     }
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     { \
568       PetscInt __i; \
569       for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \
570     }
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   __mmask8 mask;
584   PetscInt j;
585 
586   vec_y = _mm512_setzero_pd();
587   for (j = 0; j < (n >> 3); j++) {
588     vec_idx  = _mm256_loadu_si256((__m256i const *)aj);
589     vec_vals = _mm512_loadu_pd(aa);
590     vec_x    = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
591     vec_y    = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
592     aj += 8;
593     aa += 8;
594   }
595   /* masked load does not work on KNL, it requires avx512vl */
596   if ((n & 0x07) > 2) {
597     mask     = (__mmask8)(0xff >> (8 - (n & 0x07)));
598     vec_idx  = _mm256_loadu_si256((__m256i const *)aj);
599     vec_vals = _mm512_loadu_pd(aa);
600     vec_x    = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8);
601     vec_y    = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask);
602   } else if ((n & 0x07) == 2) {
603     *sum += aa[0] * x[aj[0]];
604     *sum += aa[1] * x[aj[1]];
605   } else if ((n & 0x07) == 1) {
606     *sum += aa[0] * x[aj[0]];
607   }
608   if (n > 2) *sum += _mm512_reduce_add_pd(vec_y);
609   /*
610   for (j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
611 */
612 }
613 #endif
614 
615 /*
616     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
617 
618   Input Parameters:
619 +  nnz - the number of entries
620 .  r - the array of vector values
621 .  xv - the matrix values for the row
622 -  xi - the column indices of the nonzeros in the row
623 
624   Output Parameter:
625 .  max - the max of results
626 
627 .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()`
628 */
629 #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \
630   do { \
631     for (PetscInt __i = 0; __i < (nnz); __i++) { max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); } \
632   } while (0)
633 
634 /*
635  Add column indices into table for counting the max nonzeros of merged rows
636  */
637 #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \
638   do { \
639     if ((mat)) { \
640       for (PetscInt _row = 0; _row < (nrows); _row++) { \
641         const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \
642         for (PetscInt _j = 0; _j < _nz; _j++) { \
643           PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
644           PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \
645         } \
646       } \
647     } \
648   } while (0)
649 
650 /*
651  Add column indices into table for counting the nonzeros of merged rows
652  */
653 #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \
654   do { \
655     for (PetscInt _i = 0; _i < (nrows); _i++) { \
656       const PetscInt _row = (rows)[_i]; \
657       const PetscInt _nz  = (mat)->i[_row + 1] - (mat)->i[_row]; \
658       for (PetscInt _j = 0; _j < _nz; _j++) { \
659         PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
660         PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \
661       } \
662     } \
663   } while (0)
664 
665 #endif
666