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