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