xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 71fd2e92b2fee25f70818ffc1dca9fc57f8584d7)
1 
2 #if !defined(__AIJ_H)
3 #define __AIJ_H
4 
5 #include "private/matimpl.h"
6 
7 /*
8     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
9 */
10 #define SEQAIJHEADER(datatype)	\
11   PetscTruth        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   PetscTruth        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   PetscTruth        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   PetscTruth        keepnonzeropattern;   /* keeps matrix structure same in calls to MatZeroRows()*/\
23   PetscTruth        ignorezeroentries; \
24   PetscInt          *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */\
25   Mat               XtoY;             /* used by MatAXPY() */\
26   PetscTruth        free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
27   PetscTruth        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   PetscTruth        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   PetscTruth        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 */
40 
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 */
47 
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   PetscTruth  ibdiagvalid;                     /* do ibdiag[] and bdiag[] contain the most recent values */
53 
54   PetscTruth 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   PetscTruth checked;                       /* if inodes have been checked for */
60 } Mat_SeqAIJ_Inode;
61 
62 EXTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
63 EXTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
64 EXTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
65 EXTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
66 EXTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscTruth);
67 EXTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
68 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ_Inode(Mat,IS,IS,const MatFactorInfo*,Mat*);
69 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_Inode(Mat,Mat,IS,IS,const MatFactorInfo*);
70 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_Inode(Mat,Mat,IS,IS,const MatFactorInfo*);
71 
72 typedef struct {
73   SEQAIJHEADER(MatScalar);
74   Mat_SeqAIJ_Inode inode;
75   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
76 
77   PetscScalar      *idiag,*mdiag,*ssor_work;  /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
78   PetscTruth       idiagvalid;                     /* current idiag[] and mdiag[] are valid */
79   PetscScalar      fshift,omega;                   /* last used omega and fshift */
80 
81   ISColoring       coloring;                  /* set with MatADSetColoring() used by MatADSetValues() */
82 } Mat_SeqAIJ;
83 
84 /*
85     Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
86 */
87 #undef __FUNCT__
88 #define __FUNCT__ "MatSeqXAIJFreeAIJ"
89 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
90 {
91                                      PetscErrorCode ierr;
92                                      Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
93                                      if (A->singlemalloc) {
94                                        ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
95                                      } else {
96                                        if (A->free_a  && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);}
97                                        if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);}
98                                        if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);}
99                                      }
100                                      *a = 0; *j = 0; *i = 0;
101                                      return 0;
102                                    }
103 
104 /*
105     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
106 */
107 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
108   if (NROW >= RMAX) {\
109 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
110         /* there is no extra room in row, therefore enlarge */ \
111         PetscInt   new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
112         datatype   *new_a; \
113  \
114         if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
115         /* malloc new storage space */ \
116         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
117  \
118         /* copy over old data into new slots */ \
119         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
120         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
121         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
122         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
123         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
124         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
125         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
126         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
127         /* free up old matrix storage */ \
128         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
129         AA = new_a; \
130         Ain->a = (MatScalar*) new_a;		   \
131         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
132         Ain->singlemalloc = PETSC_TRUE; \
133  \
134         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
135         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
136         Ain->maxnz += CHUNKSIZE; \
137         Ain->reallocs++; \
138       } \
139 
140 EXTERN_C_BEGIN
141 EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*);
142 EXTERN_C_END
143 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
144 EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
145 EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
146 EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
147 EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
148 EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*);
149 EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
150 
151 EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
152 EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
153 EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
154 EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
155 EXTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
156 
157 EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
158 EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
159 EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
160 
161 EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
162 EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
163 EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
164 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
165 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
166 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_newdatastruct(Mat,Mat,IS,IS,const MatFactorInfo*);
167 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
168 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_newdatastruct(Mat,Mat,const MatFactorInfo*);
169 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
170 EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
171 EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
172 EXTERN PetscErrorCode MatSolve_SeqAIJ_newdatastruct(Mat,Vec,Vec);
173 EXTERN PetscErrorCode MatSolve_SeqAIJ_newdatastruct_v2(Mat,Vec,Vec);
174 EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
175 EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_newdatastruct(Mat,Vec,Vec);
176 EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_newdatastruct_v2(Mat,Vec,Vec);
177 EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
178 EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
179 EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
180 EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
181 EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
182 EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
183 EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
184 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*,Mat*);
185 EXTERN PetscErrorCode MatILUDTFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
186 EXTERN PetscErrorCode MatILUDTFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
187 EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*);
188 EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
189 EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
190 EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
191 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
192 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
193 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
194 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
195 EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
196 EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
197 EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
198 EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
199 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
200 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
201 EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
202 EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
203 EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
204 EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
205 EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
206 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
207 EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
208 
209 EXTERN_C_BEGIN
210 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*);
211 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*);
212 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*);
213 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
214 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
215 EXTERN_C_END
216 
217 /*
218     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
219 
220   Input Parameters:
221 +  nnz - the number of entries
222 .  r - the array of vector values
223 .  xv - the matrix values for the row
224 -  xi - the column indices of the nonzeros in the row
225 
226   Output Parameter:
227 .  sum - negative the sum of results
228 
229   PETSc compile flags:
230 +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
231 -   PETSC_KERNEL_USE_UNROLL_2 -
232 
233 .seealso: PetscSparseDensePlusDot()
234 
235 */
236 #ifdef PETSC_KERNEL_USE_UNROLL_4
237 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
238 if (nnz > 0) {\
239 switch (nnz & 0x3) {\
240 case 3: sum -= *xv++ * r[*xi++];\
241 case 2: sum -= *xv++ * r[*xi++];\
242 case 1: sum -= *xv++ * r[*xi++];\
243 nnz -= 4;}\
244 while (nnz > 0) {\
245 sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\
246 	xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\
247 xv  += 4; xi += 4; nnz -= 4; }}}
248 
249 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
250 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
251 PetscInt __i,__i1,__i2;\
252 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
253 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
254 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
255 
256 #else
257 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
258 PetscInt __i;\
259 for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];}
260 #endif
261 
262 
263 
264 /*
265     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
266 
267   Input Parameters:
268 +  nnz - the number of entries
269 .  r - the array of vector values
270 .  xv - the matrix values for the row
271 -  xi - the column indices of the nonzeros in the row
272 
273   Output Parameter:
274 .  sum - the sum of results
275 
276   PETSc compile flags:
277 +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
278 -   PETSC_KERNEL_USE_UNROLL_2 -
279 
280 .seealso: PetscSparseDenseMinusDot()
281 
282 */
283 #ifdef PETSC_KERNEL_USE_UNROLL_4
284 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
285 if (nnz > 0) {\
286 switch (nnz & 0x3) {\
287 case 3: sum += *xv++ * r[*xi++];\
288 case 2: sum += *xv++ * r[*xi++];\
289 case 1: sum += *xv++ * r[*xi++];\
290 nnz -= 4;}\
291 while (nnz > 0) {\
292 sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\
293 	xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\
294 xv  += 4; xi += 4; nnz -= 4; }}}
295 
296 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
297 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
298 PetscInt __i,__i1,__i2;\
299 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
300 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
301 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
302 
303 #else
304 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
305  PetscInt __i;\
306 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];}
307 #endif
308 
309 #endif
310