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