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