Actual source code: aij.h

  1: #pragma once

  3: #include <petsc/private/matimpl.h>
  4: #include <petsc/private/hashmapi.h>
  5: #include <petsc/private/hashmapijv.h>

  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;

 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;

 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:   PetscInt          maxnz;       /* allocated nonzeros */ \
 51:   PetscInt         *imax;        /* maximum space allocated for each row */ \
 52:   PetscInt         *ilen;        /* actual length of each row */ \
 53:   PetscInt         *ipre;        /* space preallocated for each row by user */ \
 54:   PetscBool         free_imax_ilen; \
 55:   PetscInt          reallocs;           /* number of mallocs done during MatSetValues() \
 56:                                         as more values are set than were prealloced */ \
 57:   PetscInt          rmax;               /* max nonzeros in any row */ \
 58:   PetscBool         keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \
 59:   PetscBool         ignorezeroentries; \
 60:   PetscBool         free_ij;       /* free the column indices j and row offsets i when the matrix is destroyed */ \
 61:   PetscBool         free_a;        /* free the numerical values when matrix is destroy */ \
 62:   Mat_CompressedRow compressedrow; /* use compressed row format */ \
 63:   PetscInt          nz;            /* nonzeros */ \
 64:   PetscInt         *i;             /* pointer to beginning of each row */ \
 65:   PetscInt         *j;             /* column values: j + i[k] - 1 is start of row k */ \
 66:   PetscInt         *diag;          /* pointers to diagonal elements */ \
 67:   PetscInt          nonzerorowcnt; /* how many rows have nonzero entries */ \
 68:   PetscBool         free_diag; \
 69:   datatype         *a;              /* nonzero elements */ \
 70:   PetscScalar      *solve_work;     /* work space used in MatSolve */ \
 71:   IS                row, col, icol; /* index sets, used for reorderings */ \
 72:   PetscBool         pivotinblocks;  /* pivot inside factorization of each diagonal block */ \
 73:   Mat               parent;         /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
 74:                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \
 75:   Mat_SubSppt      *submatis1;      /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \
 76:   Mat_SeqAIJOps     ops[1]          /* operations for SeqAIJ and its subclasses */

 78: typedef struct {
 79:   MatTransposeColoring matcoloring;
 80:   Mat                  Bt_den;  /* dense matrix of B^T */
 81:   Mat                  ABt_den; /* dense matrix of A*B^T */
 82:   PetscBool            usecoloring;
 83: } Mat_MatMatTransMult;

 85: typedef struct { /* used by MatTransposeMatMult() */
 86:   Mat At;        /* transpose of the first matrix */
 87:   Mat mA;        /* maij matrix of A */
 88:   Vec bt, ct;    /* vectors to hold locally transposed arrays of B and C */
 89:   /* used by PtAP */
 90:   void *data;
 91:   PetscErrorCode (*destroy)(void *);
 92: } Mat_MatTransMatMult;

 94: typedef struct {
 95:   PetscInt    *api, *apj; /* symbolic structure of A*P */
 96:   PetscScalar *apa;       /* temporary array for storing one row of A*P */
 97: } Mat_AP;

 99: typedef struct {
100:   MatTransposeColoring matcoloring;
101:   Mat                  Rt;   /* sparse or dense matrix of R^T */
102:   Mat                  RARt; /* dense matrix of R*A*R^T */
103:   Mat                  ARt;  /* A*R^T used for the case -matrart_color_art */
104:   MatScalar           *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
105:   /* free intermediate products needed for PtAP */
106:   void *data;
107:   PetscErrorCode (*destroy)(void *);
108: } Mat_RARt;

110: typedef struct {
111:   Mat BC; /* temp matrix for storing B*C */
112: } Mat_MatMatMatMult;

114: /*
115:   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
116:   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
117:   j[i[k]+p] is the pth column in row k.  Note that the diagonal
118:   matrix elements are stored with the rest of the nonzeros (not separately).
119: */

