Actual source code: aij.h
1: #ifndef PETSC_MATAIJ_IMPL_H
2: #define PETSC_MATAIJ_IMPL_H
4: #include <petsc/private/matimpl.h>
5: #include <petsc/private/hashmapi.h>
6: #include <petsc/private/hashmapijv.h>
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;
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;
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 */
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;
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;
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;
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;
112: typedef struct {
113: Mat BC; /* temp matrix for storing B*C */
114: } Mat_MatMatMatMult;
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: */
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 */
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;
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(Mat, Mat, const MatFactorInfo *);
146: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **);
147: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **);
149: typedef struct {
150: SEQAIJHEADER(MatScalar);
151: Mat_SeqAIJ_Inode inode;
152: MatScalar *saved_values; /* location for stashing nonzero values of matrix */
154: PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
155: PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */
156: PetscScalar *ibdiag; /* inverses of block diagonals */
157: PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */
158: PetscBool diagonaldense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
159: PetscScalar fshift, omega; /* last used omega and fshift */
161: /* MatSetValuesCOO() related fields on host */
162: PetscCount coo_n; /* Number of entries in MatSetPreallocationCOO() */
163: PetscCount Atot; /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */
164: PetscCount *jmap; /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */
165: PetscCount *perm; /* The permutation array in sorting (i,j) by row and then by col */
167: /* MatSetValues() via hash related fields */
168: PetscHMapIJV ht;
169: PetscInt *dnz;
170: struct _MatOps cops;
171: } Mat_SeqAIJ;
173: /*
174: Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
175: */
176: static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i)
177: {
178: Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data;
180: PetscFunctionBegin;
181: if (A->singlemalloc) {
182: PetscCall(PetscFree3(*a, *j, *i));
183: } else {
184: if (A->free_a) PetscCall(PetscFree(*a));
185: if (A->free_ij) PetscCall(PetscFree(*j));
186: if (A->free_ij) PetscCall(PetscFree(*i));
187: }
188: PetscFunctionReturn(PETSC_SUCCESS);
189: }
190: /*
191: Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
192: This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
193: */
194: #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \
195: if (NROW >= RMAX) { \
196: Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
197: /* there is no extra room in row, therefore enlarge */ \
198: PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
199: datatype *new_a; \
200: \
201: 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); \
202: /* malloc new storage space */ \
203: PetscCall(PetscMalloc3(BS2 *new_nz, &new_a, new_nz, &new_j, AM + 1, &new_i)); \
204: \
205: /* copy over old data into new slots */ \
206: for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
207: for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
208: PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
209: len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
210: PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
211: PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \
212: PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \
213: PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), AA + BS2 * (AI[ROW] + NROW), BS2 * len)); \
214: /* free up old matrix storage */ \
215: PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
216: AA = new_a; \
217: Ain->a = (MatScalar *)new_a; \
218: AI = Ain->i = new_i; \
219: AJ = Ain->j = new_j; \
220: Ain->singlemalloc = PETSC_TRUE; \
221: \
222: RP = AJ + AI[ROW]; \
223: AP = AA + BS2 * AI[ROW]; \
224: RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
225: Ain->maxnz += BS2 * CHUNKSIZE; \
226: Ain->reallocs++; \
227: }
229: #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \
230: if (NROW >= RMAX) { \
231: Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \
232: /* there is no extra room in row, therefore enlarge */ \
233: PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \
234: \
235: 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); \
236: /* malloc new storage space */ \
237: PetscCall(PetscMalloc1(new_nz, &new_j)); \
238: PetscCall(PetscMalloc1(AM + 1, &new_i)); \
239: \
240: /* copy over old data into new slots */ \
241: for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \
242: for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \
243: PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \
244: len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
245: PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \
246: \
247: /* free up old matrix storage */ \
248: PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \
249: Ain->a = NULL; \
250: AI = Ain->i = new_i; \
251: AJ = Ain->j = new_j; \
252: Ain->singlemalloc = PETSC_FALSE; \
253: Ain->free_a = PETSC_FALSE; \
254: \
255: RP = AJ + AI[ROW]; \
256: RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
257: Ain->maxnz += BS2 * CHUNKSIZE; \
258: Ain->reallocs++; \
259: }
261: PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *);
262: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]);
263: PETSC_INTERN PetscErrorCode MatResetPreallocationCOO_SeqAIJ(Mat);
265: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
266: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *);
268: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
269: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *);
270: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
271: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
272: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *);
273: PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure);
274: PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *);
275: PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
276: PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **);
278: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec);
279: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec);
280: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec);
281: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec);
282: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec);
283: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
284: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
285: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
287: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool);
289: PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
290: PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]);
291: PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]);
292: PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *);
293: PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *);
295: PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **);
296: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *);
297: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *);
298: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *);
299: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *);
300: PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *);
301: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec);
302: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec);
303: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec);
304: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec);
305: PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec);
306: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec);
307: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec);
308: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec);
309: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec);
310: PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat);
311: PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *);
312: PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring);
313: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring);
314: PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt);
315: PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer);
316: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer);
317: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer);
318: PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
320: #if defined(PETSC_HAVE_HYPRE)
321: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
322: #endif
323: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);
325: PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
326: PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
327: PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);
329: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
330: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat);
331: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat);
332: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat);
333: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat);
334: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat);
335: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat);
336: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat);
337: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat);
338: #if defined(PETSC_HAVE_HYPRE)
339: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
340: #endif
342: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
343: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat);
345: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat);
346: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat);
348: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat);
349: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
350: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat);
352: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
353: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat);
354: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat);
355: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
356: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat);
357: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat);
359: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
360: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
361: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *);
363: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat);
364: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat);
365: PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring);
366: PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat);
367: PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat);
369: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat);
370: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat);
372: PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom);
373: PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
374: PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
375: PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
376: PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar);
377: PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec);
378: PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode);
379: PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure);
380: PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
381: PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
382: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
383: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
384: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
385: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *);
386: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
387: PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer);
389: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
390: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
391: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
392: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
394: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *);
396: #if defined(PETSC_HAVE_MATLAB)
397: PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *);
398: PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *);
399: #endif
400: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *);
401: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *);
402: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *);
403: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *);
404: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *);
405: #if defined(PETSC_HAVE_SCALAPACK)
406: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
407: #endif
408: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *);
409: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *);
410: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *);
411: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *);
412: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *);
413: PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS);
414: PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *);
415: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
416: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType);
417: PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat);
419: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *);
420: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
421: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
423: PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat);
424: PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqAIJ(Mat A);
425: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *);
426: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
427: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
428: PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]);
429: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *);
431: PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *);
432: PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat);
434: /*
435: PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
437: Input Parameters:
438: + nnz - the number of entries
439: . r - the array of vector values
440: . xv - the matrix values for the row
441: - xi - the column indices of the nonzeros in the row
443: Output Parameter:
444: . sum - negative the sum of results
446: PETSc compile flags:
447: + PETSC_KERNEL_USE_UNROLL_4
448: - PETSC_KERNEL_USE_UNROLL_2
450: Developer Note:
451: The macro changes sum but not other parameters
453: .seealso: `PetscSparseDensePlusDot()`
454: */
455: #if defined(PETSC_KERNEL_USE_UNROLL_4)
456: #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
457: { \
458: if (nnz > 0) { \
459: PetscInt nnz2 = nnz, rem = nnz & 0x3; \
460: switch (rem) { \
461: case 3: \
462: sum -= *xv++ * r[*xi++]; \
463: case 2: \
464: sum -= *xv++ * r[*xi++]; \
465: case 1: \
466: sum -= *xv++ * r[*xi++]; \
467: nnz2 -= rem; \
468: } \
469: while (nnz2 > 0) { \
470: sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
471: xv += 4; \
472: xi += 4; \
473: nnz2 -= 4; \
474: } \
475: xv -= nnz; \
476: xi -= nnz; \
477: } \
478: }
480: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
481: #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
482: { \
483: PetscInt __i, __i1, __i2; \
484: for (__i = 0; __i < nnz - 1; __i += 2) { \
485: __i1 = xi[__i]; \
486: __i2 = xi[__i + 1]; \
487: sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
488: } \
489: if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \
490: }
492: #else
493: #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \
494: { \
495: PetscInt __i; \
496: for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \
497: }
498: #endif
500: /*
501: PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
503: Input Parameters:
504: + nnz - the number of entries
505: . r - the array of vector values
506: . xv - the matrix values for the row
507: - xi - the column indices of the nonzeros in the row
509: Output Parameter:
510: . sum - the sum of results
512: PETSc compile flags:
513: + PETSC_KERNEL_USE_UNROLL_4
514: - PETSC_KERNEL_USE_UNROLL_2
516: Developer Note:
517: The macro changes sum but not other parameters
519: .seealso: `PetscSparseDenseMinusDot()`
520: */
521: #if defined(PETSC_KERNEL_USE_UNROLL_4)
522: #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
523: { \
524: if (nnz > 0) { \
525: PetscInt nnz2 = nnz, rem = nnz & 0x3; \
526: switch (rem) { \
527: case 3: \
528: sum += *xv++ * r[*xi++]; \
529: case 2: \
530: sum += *xv++ * r[*xi++]; \
531: case 1: \
532: sum += *xv++ * r[*xi++]; \
533: nnz2 -= rem; \
534: } \
535: while (nnz2 > 0) { \
536: sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
537: xv += 4; \
538: xi += 4; \
539: nnz2 -= 4; \
540: } \
541: xv -= nnz; \
542: xi -= nnz; \
543: } \
544: }
546: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
547: #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
548: { \
549: PetscInt __i, __i1, __i2; \
550: for (__i = 0; __i < nnz - 1; __i += 2) { \
551: __i1 = xi[__i]; \
552: __i2 = xi[__i + 1]; \
553: sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \
554: } \
555: if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \
556: }
558: #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)
559: #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz))
561: #else
562: #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \
563: { \
564: PetscInt __i; \
565: for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \
566: }
567: #endif
569: #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)
570: #include <immintrin.h>
571: #if !defined(_MM_SCALE_8)
572: #define _MM_SCALE_8 8
573: #endif
575: static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n)
576: {
577: __m512d vec_x, vec_y, vec_vals;
578: __m256i vec_idx;
579: PetscInt j;
581: vec_y = _mm512_setzero_pd();
582: for (j = 0; j < (n >> 3); j++) {
583: vec_idx = _mm256_loadu_si256((__m256i const *)aj);
584: vec_vals = _mm512_loadu_pd(aa);
585: vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
586: vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
587: aj += 8;
588: aa += 8;
589: }
590: #if defined(__AVX512VL__)
591: /* masked load requires avx512vl, which is not supported by KNL */
592: if (n & 0x07) {
593: __mmask8 mask;
594: mask = (__mmask8)(0xff >> (8 - (n & 0x07)));
595: vec_idx = _mm256_mask_loadu_epi32(vec_idx, mask, aj);
596: vec_vals = _mm512_mask_loadu_pd(vec_vals, mask, aa);
597: vec_x = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8);
598: vec_y = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask);
599: }
600: *sum += _mm512_reduce_add_pd(vec_y);
601: #else
602: *sum += _mm512_reduce_add_pd(vec_y);
603: for (j = 0; j < (n & 0x07); j++) *sum += aa[j] * x[aj[j]];
604: #endif
605: }
606: #endif
608: /*
609: PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
611: Input Parameters:
612: + nnz - the number of entries
613: . r - the array of vector values
614: . xv - the matrix values for the row
615: - xi - the column indices of the nonzeros in the row
617: Output Parameter:
618: . max - the max of results
620: .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()`
621: */
622: #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \
623: do { \
624: for (PetscInt __i = 0; __i < (nnz); __i++) { max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); } \
625: } while (0)
627: /*
628: Add column indices into table for counting the max nonzeros of merged rows
629: */
630: #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \
631: do { \
632: if ((mat)) { \
633: for (PetscInt _row = 0; _row < (nrows); _row++) { \
634: const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \
635: for (PetscInt _j = 0; _j < _nz; _j++) { \
636: PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
637: PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \
638: } \
639: } \
640: } \
641: } while (0)
643: /*
644: Add column indices into table for counting the nonzeros of merged rows
645: */
646: #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \
647: do { \
648: for (PetscInt _i = 0; _i < (nrows); _i++) { \
649: const PetscInt _row = (rows)[_i]; \
650: const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \
651: for (PetscInt _j = 0; _j < _nz; _j++) { \
652: PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \
653: PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \
654: } \
655: } \
656: } while (0)
658: #endif