Actual source code: aij.h
5: #include include/private/matimpl.h
6: /*
7: Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
8: */
9: #define SEQAIJHEADER(datatype) \
10: PetscTruth sorted; /* if true, rows are sorted by increasing columns */\
11: PetscTruth roworiented; /* if true, row-oriented input, default */\
12: PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */\
13: PetscTruth singlemalloc; /* if true a, i, and j have been obtained with one big malloc */\
14: PetscInt maxnz; /* allocated nonzeros */\
15: PetscInt *imax; /* maximum space allocated for each row */\
16: PetscInt *ilen; /* actual length of each row */\
17: PetscInt reallocs; /* number of mallocs done during MatSetValues() \
18: as more values are set than were prealloced */\
19: PetscInt rmax; /* max nonzeros in any row */\
20: PetscTruth keepzeroedrows; /* keeps matrix structure same in calls to MatZeroRows()*/\
21: PetscTruth ignorezeroentries; \
22: PetscInt *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */\
23: Mat XtoY; /* used by MatAXPY() */\
24: PetscTruth free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \
25: PetscTruth free_a; /* free the numerical values when matrix is destroy */ \
26: Mat_CompressedRow compressedrow; /* use compressed row format */ \
27: PetscInt nz; /* nonzeros */ \
28: PetscInt *i; /* pointer to beginning of each row */ \
29: PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \
30: PetscInt *diag; /* pointers to diagonal elements */ \
31: datatype *a; /* nonzero elements */ \
32: PetscScalar *solve_work; /* work space used in MatSolve */ \
33: IS row, col, icol /* index sets, used for reorderings */
35: /*
36: MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
37: format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example,
38: j[i[k]+p] is the pth column in row k. Note that the diagonal
39: matrix elements are stored with the rest of the nonzeros (not separately).
40: */
42: /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
43: typedef struct {
44: PetscTruth use;
45: PetscInt node_count; /* number of inodes */
46: PetscInt *size; /* size of each inode */
47: PetscInt limit; /* inode limit */
48: PetscInt max_limit; /* maximum supported inode limit */
49: PetscTruth checked; /* if inodes have been checked for */
50: } Mat_Inode;
52: EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer);
53: EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType);
54: EXTERN PetscErrorCode MatDestroy_Inode(Mat);
55: EXTERN PetscErrorCode MatCreate_Inode(Mat);
56: EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption);
57: EXTERN PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *B);
58: EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
59: EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
60: EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
63: typedef struct {
64: SEQAIJHEADER(PetscScalar);
65: Mat_Inode inode;
66: PetscScalar *saved_values; /* location for stashing nonzero values of matrix */
67: PetscScalar *idiag,*ssor; /* inverse of diagonal entries; space for eisen */
68: ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */
69: } Mat_SeqAIJ;
71: /*
72: Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
73: */
76: PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,PetscScalar **a,PetscInt **j,PetscInt **i) {
78: Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data;
79: if (A->singlemalloc) {
80: PetscFree3(*a,*j,*i);
81: } else {
82: if (A->free_a && *a) {PetscFree(*a);}
83: if (A->free_ij && *j) {PetscFree(*j);}
84: if (A->free_ij && *i) {PetscFree(*i);}
85: }
86: *a = 0; *j = 0; *i = 0;
87: return 0;
88: }
90: /*
91: Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
92: */
93: #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
94: if (NROW >= RMAX) {\
95: Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
96: /* there is no extra room in row, therefore enlarge */ \
97: PetscInt new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
98: datatype *new_a; \
99: \
100: if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
101: /* malloc new storage space */ \
102: PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);\
103: \
104: /* copy over old data into new slots */ \
105: for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
106: for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
107: PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt)); \
108: len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
109: PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt)); \
110: PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype)); \
111: PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));\
112: PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype)); \
113: /* free up old matrix storage */ \
114: MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);\
115: AA = new_a; \
116: Ain->a = (MatScalar*) new_a; \
117: AI = Ain->i = new_i; AJ = Ain->j = new_j; \
118: Ain->singlemalloc = PETSC_TRUE; \
119: \
120: RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
121: RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
122: Ain->maxnz += CHUNKSIZE; \
123: Ain->reallocs++; \
124: } \
127: EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*);
129: EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
130: EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat *);
131: EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat*);
132: EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*);
133: EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
134: EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*);
135: EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
137: EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
138: EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
139: EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
140: EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
141: EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
143: EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
144: EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
145: EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
147: EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
148: EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
149: EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
150: EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
151: EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
152: EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*);
153: EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,MatFactorInfo*,Mat*);
154: EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*);
155: EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
156: EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
157: EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
158: EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
159: EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
160: EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
161: EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
162: EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
163: EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
164: EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, MatType,Mat*);
165: EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
166: EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
167: EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
168: EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
169: EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
170: EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
171: EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
172: EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
173: EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
174: EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
175: EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
176: EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
177: EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
178: EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
179: EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
182: EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType,MatReuse,Mat*);
183: EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
184: EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqCSRPERM(Mat,MatType,MatReuse,Mat*);
185: EXTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
188: #endif