xref: /petsc/src/mat/impls/aij/seq/aij.h (revision bcaeba4d41d6ca6c6dc4189db20683073a9959ce)
1 
2 #if !defined(__AIJ_H)
3 #define __AIJ_H
4 
5 #include <petsc-private/matimpl.h>
6 
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 */
40 
41 typedef struct {
42   MatTransposeColoring      matcoloring;
43   Mat                       Bt_den;  /* dense matrix of B^T */
44   Mat                       ABt_den; /* dense matrix of A*B^T */
45   PetscBool                 usecoloring;
46   PetscErrorCode (*destroy)(Mat);
47 } Mat_MatMatTransMult;
48 
49 typedef struct {
50   PetscInt       *api,*apj;    /* symbolic structure of A*P */
51   PetscScalar    *apa;         /* temporary array for storing one row of A*P */
52   PetscErrorCode (*destroy)(Mat);
53 } Mat_PtAP;
54 
55 typedef struct {
56   MatTransposeColoring matcoloring;
57   Mat                  Rt;    /* dense matrix of R^T */
58   Mat                  RARt;  /* dense matrix of R*A*R^T */
59   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
60   PetscErrorCode (*destroy)(Mat);
61 } Mat_RARt;
62 
63 typedef struct {
64   Mat            BC;    /* temp matrix for storing B*C */
65   PetscErrorCode (*destroy)(Mat);
66 } Mat_MatMatMatMult;
67 
68 /*
69   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
70   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
71   j[i[k]+p] is the pth column in row k.  Note that the diagonal
72   matrix elements are stored with the rest of the nonzeros (not separately).
73 */
74 
75 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
76 typedef struct {
77   MatScalar   *bdiag,*ibdiag,*ssor_work;      /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
78   PetscInt    bdiagsize;                       /* length of bdiag and ibdiag */
79   PetscBool   ibdiagvalid;                     /* do ibdiag[] and bdiag[] contain the most recent values */
80 
81   PetscBool  use;
82   PetscInt   node_count;                    /* number of inodes */
83   PetscInt   *size;                         /* size of each inode */
84   PetscInt   limit;                         /* inode limit */
85   PetscInt   max_limit;                     /* maximum supported inode limit */
86   PetscBool  checked;                       /* if inodes have been checked for */
87 } Mat_SeqAIJ_Inode;
88 
89 extern PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
90 extern PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
91 extern PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
92 extern PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
93 extern PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool );
94 extern PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
95 extern PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool );
96 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
97 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
98 
99 typedef struct {
100   SEQAIJHEADER(MatScalar);
101   Mat_SeqAIJ_Inode inode;
102   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
103 
104   PetscScalar      *idiag,*mdiag,*ssor_work;  /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
105   PetscBool        idiagvalid;                     /* current idiag[] and mdiag[] are valid */
106   PetscScalar      *ibdiag;                   /* inverses of block diagonals */
107   PetscBool        ibdiagvalid;               /* inverses of block diagonals are valid. */
108   PetscScalar      fshift,omega;                   /* last used omega and fshift */
109 
110   ISColoring       coloring;                  /* set with MatADSetColoring() used by MatADSetValues() */
111 
112   PetscScalar      *matmult_abdense;     /* used by MatMatMult() */
113   Mat_PtAP         *ptap;                /* used by MatPtAP() */
114   Mat_MatMatMatMult *matmatmatmult;      /* used by MatMatMatMult() */
115 } Mat_SeqAIJ;
116 
117 /*
118   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
119 */
120 #undef __FUNCT__
121 #define __FUNCT__ "MatSeqXAIJFreeAIJ"
122 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
123 {
124   PetscErrorCode ierr;
125   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
126   if (A->singlemalloc) {
127     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
128   } else {
129     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
130     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
131     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
132   }
133   return 0;
134 }
135 /*
136     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
137     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
138 */
139 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
140   if (NROW >= RMAX) {\
141 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
142         /* there is no extra room in row, therefore enlarge */ \
143         PetscInt   CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
144         datatype   *new_a; \
145  \
146         if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
147         /* malloc new storage space */ \
148         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
149  \
150         /* copy over old data into new slots */ \
151         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
152         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
153         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
154         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
155         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
156         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
157         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
158         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
159         /* free up old matrix storage */ \
160         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
161         AA = new_a; \
162         Ain->a = (MatScalar*) new_a;		   \
163         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
164         Ain->singlemalloc = PETSC_TRUE; \
165  \
166         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
167         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
168         Ain->maxnz += BS2*CHUNKSIZE; \
169         Ain->reallocs++; \
170       } \
171 
172 
173 EXTERN_C_BEGIN
174 extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
175 EXTERN_C_END
176 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
177 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
178 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
179 
180 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
181 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
182 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
183 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
184 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
185 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
186 extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
187 extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
188 extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*);
189 extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
190 extern PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
191 
192 extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
193 extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
194 extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
195 extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
196 extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
197 
198 extern PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
199 extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
200 extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
201 
202 extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
203 extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
204 extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
205 extern PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
206 extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
207 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
208 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
209 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
210 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
211 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
212 extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
213 extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
214 extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
215 extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
216 extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
217 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
218 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
219 extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
220 extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
221 extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
222 extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
223 extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
224 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
225 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
226 extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
227 extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
228 extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg);
229 extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
230 extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
231 extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
232 
233 extern PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
234 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
235 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
236 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
237 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
238 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
239 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
240 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
241 
242 extern PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
243 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
244 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
245 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
246 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
247 
248 extern PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
249 extern PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
250 
251 extern PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
252 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
253 extern PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
254 extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
255 extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
256 extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
257 extern PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
258 extern PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
259 extern PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
260 
261 extern PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
262 extern PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
263 extern PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
264 
265 extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
266 extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
267 extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
268 extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
269 extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
270 extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
271 extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
272 extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
273 extern PetscErrorCode MatDestroy_SeqAIJ(Mat);
274 extern PetscErrorCode MatSetUp_SeqAIJ(Mat);
275 extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
276 
277 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool );
278 extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool );
279 
280 extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
281 
282 EXTERN_C_BEGIN
283 extern PetscErrorCode  MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
284 extern PetscErrorCode  MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
285 extern PetscErrorCode  MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
286 extern PetscErrorCode  MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
287 extern PetscErrorCode  MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
288 extern PetscErrorCode  MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
289 extern PetscErrorCode  MatCreate_SeqAIJ(Mat);
290 EXTERN_C_END
291 extern PetscErrorCode  MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
292 extern PetscErrorCode  MatDestroy_SeqAIJ(Mat);
293 
294 
295 /*
296     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
297 
298   Input Parameters:
299 +  nnz - the number of entries
300 .  r - the array of vector values
301 .  xv - the matrix values for the row
302 -  xi - the column indices of the nonzeros in the row
303 
304   Output Parameter:
305 .  sum - negative the sum of results
306 
307   PETSc compile flags:
308 +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
309 -   PETSC_KERNEL_USE_UNROLL_2 -
310 
311 .seealso: PetscSparseDensePlusDot()
312 
313 */
314 #ifdef PETSC_KERNEL_USE_UNROLL_4
315 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
316 if (nnz > 0) {\
317 switch (nnz & 0x3) {\
318 case 3: sum -= *xv++ * r[*xi++];\
319 case 2: sum -= *xv++ * r[*xi++];\
320 case 1: sum -= *xv++ * r[*xi++];\
321 nnz -= 4;}\
322 while (nnz > 0) {\
323 sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\
324 	xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\
325 xv  += 4; xi += 4; nnz -= 4; }}}
326 
327 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
328 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
329 PetscInt __i,__i1,__i2;\
330 for (__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
331 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
332 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
333 
334 #else
335 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
336 PetscInt __i;\
337 for (__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];}
338 #endif
339 
340 
341 
342 /*
343     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
344 
345   Input Parameters:
346 +  nnz - the number of entries
347 .  r - the array of vector values
348 .  xv - the matrix values for the row
349 -  xi - the column indices of the nonzeros in the row
350 
351   Output Parameter:
352 .  sum - the sum of results
353 
354   PETSc compile flags:
355 +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
356 -   PETSC_KERNEL_USE_UNROLL_2 -
357 
358 .seealso: PetscSparseDenseMinusDot()
359 
360 */
361 #ifdef PETSC_KERNEL_USE_UNROLL_4
362 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
363 if (nnz > 0) {\
364 switch (nnz & 0x3) {\
365 case 3: sum += *xv++ * r[*xi++];\
366 case 2: sum += *xv++ * r[*xi++];\
367 case 1: sum += *xv++ * r[*xi++];\
368 nnz -= 4;}\
369 while (nnz > 0) {\
370 sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\
371 	xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\
372 xv  += 4; xi += 4; nnz -= 4; }}}
373 
374 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
375 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
376 PetscInt __i,__i1,__i2;\
377 for (__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
378 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
379 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
380 
381 #else
382 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
383  PetscInt __i;\
384 for (__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];}
385 #endif
386 
387 #endif
388