Actual source code: aij.h
petsc-3.9.4 2018-09-11
5: #include <petsc/private/matimpl.h>
6: #include <petscctable.h>
8: /*
9: Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
10: */
11: #define SEQAIJHEADER(datatype) \
12: PetscBool roworiented; /* if true, row-oriented input, default */ \
13: PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \
14: PetscInt nounused; /* -1 generate error on unused space */ \
15: PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \
16: PetscInt maxnz; /* allocated nonzeros */ \
17: PetscInt *imax; /* maximum space allocated for each row */ \
18: PetscInt *ilen; /* actual length of each row */ \
19: PetscInt *ipre; /* space preallocated for each row by user */ \
20: PetscBool free_imax_ilen; \
21: PetscInt reallocs; /* number of mallocs done during MatSetValues() \
22: as more values are set than were prealloced */\
23: PetscInt rmax; /* max nonzeros in any row */ \
24: PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \
25: PetscBool ignorezeroentries; \
26: PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \
27: PetscBool free_a; /* free the numerical values when matrix is destroy */ \
28: Mat_CompressedRow compressedrow; /* use compressed row format */ \
29: PetscInt nz; /* nonzeros */ \
30: PetscInt *i; /* pointer to beginning of each row */ \
31: PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \
32: PetscInt *diag; /* pointers to diagonal elements */ \
33: PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \
34: PetscBool free_diag; \
35: datatype *a; /* nonzero elements */ \
36: PetscScalar *solve_work; /* work space used in MatSolve */ \
37: IS row, col, icol; /* index sets, used for reorderings */ \
38: PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \
39: Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \
40: means that this shares some data structures with the parent including diag, ilen, imax, i, j */\
41: Mat_SubSppt *submatis1 /* used by MatCreateSubMatrices_MPIXAIJ_Local */
43: typedef struct {
44: MatTransposeColoring matcoloring;
45: Mat Bt_den; /* dense matrix of B^T */
46: Mat ABt_den; /* dense matrix of A*B^T */
47: PetscBool usecoloring;
48: PetscErrorCode (*destroy)(Mat);
49: } Mat_MatMatTransMult;
51: typedef struct { /* used by MatTransposeMatMult() */
52: Mat At; /* transpose of the first matrix */
53: Mat mA; /* maij matrix of A */
54: Vec bt,ct; /* vectors to hold locally transposed arrays of B and C */
55: PetscErrorCode (*destroy)(Mat);
56: } Mat_MatTransMatMult;
58: typedef struct {
59: PetscInt *api,*apj; /* symbolic structure of A*P */
60: PetscScalar *apa; /* temporary array for storing one row of A*P */
61: PetscErrorCode (*destroy)(Mat);
62: } Mat_PtAP;
64: typedef struct {
65: MatTransposeColoring matcoloring;
66: Mat Rt; /* sparse or dense matrix of R^T */
67: Mat RARt; /* dense matrix of R*A*R^T */
68: Mat ARt; /* A*R^T used for the case -matrart_color_art */
69: MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
70: PetscErrorCode (*destroy)(Mat);
71: } Mat_RARt;
73: typedef struct {
74: Mat BC; /* temp matrix for storing B*C */
75: PetscErrorCode (*destroy)(Mat);
76: } Mat_MatMatMatMult;
78: /*
79: MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
80: format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example,
81: j[i[k]+p] is the pth column in row k. Note that the diagonal
82: matrix elements are stored with the rest of the nonzeros (not separately).
83: */
85: /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
86: typedef struct {
87: MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
88: PetscInt bdiagsize; /* length of bdiag and ibdiag */
89: PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */
91: PetscBool use;
92: PetscInt node_count; /* number of inodes */
93: PetscInt *size; /* size of each inode */
94: PetscInt limit; /* inode limit */
95: PetscInt max_limit; /* maximum supported inode limit */
96: PetscBool checked; /* if inodes have been checked for */
97: PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */
98: } Mat_SeqAIJ_Inode;
100: PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
101: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
102: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
103: PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
104: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
105: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
106: PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
107: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
108: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
110: typedef struct {
111: SEQAIJHEADER(MatScalar);
112: Mat_SeqAIJ_Inode inode;
113: MatScalar *saved_values; /* location for stashing nonzero values of matrix */
115: PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
116: PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */
117: PetscScalar *ibdiag; /* inverses of block diagonals */
118: PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */
119: PetscScalar fshift,omega; /* last used omega and fshift */
121: ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */
123: PetscScalar *matmult_abdense; /* used by MatMatMult() */
124: Mat_PtAP *ptap; /* used by MatPtAP() */
125: Mat_MatMatMatMult *matmatmatmult; /* used by MatMatMatMult() */
126: Mat_RARt *rart; /* used by MatRARt() */
127: Mat_MatMatTransMult *abt; /* used by MatMatTransposeMult() */
128: Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult() */
129: } Mat_SeqAIJ;
131: /*
132: Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
133: */
134: PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
135: {
137: Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data;
138: if (A->singlemalloc) {
139: PetscFree3(*a,*j,*i);
140: } else {
141: if (A->free_a) {PetscFree(*a);}
142: if (A->free_ij) {PetscFree(*j);}
143: if (A->free_ij) {PetscFree(*i);}
144: }
145: return 0;
146: }
147: /*
148: Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
149: This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
150: */
151: #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
152: if (NROW >= RMAX) { \
153: Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
154: /* there is no extra room in row, therefore enlarge */ \
155: PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
156: datatype *new_a; \
157: \
158: if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \
159: /* malloc new storage space */ \
160: PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i); \
161: \
162: /* copy over old data into new slots */ \
163: for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
164: for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
165: PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt)); \
166: len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
167: PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt)); \
168: PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype)); \
169: PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype)); \
170: PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype)); \
171: /* free up old matrix storage */ \
172: MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i); \
173: AA = new_a; \
174: Ain->a = (MatScalar*) new_a; \
175: AI = Ain->i = new_i; AJ = Ain->j = new_j; \
176: Ain->singlemalloc = PETSC_TRUE; \
177: \
178: RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
179: RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
180: Ain->maxnz += BS2*CHUNKSIZE; \
181: Ain->reallocs++; \
182: } \
184: #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \
185: if (NROW >= RMAX) { \
186: Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
187: /* there is no extra room in row, therefore enlarge */ \
188: PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
189: \
190: if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \
191: /* malloc new storage space */ \
192: PetscMalloc1(new_nz,&new_j); \
193: PetscMalloc1(AM+1,&new_i);\
194: \
195: /* copy over old data into new slots */ \
196: for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
197: for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
198: PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt)); \
199: len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
200: PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt)); \
201: \
202: /* free up old matrix storage */ \
203: MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i); \
204: Ain->a = NULL; \
205: AI = Ain->i = new_i; AJ = Ain->j = new_j; \
206: Ain->singlemalloc = PETSC_FALSE; \
207: Ain->free_a = PETSC_FALSE; \
208: \
209: RP = AJ + AI[ROW]; \
210: RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
211: Ain->maxnz += BS2*CHUNKSIZE; \
212: Ain->reallocs++; \
213: } \
215: PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
216: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
217: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
218: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
220: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
221: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
222: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
223: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
224: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
225: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
226: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
227: PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
228: PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
229: PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
230: PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
232: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
233: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
234: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
235: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
236: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
238: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
240: PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
241: PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
242: PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
243: PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
244: PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
245: PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
246: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
247: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
248: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
249: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
250: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
251: PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
252: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
253: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
254: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
255: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
256: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
257: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
258: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
259: PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
260: PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
261: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
262: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
263: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
264: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
265: PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
266: PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
267: PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
268: PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
269: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
270: PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
271: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
272: PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
274: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
275: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
276: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat*);
277: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
278: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
279: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
280: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
281: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
282: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
283: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
285: PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
286: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
287: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
288: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
290: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
291: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
292: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
293: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
294: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
295: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
297: PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
298: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
299: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
300: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat);
302: PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
303: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
304: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
306: PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
307: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
308: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
309: PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
310: PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
311: PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
313: PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
314: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
315: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
317: PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
318: PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
319: PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
320: PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
321: PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
322: PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
323: PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
324: PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
325: PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
326: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
327: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
328: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
329: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
330: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
331: PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
332: PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
334: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
335: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
336: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
337: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
339: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
341: PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
342: PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
343: PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
344: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
345: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
346: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
347: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
348: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
349: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
350: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
351: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
352: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
353: PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
354: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
355: PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
356: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
357: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
359: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
360: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
361: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
363: PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
364: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
365: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
366: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
367: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
369: /*
370: PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
372: Input Parameters:
373: + nnz - the number of entries
374: . r - the array of vector values
375: . xv - the matrix values for the row
376: - xi - the column indices of the nonzeros in the row
378: Output Parameter:
379: . sum - negative the sum of results
381: PETSc compile flags:
382: + PETSC_KERNEL_USE_UNROLL_4
383: - PETSC_KERNEL_USE_UNROLL_2
385: Developer Notes:
386: The macro changes sum but not other parameters
388: .seealso: PetscSparseDensePlusDot()
390: */
391: #if defined(PETSC_KERNEL_USE_UNROLL_4)
392: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
393: if (nnz > 0) { \
394: PetscInt nnz2=nnz,rem=nnz&0x3; \
395: switch (rem) { \
396: case 3: sum -= *xv++ *r[*xi++]; \
397: case 2: sum -= *xv++ *r[*xi++]; \
398: case 1: sum -= *xv++ *r[*xi++]; \
399: nnz2 -= rem;} \
400: while (nnz2 > 0) { \
401: sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
402: xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
403: xv += 4; xi += 4; nnz2 -= 4; \
404: } \
405: xv -= nnz; xi -= nnz; \
406: } \
407: }
409: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
410: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
411: PetscInt __i,__i1,__i2; \
412: for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
413: sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
414: if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
416: #else
417: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
418: PetscInt __i; \
419: for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
420: #endif
424: /*
425: PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
427: Input Parameters:
428: + nnz - the number of entries
429: . r - the array of vector values
430: . xv - the matrix values for the row
431: - xi - the column indices of the nonzeros in the row
433: Output Parameter:
434: . sum - the sum of results
436: PETSc compile flags:
437: + PETSC_KERNEL_USE_UNROLL_4
438: - PETSC_KERNEL_USE_UNROLL_2
440: Developer Notes:
441: The macro changes sum but not other parameters
443: .seealso: PetscSparseDenseMinusDot()
445: */
446: #if defined(PETSC_KERNEL_USE_UNROLL_4)
447: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
448: if (nnz > 0) { \
449: PetscInt nnz2=nnz,rem=nnz&0x3; \
450: switch (rem) { \
451: case 3: sum += *xv++ *r[*xi++]; \
452: case 2: sum += *xv++ *r[*xi++]; \
453: case 1: sum += *xv++ *r[*xi++]; \
454: nnz2 -= rem;} \
455: while (nnz2 > 0) { \
456: sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
457: xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
458: xv += 4; xi += 4; nnz2 -= 4; \
459: } \
460: xv -= nnz; xi -= nnz; \
461: } \
462: }
464: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
465: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
466: PetscInt __i,__i1,__i2; \
467: for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
468: sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
469: if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
471: #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)
472: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))
474: #else
475: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
476: PetscInt __i; \
477: for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
478: #endif
480: #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)
481: #include <immintrin.h>
482: #if !defined(_MM_SCALE_8)
483: #define _MM_SCALE_8 8
484: #endif
486: PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
487: {
488: __m512d vec_x,vec_y,vec_vals;
489: __m256i vec_idx;
490: __mmask8 mask;
491: PetscInt j;
493: vec_y = _mm512_setzero_pd();
494: for (j=0; j<(n>>3); j++) {
495: vec_idx = _mm256_loadu_si256((__m256i const*)aj);
496: vec_vals = _mm512_loadu_pd(aa);
497: vec_x = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8);
498: vec_y = _mm512_fmadd_pd(vec_x,vec_vals,vec_y);
499: aj += 8; aa += 8;
500: }
501: /* masked load does not work on KNL, it requires avx512vl */
502: if ((n&0x07)>2) {
503: mask = (__mmask8)(0xff >> (8-(n&0x07)));
504: vec_idx = _mm256_loadu_si256((__m256i const*)aj);
505: vec_vals = _mm512_loadu_pd(aa);
506: vec_x = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8);
507: vec_y = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask);
508: } else if ((n&0x07)==2) {
509: *sum += aa[0]*x[aj[0]];
510: *sum += aa[1]*x[aj[1]];
511: } else if ((n&0x07)==1) {
512: *sum += aa[0]*x[aj[0]];
513: }
514: if (n>2) *sum += _mm512_reduce_add_pd(vec_y);
515: /*
516: for(j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
517: */
518: }
519: #endif
521: /*
522: PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
524: Input Parameters:
525: + nnz - the number of entries
526: . r - the array of vector values
527: . xv - the matrix values for the row
528: - xi - the column indices of the nonzeros in the row
530: Output Parameter:
531: . max - the max of results
533: .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
535: */
536: #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
537: PetscInt __i; \
538: for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
540: /*
541: Add column indices into table for counting the max nonzeros of merged rows
542: */
543: #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) { \
544: PetscInt _j,_row,_nz,*_col; \
545: if (mat) { \
546: for (_row=0; _row<nrows; _row++) { \
547: _nz = mat->i[_row+1] - mat->i[_row]; \
548: for (_j=0; _j<_nz; _j++) { \
549: _col = _j + mat->j + mat->i[_row]; \
550: PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
551: } \
552: } \
553: } \
554: }
556: /*
557: Add column indices into table for counting the nonzeros of merged rows
558: */
559: #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) { \
560: PetscInt _j,_row,_nz,*_col,_i; \
561: for (_i=0; _i<nrows; _i++) {\
562: _row = rows[_i]; \
563: _nz = mat->i[_row+1] - mat->i[_row]; \
564: for (_j=0; _j<_nz; _j++) { \
565: _col = _j + mat->j + mat->i[_row]; \
566: PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
567: } \
568: } \
569: }
571: #endif