121: /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
122: typedef struct {
123:   MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
124:   PetscInt   bdiagsize;                  /* length of bdiag and ibdiag */
125:   PetscBool  ibdiagvalid;                /* do ibdiag[] and bdiag[] contain the most recent values */

127:   PetscBool        use;
128:   PetscInt         node_count;       /* number of inodes */
129:   PetscInt        *size;             /* size of each inode */
130:   PetscInt         limit;            /* inode limit */
131:   PetscInt         max_limit;        /* maximum supported inode limit */
132:   PetscBool        checked;          /* if inodes have been checked for */
133:   PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */
134: } Mat_SeqAIJ_Inode;

136: PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer);
137: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType);
138: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
139: PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
140: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool);
141: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *);
142: PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool);
143: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *);
144: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **);
145: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **);

147: typedef struct {
148:   SEQAIJHEADER(MatScalar);
149:   Mat_SeqAIJ_Inode inode;
150:   MatScalar       *saved_values; /* location for stashing nonzero values of matrix */

152:   PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
153:   PetscBool    idiagvalid;                /* current idiag[] and mdiag[] are valid */
154:   PetscScalar *ibdiag;                    /* inverses of block diagonals */
155:   PetscBool    ibdiagvalid;               /* inverses of block diagonals are valid. */
156:   PetscBool    diagonaldense;             /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
157:   PetscScalar  fshift, omega;             /* last used omega and fshift */

159:   /* MatSetValues() via hash related fields */
160:   PetscHMapIJV   ht;
161:   PetscInt      *dnz;
162:   struct _MatOps cops;
163: } Mat_SeqAIJ;

165: typedef struct {
166:   PetscInt    nz;   /* nz of the matrix after assembly */
167:   PetscCount  n;    /* Number of entries in MatSetPreallocationCOO() */
168:   PetscCount  Atot; /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */
169:   PetscCount *jmap; /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */
170:   PetscCount *perm; /* The permutation array in sorting (i,j) by row and then by col */
171: } MatCOOStruct_SeqAIJ;

173: #define MatSeqXAIJGetOptions_Private(A) \
174:   { \
175:     const PetscBool oldvalues = (PetscBool)(A != PETSC_NULLPTR); \
176:     PetscInt        nonew = 0, nounused = 0; \
177:     PetscBool       roworiented = PETSC_FALSE; \
178:     if (oldvalues) { \
179:       nonew       = ((Mat_SeqAIJ *)A->data)->nonew; \
180:       nounused    = ((Mat_SeqAIJ *)A->data)->nounused; \
181:       roworiented = ((Mat_SeqAIJ *)A->data)->roworiented; \
182:     } \
183:     (void)0

185: #define MatSeqXAIJRestoreOptions_Private(A) \
186:   if (oldvalues) { \
187:     ((Mat_SeqAIJ *)A->data)->nonew       = nonew; \
188:     ((Mat_SeqAIJ *)A->data)->nounused    = nounused; \
189:     ((Mat_SeqAIJ *)A->data)->roworiented = roworiented; \
190:   } \
191:   } \
192:   (void)0

194: static inline PetscErrorCode MatXAIJAllocatea(Mat A, PetscInt nz, PetscScalar **array)
195: {
196:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

198:   PetscFunctionBegin;
199:   PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscScalar), (void **)array));
200:   a->free_a = PETSC_TRUE;
201:   PetscFunctionReturn(PETSC_SUCCESS);
202: }

