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