Actual source code: aij.h
5: #include <private/matimpl.h>
7: /*
8: Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
9: */
10: #define SEQAIJHEADER(datatype) \
11: PetscBool roworiented; /* if true, row-oriented input, default */\
12: PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */\
13: PetscInt nounused; /* -1 generate error on unused space */\
14: PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\
15: PetscInt maxnz; /* allocated nonzeros */\
16: PetscInt *imax; /* maximum space allocated for each row */\
17: PetscInt *ilen; /* actual length of each row */\
18: PetscBool free_imax_ilen; \
19: PetscInt reallocs; /* number of mallocs done during MatSetValues() \
20: as more values are set than were prealloced */\
21: PetscInt rmax; /* max nonzeros in any row */\
22: PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/\
23: PetscBool ignorezeroentries; \
24: PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\
25: Mat XtoY; /* used by MatAXPY() */\
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: PetscBool free_diag; \
34: datatype *a; /* nonzero elements */ \
35: PetscScalar *solve_work; /* work space used in MatSolve */ \
36: IS row, col, icol; /* index sets, used for reorderings */ \
37: PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \
38: Mat parent /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....);
39: means that this shares some data structures with the parent including diag, ilen, imax, i, j */
41: /*
42: MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
43: format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example,
44: j[i[k]+p] is the pth column in row k. Note that the diagonal
45: matrix elements are stored with the rest of the nonzeros (not separately).
46: */
48: /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
49: typedef struct {
50: MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
51: PetscInt bdiagsize; /* length of bdiag and ibdiag */
52: PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */
54: PetscBool use;
55: PetscInt node_count; /* number of inodes */
56: PetscInt *size; /* size of each inode */
57: PetscInt limit; /* inode limit */
58: PetscInt max_limit; /* maximum supported inode limit */
59: PetscBool checked; /* if inodes have been checked for */
60: } Mat_SeqAIJ_Inode;
72: typedef struct {
73: SEQAIJHEADER(MatScalar);
74: Mat_SeqAIJ_Inode inode;
75: MatScalar *saved_values; /* location for stashing nonzero values of matrix */
77: PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
78: PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */
79: PetscScalar *ibdiag; /* inverses of block diagonals */
80: PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */
81: PetscScalar fshift,omega; /* last used omega and fshift */
83: ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */
84: } Mat_SeqAIJ;
86: /*
87: Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
88: */
91: PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
92: {
94: Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data;
95: if (A->singlemalloc) {
96: PetscFree3(*a,*j,*i);
97: } else {
98: if (A->free_a) {PetscFree(*a);}
99: if (A->free_ij) {PetscFree(*j);}
100: if (A->free_ij) {PetscFree(*i);}
101: }
102: return 0;
103: }
104: /*
105: Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
106: This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
107: */
108: #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
109: if (NROW >= RMAX) {\
110: Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
111: /* there is no extra room in row, therefore enlarge */ \
112: PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
113: datatype *new_a; \
114: \
115: if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
116: /* malloc new storage space */ \
117: PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);\
118: \
119: /* copy over old data into new slots */ \
120: for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
121: for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
122: PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt)); \
123: len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
124: PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt)); \
125: PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype)); \
126: PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));\
127: PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype)); \
128: /* free up old matrix storage */ \
129: MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);\
130: AA = new_a; \
131: Ain->a = (MatScalar*) new_a; \
132: AI = Ain->i = new_i; AJ = Ain->j = new_j; \
133: Ain->singlemalloc = PETSC_TRUE; \
134: \
135: RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
136: RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
137: Ain->maxnz += BS2*CHUNKSIZE; \
138: Ain->reallocs++; \
139: } \
236: /*
237: PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
239: Input Parameters:
240: + nnz - the number of entries
241: . r - the array of vector values
242: . xv - the matrix values for the row
243: - xi - the column indices of the nonzeros in the row
245: Output Parameter:
246: . sum - negative the sum of results
248: PETSc compile flags:
249: + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG
250: - PETSC_KERNEL_USE_UNROLL_2 -
252: .seealso: PetscSparseDensePlusDot()
254: */
255: #ifdef PETSC_KERNEL_USE_UNROLL_4
256: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
257: if (nnz > 0) {\
258: switch (nnz & 0x3) {\
259: case 3: sum -= *xv++ * r[*xi++];\
260: case 2: sum -= *xv++ * r[*xi++];\
261: case 1: sum -= *xv++ * r[*xi++];\
262: nnz -= 4;}\
263: while (nnz > 0) {\
264: sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\
265: xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\
266: xv += 4; xi += 4; nnz -= 4; }}}
268: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
269: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
270: PetscInt __i,__i1,__i2;\
271: for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
272: sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
273: if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
275: #else
276: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
277: PetscInt __i;\
278: for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];}
279: #endif
283: /*
284: PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
286: Input Parameters:
287: + nnz - the number of entries
288: . r - the array of vector values
289: . xv - the matrix values for the row
290: - xi - the column indices of the nonzeros in the row
292: Output Parameter:
293: . sum - the sum of results
295: PETSc compile flags:
296: + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG
297: - PETSC_KERNEL_USE_UNROLL_2 -
299: .seealso: PetscSparseDenseMinusDot()
301: */
302: #ifdef PETSC_KERNEL_USE_UNROLL_4
303: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
304: if (nnz > 0) {\
305: switch (nnz & 0x3) {\
306: case 3: sum += *xv++ * r[*xi++];\
307: case 2: sum += *xv++ * r[*xi++];\
308: case 1: sum += *xv++ * r[*xi++];\
309: nnz -= 4;}\
310: while (nnz > 0) {\
311: sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\
312: xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\
313: xv += 4; xi += 4; nnz -= 4; }}}
315: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
316: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
317: PetscInt __i,__i1,__i2;\
318: for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
319: sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
320: if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
322: #else
323: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
324: PetscInt __i;\
325: for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];}
326: #endif
328: #endif