204: static inline PetscErrorCode MatXAIJDeallocatea(Mat A, PetscScalar **array)
205: {
206:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

208:   PetscFunctionBegin;
209:   if (a->free_a) PetscCall(PetscShmgetDeallocateArray((void **)array));
210:   a->free_a = PETSC_FALSE;
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /*
215:   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
216: */
217: static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i)
218: {
219:   Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data;

221:   PetscFunctionBegin;
222:   if (A->free_a) PetscCall(PetscShmgetDeallocateArray((void **)a));
223:   if (A->free_ij) PetscCall(PetscShmgetDeallocateArray((void **)j));
224:   if (A->free_ij) PetscCall(PetscShmgetDeallocateArray((void **)i));
225:   PetscFunctionReturn(PETSC_SUCCESS);
226: }
227: /*
228:     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
229:     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
230: */
231: #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \
232:   do { \
233:     if (NROW >= RMAX) { \
234:       Mat_SeqAIJ *Ain       = (Mat_SeqAIJ *)Amat->data; \
235:       PetscInt    CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
236:       datatype   *new_a; \
237: \
238:       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); \
239:       /* malloc new storage space */ \
240:       PetscCall(PetscShmgetAllocateArray(BS2 *new_nz, sizeof(PetscScalar), (void **)&new_a)); \
241:       PetscCall(PetscShmgetAllocateArray(new_nz, sizeof(PetscInt), (void **)&new_j)); \
242:       PetscCall(PetscShmgetAllocateArray(AM + 1, sizeof(PetscInt), (void **)&new_i)); \
243:       Ain->free_a  = PETSC_TRUE; \
244:       Ain->free_ij = PETSC_TRUE; \
245:       /* copy over old data into new slots */ \
246:       for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
247:       for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
248:       PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
249:       len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
250:       PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, PetscSafePointerPlusOffset(AJ, AI[ROW] + NROW), len)); \
251:       PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \
252:       PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \
253:       PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), PetscSafePointerPlusOffset(AA, BS2 * (AI[ROW] + NROW)), BS2 * len)); \
254:       /* free up old matrix storage */ \
255:       PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
256:       AA     = new_a; \
257:       Ain->a = (MatScalar *)new_a; \
258:       AI = Ain->i = new_i; \
259:       AJ = Ain->j = new_j; \
260: \
261:       RP   = AJ + AI[ROW]; \
262:       AP   = AA + BS2 * AI[ROW]; \
263:       RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
264:       Ain->maxnz += BS2 * CHUNKSIZE; \
265:       Ain->reallocs++; \
266:       Amat->nonzerostate++; \
267:     } \
268:   } while (0)

270: #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \
271:   do { \
272:     if (NROW >= RMAX) { \
273:       Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
274:       /* there is no extra room in row, therefore enlarge */ \
275:       PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
276: \
277:       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); \
278:       /* malloc new storage space */ \
279:       PetscCall(PetscShmgetAllocateArray(new_nz, sizeof(PetscInt), (void **)&new_j)); \
280:       PetscCall(PetscShmgetAllocateArray(AM + 1, sizeof(PetscInt), (void **)&new_i)); \
281:       Ain->free_a  = PETSC_FALSE; \
282:       Ain->free_ij = PETSC_TRUE; \
283: \
284:       /* copy over old data into new slots */ \
285:       for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
286:       for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
287:       PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
288:       len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
289:       PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
290: \
291:       /* free up old matrix storage */ \
292:       PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
293:       Ain->a = NULL; \
294:       AI = Ain->i = new_i; \
295:       AJ = Ain->j = new_j; \
296: \
297:       RP   = AJ + AI[ROW]; \
298:       RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
299:       Ain->maxnz += BS2 * CHUNKSIZE; \
300:       Ain->reallocs++; \
301:       Amat->nonzerostate++; \
302:     } \
303:   } while (0)

305: PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *);
306: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);

308: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
309: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *);

311: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
312: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
313: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
314: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
315: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *);
316: PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure);
317: PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *);
318: PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
319: PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **);

321: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec);
322: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec);
323: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec);
324: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec);
325: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec);
326: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
327: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
328: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);

330: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool);

332: PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
333: PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
334: PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]);
335: PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *);
336: PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *);

338: PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **);
339: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
340: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
341: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
342: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *);
343: PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *);
344: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec);
345: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec);
346: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec);
347: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec);
348: PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec);
349: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec);
350: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec);
351: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
352: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
353: PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat);
354: PETSC_INTERN PetscErrorCode MatMatSolveTranspose_SeqAIJ(Mat, Mat, Mat);
355: PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *);
356: PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring);
357: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring);
358: PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt);
359: PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer);
360: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer);
361: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer);
362: PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);

364: #if defined(PETSC_HAVE_HYPRE)
365: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
366: #endif
367: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);

