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