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