xref: /petsc/src/mat/impls/aij/seq/aij.h (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
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 /*
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   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
94 } Mat_SeqAIJ_Inode;
95 
96 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
97 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
98 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
99 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
100 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
101 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
102 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
103 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
104 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
105 
106 typedef struct {
107   SEQAIJHEADER(MatScalar);
108   Mat_SeqAIJ_Inode inode;
109   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
110 
111   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
112   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
113   PetscScalar *ibdiag;                        /* inverses of block diagonals */
114   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
115   PetscScalar fshift,omega;                   /* last used omega and fshift */
116 
117   ISColoring coloring;                        /* set with MatADSetColoring() used by MatADSetValues() */
118 
119   PetscScalar       *matmult_abdense;    /* used by MatMatMult() */
120   Mat_PtAP          *ptap;               /* used by MatPtAP() */
121   Mat_MatMatMatMult *matmatmatmult;      /* used by MatMatMatMult() */
122   Mat_RARt          *rart;               /* used by MatRARt() */
123   Mat_MatMatTransMult *abt;              /* used by MatMatTransposeMult() */
124 } Mat_SeqAIJ;
125 
126 /*
127   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
128 */
129 #undef __FUNCT__
130 #define __FUNCT__ "MatSeqXAIJFreeAIJ"
131 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
132 {
133   PetscErrorCode ierr;
134   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
135   if (A->singlemalloc) {
136     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
137   } else {
138     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
139     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
140     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
141   }
142   return 0;
143 }
144 /*
145     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
146     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
147 */
148 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
149   if (NROW >= RMAX) { \
150     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
151     /* there is no extra room in row, therefore enlarge */ \
152     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
153     datatype *new_a; \
154  \
155     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); \
156     /* malloc new storage space */ \
157     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
158  \
159     /* copy over old data into new slots */ \
160     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
161     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
162     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
163     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
164     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
165     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
166     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
167     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
168     /* free up old matrix storage */ \
169     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
170     AA                = new_a; \
171     Ain->a            = (MatScalar*) new_a;                   \
172     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
173     Ain->singlemalloc = PETSC_TRUE; \
174  \
175     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
176     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
177     Ain->maxnz += BS2*CHUNKSIZE; \
178     Ain->reallocs++; \
179   } \
180 
181 
182 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
183 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
184 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
185 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
186 
187 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
188 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
189 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
190 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
191 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
192 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
193 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
194 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
195 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
196 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
197 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
198 
199 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
200 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
201 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
202 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
203 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
204 
205 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
206 PETSC_INTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
207 PETSC_INTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
208 
209 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
210 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
211 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
212 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
213 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
214 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
215 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
216 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
217 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
218 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
219 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
220 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
221 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
222 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
223 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
224 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
225 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
226 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
227 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
228 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
229 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
230 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
231 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
232 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
233 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
234 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
235 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
236 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
237 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
238 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
239 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
240 PETSC_INTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*);
241 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
242 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
243 
244 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
245 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
246 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
247 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
248 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
249 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
250 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
251 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
252 
253 PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
254 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat,Mat,PetscReal,Mat*);
255 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
256 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
257 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
258 
259 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
260 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
261 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
262 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
263 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
264 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
265 
266 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
267 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
268 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
269 
270 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
271 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
272 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
273 
274 PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
275 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
276 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
277 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
278 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
279 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
280 
281 PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
282 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
283 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
284 
285 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
286 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
287 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
288 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
289 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
290 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
291 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
292 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
293 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
294 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
295 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
296 PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
297 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
298 
299 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
300 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
301 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
302 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
303 
304 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
305 
306 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
307 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
308 PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
309 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
310 PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
311 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
312 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
313 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
314 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
315 
316 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
317 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
318 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
319 
320 /*
321     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
322 
323   Input Parameters:
324 +  nnz - the number of entries
325 .  r - the array of vector values
326 .  xv - the matrix values for the row
327 -  xi - the column indices of the nonzeros in the row
328 
329   Output Parameter:
330 .  sum - negative the sum of results
331 
332   PETSc compile flags:
333 +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
334 -   PETSC_KERNEL_USE_UNROLL_2 -
335 
336 .seealso: PetscSparseDensePlusDot()
337 
338 */
339 #if defined(PETSC_KERNEL_USE_UNROLL_4)
340 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
341     if (nnz > 0) { \
342       switch (nnz & 0x3) { \
343       case 3: sum -= *xv++ *r[*xi++]; \
344       case 2: sum -= *xv++ *r[*xi++]; \
345       case 1: sum -= *xv++ *r[*xi++]; \
346         nnz       -= 4;} \
347       while (nnz > 0) { \
348         sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \
349                xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \
350         xv += 4; xi += 4; nnz -= 4; }}}
351 
352 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
353 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
354     PetscInt __i,__i1,__i2; \
355     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
356                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
357     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
358 
359 #else
360 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
361     PetscInt __i; \
362     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
363 #endif
364 
365 
366 
367 /*
368     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
369 
370   Input Parameters:
371 +  nnz - the number of entries
372 .  r - the array of vector values
373 .  xv - the matrix values for the row
374 -  xi - the column indices of the nonzeros in the row
375 
376   Output Parameter:
377 .  sum - the sum of results
378 
379   PETSc compile flags:
380 +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
381 -   PETSC_KERNEL_USE_UNROLL_2 -
382 
383 .seealso: PetscSparseDenseMinusDot()
384 
385 */
386 #if defined(PETSC_KERNEL_USE_UNROLL_4)
387 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
388     if (nnz > 0) { \
389       switch (nnz & 0x3) { \
390       case 3: sum += *xv++ *r[*xi++]; \
391       case 2: sum += *xv++ *r[*xi++]; \
392       case 1: sum += *xv++ *r[*xi++]; \
393         nnz       -= 4;} \
394       while (nnz > 0) { \
395         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
396                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
397         xv += 4; xi += 4; nnz -= 4; }}}
398 
399 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
400 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
401     PetscInt __i,__i1,__i2; \
402     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
403                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
404     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
405 
406 #else
407 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
408     PetscInt __i; \
409     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
410 #endif
411 
412 
413 /*
414     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
415 
416   Input Parameters:
417 +  nnz - the number of entries
418 .  r - the array of vector values
419 .  xv - the matrix values for the row
420 -  xi - the column indices of the nonzeros in the row
421 
422   Output Parameter:
423 .  max - the max of results
424 
425 .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
426 
427 */
428 #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
429     PetscInt __i; \
430     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
431 
432 #endif
433