369: PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
370: PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
371: PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);

373: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
374: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat);
375: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat);
376: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat);
377: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat);
378: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat);
379: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat);
380: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat);
381: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat);
382: #if defined(PETSC_HAVE_HYPRE)
383: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
384: #endif

386: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
387: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat);

389: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat);
390: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat);

392: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat);
393: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
394: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat);

396: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
397: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat);
398: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat);
399: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
400: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat);
401: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat);

403: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
404: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
405: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *);

407: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
408: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
409: PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring);
410: PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat);
411: PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat);

413: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat);
414: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat);

416: PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom);
417: PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
418: PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
419: PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
420: PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar);
421: PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec);
422: PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode);
423: PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure);
424: PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
425: PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
426: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
427: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
428: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
429: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
430: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
431: PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer);

433: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
434: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
435: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
436: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);

438: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *);

440: #if defined(PETSC_HAVE_MATLAB)
441: PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *);
442: PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *);
443: #endif
444: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *);
445: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *);
446: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *);
447: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
448: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
449: #if defined(PETSC_HAVE_SCALAPACK)
450: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
451: #endif
452: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *);
453: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *);
454: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *);
455: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *);
456: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *);
457: PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS);
458: PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *);
459: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
460: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType);
461: PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat);

463: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *);
464: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
465: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);

467: PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat);
468: PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqAIJ(Mat, PetscBool);
469: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *);
470: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
471: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
472: PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]);
473: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *);

475: PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat);

477: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *);

479: /*
480:     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage

482:   Input Parameters:
483: +  nnz - the number of entries
484: .  r - the array of vector values
485: .  xv - the matrix values for the row
486: -  xi - the column indices of the nonzeros in the row

488:   Output Parameter:
489: .  sum - negative the sum of results

491:   PETSc compile flags:
492: +   PETSC_KERNEL_USE_UNROLL_4
493: -   PETSC_KERNEL_USE_UNROLL_2

495:   Developer Note:
496:     The macro changes sum but not other parameters

498: .seealso: `PetscSparseDensePlusDot()`
499: */
500: #if defined(PETSC_KERNEL_USE_UNROLL_4)
501:   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
502:     do { \
503:       if (nnz > 0) { \
504:         PetscInt nnz2 = nnz, rem = nnz & 0x3; \
505:         switch (rem) { \
506:         case 3: \
507:           sum -= *xv++ * r[*xi++]; \
508:         case 2: \
509:           sum -= *xv++ * r[*xi++]; \
510:         case 1: \
511:           sum -= *xv++ * r[*xi++]; \
512:           nnz2 -= rem; \
513:         } \
514:         while (nnz2 > 0) { \
515:           sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
516:           xv += 4; \
517:           xi += 4; \
518:           nnz2 -= 4; \
519:         } \
520:         xv -= nnz; \
521:         xi -= nnz; \
522:       } \
523:     } while (0)

525: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
526:   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
527:     do { \
528:       PetscInt __i, __i1, __i2; \
529:       for (__i = 0; __i < nnz - 1; __i += 2) { \
530:         __i1 = xi[__i]; \
531:         __i2 = xi[__i + 1]; \
532:         sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
533:       } \
534:       if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \
535:     } while (0)

537: #else
538:   #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
539:     do { \
540:       PetscInt __i; \
541:       for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \
542:     } while (0)
543: #endif

545: /*
546:     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage

548:   Input Parameters:
549: +  nnz - the number of entries
550: .  r - the array of vector values
551: .  xv - the matrix values for the row
552: -  xi - the column indices of the nonzeros in the row

554:   Output Parameter:
555: .  sum - the sum of results

557:   PETSc compile flags:
558: +   PETSC_KERNEL_USE_UNROLL_4
559: -   PETSC_KERNEL_USE_UNROLL_2

561:   Developer Note:
562:     The macro changes sum but not other parameters

564: .seealso: `PetscSparseDenseMinusDot()`
565: */
566: #if defined(PETSC_KERNEL_USE_UNROLL_4)
567:   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
568:     do { \
569:       if (nnz > 0) { \
570:         PetscInt nnz2 = nnz, rem = nnz & 0x3; \
571:         switch (rem) { \
572:         case 3: \
573:           sum += *xv++ * r[*xi++]; \
574:         case 2: \
575:           sum += *xv++ * r[*xi++]; \
576:         case 1: \
577:           sum += *xv++ * r[*xi++]; \
578:           nnz2 -= rem; \
579:         } \
580:         while (nnz2 > 0) { \
581:           sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
582:           xv += 4; \
583:           xi += 4; \
584:           nnz2 -= 4; \
585:         } \
586:         xv -= nnz; \
587:         xi -= nnz; \
588:       } \
589:     } while (0)

