xref: /petsc/src/mat/impls/aij/seq/aij.h (revision e4d1774bba497c88855ca48ceec8de09b31e2d79)
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 /*
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   PetscBool   ibdiagvalid;                     /* do ibdiag[] and bdiag[] contain the most recent values */
53 
54   PetscBool  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   PetscBool  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,PetscBool );
67 extern PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
68 extern PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool );
69 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
70 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,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   PetscBool        idiagvalid;                     /* current idiag[] and mdiag[] are valid */
79   PetscScalar      *ibdiag;                   /* inverses of block diagonals */
80   PetscBool        ibdiagvalid;               /* inverses of block diagonals are valid. */
81   PetscScalar      fshift,omega;                   /* last used omega and fshift */
82 
83   ISColoring       coloring;                  /* set with MatADSetColoring() used by MatADSetValues() */
84 } Mat_SeqAIJ;
85 
86 typedef struct {
87   MatTransposeColoring      matcoloring;
88   Mat                       Bt_den;  /* dense matrix of B^T */
89   Mat                       ABt_den; /* dense matrix of A*B^T */
90   PetscBool                 usecoloring;
91   PetscErrorCode (*destroy)(Mat);
92 } Mat_MatMatMultTrans;
93 
94 /*
95   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
96 */
97 #undef __FUNCT__
98 #define __FUNCT__ "MatSeqXAIJFreeAIJ"
99 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
100 {
101   PetscErrorCode ierr;
102   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
103   if (A->singlemalloc) {
104     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
105   } else {
106     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
107     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
108     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
109   }
110   return 0;
111 }
112 /*
113     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
114     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
115 */
116 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
117   if (NROW >= RMAX) {\
118 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
119         /* there is no extra room in row, therefore enlarge */ \
120         PetscInt   CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
121         datatype   *new_a; \
122  \
123         if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
124         /* malloc new storage space */ \
125         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
126  \
127         /* copy over old data into new slots */ \
128         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
129         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
130         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
131         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
132         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
133         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
134         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
135         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
136         /* free up old matrix storage */ \
137         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
138         AA = new_a; \
139         Ain->a = (MatScalar*) new_a;		   \
140         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
141         Ain->singlemalloc = PETSC_TRUE; \
142  \
143         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
144         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
145         Ain->maxnz += BS2*CHUNKSIZE; \
146         Ain->reallocs++; \
147       } \
148 
149 
150 EXTERN_C_BEGIN
151 extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
152 EXTERN_C_END
153 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
154 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
155 extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
156 
157 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
158 extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
159 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
160 extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
161 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
162 extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
163 extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
164 extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
165 extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*);
166 extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
167 
168 extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
169 extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
170 extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
171 extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
172 extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
173 
174 extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
175 extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
176 extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
177 
178 extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
179 extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
180 extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
181 extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
182 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
183 extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
184 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
185 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
186 extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
187 extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
188 extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
189 extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
190 extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
191 extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
192 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
193 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
194 extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
195 extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
196 extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
197 extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
198 extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
199 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
200 extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
201 extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
202 extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
203 extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg);
204 extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
205 extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
206 extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
207 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
208 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
209 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
210 extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
211 extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
212 extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
213 
214 extern PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
215 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
216 extern PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
217 extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
218 extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
219 extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
220 extern PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
221 extern PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
222 extern PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
223 
224 extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
225 extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
226 extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
227 extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
228 extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
229 extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
230 extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
231 extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
232 extern PetscErrorCode MatDestroy_SeqAIJ(Mat);
233 extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
234 
235 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool );
236 extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool );
237 
238 extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
239 
240 EXTERN_C_BEGIN
241 extern PetscErrorCode  MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*);
242 extern PetscErrorCode  MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*);
243 extern PetscErrorCode  MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*);
244 extern PetscErrorCode  MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
245 extern PetscErrorCode  MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
246 extern PetscErrorCode  MatCreate_SeqAIJ(Mat);
247 EXTERN_C_END
248 extern PetscErrorCode  MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
249 extern PetscErrorCode  MatDestroy_SeqAIJ(Mat);
250 
251 
252 /*
253     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \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 - negative 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: PetscSparseDensePlusDot()
269 
270 */
271 #ifdef PETSC_KERNEL_USE_UNROLL_4
272 #define PetscSparseDenseMinusDot(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 PetscSparseDenseMinusDot(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 PetscSparseDenseMinusDot(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 
298 
299 /*
300     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
301 
302   Input Parameters:
303 +  nnz - the number of entries
304 .  r - the array of vector values
305 .  xv - the matrix values for the row
306 -  xi - the column indices of the nonzeros in the row
307 
308   Output Parameter:
309 .  sum - the sum of results
310 
311   PETSc compile flags:
312 +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
313 -   PETSC_KERNEL_USE_UNROLL_2 -
314 
315 .seealso: PetscSparseDenseMinusDot()
316 
317 */
318 #ifdef PETSC_KERNEL_USE_UNROLL_4
319 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
320 if (nnz > 0) {\
321 switch (nnz & 0x3) {\
322 case 3: sum += *xv++ * r[*xi++];\
323 case 2: sum += *xv++ * r[*xi++];\
324 case 1: sum += *xv++ * r[*xi++];\
325 nnz -= 4;}\
326 while (nnz > 0) {\
327 sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\
328 	xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\
329 xv  += 4; xi += 4; nnz -= 4; }}}
330 
331 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
332 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
333 PetscInt __i,__i1,__i2;\
334 for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
335 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
336 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
337 
338 #else
339 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
340  PetscInt __i;\
341 for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];}
342 #endif
343 
344 #endif
345