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