591: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
592:   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
593:     do { \
594:       PetscInt __i, __i1, __i2; \
595:       for (__i = 0; __i < nnz - 1; __i += 2) { \
596:         __i1 = xi[__i]; \
597:         __i2 = xi[__i + 1]; \
598:         sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
599:       } \
600:       if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \
601:     } while (0)

603: #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)
604:   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz))

606: #else
607:   #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
608:     do { \
609:       PetscInt __i; \
610:       for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \
611:     } while (0)
612: #endif

614: #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)
615:   #include <immintrin.h>
616:   #if !defined(_MM_SCALE_8)
617:     #define _MM_SCALE_8 8
618:   #endif

620: static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n)
621: {
622:   __m512d  vec_x, vec_y, vec_vals;
623:   __m256i  vec_idx;
624:   PetscInt j;

626:   vec_y = _mm512_setzero_pd();
627:   for (j = 0; j < (n >> 3); j++) {
628:     vec_idx  = _mm256_loadu_si256((__m256i const *)aj);
629:     vec_vals = _mm512_loadu_pd(aa);
630:     vec_x    = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
631:     vec_y    = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
632:     aj += 8;
633:     aa += 8;
634:   }
635:   #if defined(__AVX512VL__)
636:   /* masked load requires avx512vl, which is not supported by KNL */
637:   if (n & 0x07) {
638:     __mmask8 mask;
639:     mask     = (__mmask8)(0xff >> (8 - (n & 0x07)));
640:     vec_idx  = _mm256_mask_loadu_epi32(vec_idx, mask, aj);
641:     vec_vals = _mm512_mask_loadu_pd(vec_vals, mask, aa);
642:     vec_x    = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8);
643:     vec_y    = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask);
644:   }
645:   *sum += _mm512_reduce_add_pd(vec_y);
646:   #else
647:   *sum += _mm512_reduce_add_pd(vec_y);
648:   for (j = 0; j < (n & 0x07); j++) *sum += aa[j] * x[aj[j]];
649:   #endif
650: }
651: #endif

653: /*
654:     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage

656:   Input Parameters:
657: +  nnz - the number of entries
658: .  r - the array of vector values
659: .  xv - the matrix values for the row
660: -  xi - the column indices of the nonzeros in the row

662:   Output Parameter:
663: .  max - the max of results

665: .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()`
666: */
667: #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \
668:   do { \
669:     for (PetscInt __i = 0; __i < (nnz); __i++) { max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); } \
670:   } while (0)

672: /*
673:  Add column indices into table for counting the max nonzeros of merged rows
674:  */
675: #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \
676:   do { \
677:     if (mat) { \
678:       for (PetscInt _row = 0; _row < (nrows); _row++) { \
679:         const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \
680:         for (PetscInt _j = 0; _j < _nz; _j++) { \
681:           PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
682:           PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \
683:         } \
684:       } \
685:     } \
686:   } while (0)

688: /*
689:  Add column indices into table for counting the nonzeros of merged rows
690:  */
691: #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \
692:   do { \
693:     for (PetscInt _i = 0; _i < (nrows); _i++) { \
694:       const PetscInt _row = (rows)[_i]; \
695:       const PetscInt _nz  = (mat)->i[_row + 1] - (mat)->i[_row]; \
696:       for (PetscInt _j = 0; _j < _nz; _j++) { \
697:         PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
698:         PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \
699:       } \
700:     } \
701:   } while (0)