xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 2ec85ea456d4da44b4dd09fd939e615b8c7fb52b)
1 
2 #if !defined(__AIJ_H)
3 #define __AIJ_H
4 
5 #include <petsc/private/matimpl.h>
6 #include <petscctable.h>
7 
8 /*
9     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
10 */
11 #define SEQAIJHEADER(datatype)        \
12   PetscBool roworiented;              /* if true, row-oriented input, default */ \
13   PetscInt  nonew;                    /* 1 don't add new nonzeros, -1 generate error on new */ \
14   PetscInt  nounused;                 /* -1 generate error on unused space */ \
15   PetscBool singlemalloc;             /* if true a, i, and j have been obtained with one big malloc */ \
16   PetscInt  maxnz;                    /* allocated nonzeros */ \
17   PetscInt  *imax;                    /* maximum space allocated for each row */ \
18   PetscInt  *ilen;                    /* actual length of each row */ \
19   PetscBool free_imax_ilen;  \
20   PetscInt  reallocs;                 /* number of mallocs done during MatSetValues() \
21                                         as more values are set than were prealloced */\
22   PetscInt          rmax;             /* max nonzeros in any row */ \
23   PetscBool         keepnonzeropattern;   /* keeps matrix structure same in calls to MatZeroRows()*/ \
24   PetscBool         ignorezeroentries; \
25   PetscBool         free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
26   PetscBool         free_a;           /* free the numerical values when matrix is destroy */ \
27   Mat_CompressedRow compressedrow;    /* use compressed row format */                      \
28   PetscInt          nz;               /* nonzeros */                                       \
29   PetscInt          *i;               /* pointer to beginning of each row */               \
30   PetscInt          *j;               /* column values: j + i[k] - 1 is start of row k */  \
31   PetscInt          *diag;            /* pointers to diagonal elements */                  \
32   PetscInt          nonzerorowcnt;    /* how many rows have nonzero entries */             \
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;    /* sparse or dense matrix of R^T */
64   Mat                  RARt;  /* dense matrix of R*A*R^T */
65   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
66   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
67   PetscErrorCode (*destroy)(Mat);
68 } Mat_RARt;
69 
70 typedef struct {
71   Mat BC;               /* temp matrix for storing B*C */
72   PetscErrorCode (*destroy)(Mat);
73 } Mat_MatMatMatMult;
74 
75 /*
76   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
77   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
78   j[i[k]+p] is the pth column in row k.  Note that the diagonal
79   matrix elements are stored with the rest of the nonzeros (not separately).
80 */
81 
82 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
83 typedef struct {
84   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
85   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
86   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
87 
88   PetscBool        use;
89   PetscInt         node_count;                     /* number of inodes */
90   PetscInt         *size;                          /* size of each inode */
91   PetscInt         limit;                          /* inode limit */
92   PetscInt         max_limit;                      /* maximum supported inode limit */
93   PetscBool        checked;                        /* if inodes have been checked for */
94   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
95 } Mat_SeqAIJ_Inode;
96 
97 typedef struct { /* used by MatGetSubMatrices_MPIAIJ_SingleIS_Local() and MatGetSubMatrices_MPIAIJ_Local */
98   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
99   PetscInt   nrqs,nrqr;
100   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
101   PetscInt   **ptr;
102   PetscInt   *tmp;
103   PetscInt   *ctr;
104   PetscInt   *pa; /* proc array */
105   PetscInt   *req_size,*req_source1,*req_source2;
106   PetscBool  allcolumns;
107   PetscBool  singleis;
108   PetscInt   *row2proc; /* row to proc map */
109   PetscInt   nstages;
110 #if defined(PETSC_USE_CTABLE)
111   PetscTable     cmap,rmap;
112   PetscInt       *cmap_loc,*rmap_loc;
113 #else
114   PetscInt       *cmap,*rmap;
115 #endif
116 
117   PetscErrorCode (*destroy)(Mat);
118 } Mat_SubMat;
119 
120 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
121 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
122 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
123 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
124 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
125 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
126 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
127 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
128 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
129 
130 typedef struct {
131   SEQAIJHEADER(MatScalar);
132   Mat_SeqAIJ_Inode inode;
133   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
134 
135   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
136   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
137   PetscScalar *ibdiag;                        /* inverses of block diagonals */
138   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
139   PetscScalar fshift,omega;                   /* last used omega and fshift */
140 
141   ISColoring coloring;                        /* set with MatADSetColoring() used by MatADSetValues() */
142 
143   PetscScalar       *matmult_abdense;    /* used by MatMatMult() */
144   Mat_PtAP          *ptap;               /* used by MatPtAP() */
145   Mat_MatMatMatMult *matmatmatmult;      /* used by MatMatMatMult() */
146   Mat_RARt          *rart;               /* used by MatRARt() */
147   Mat_MatMatTransMult *abt;              /* used by MatMatTransposeMult() */
148 
149   /* used by MatGetSubMatrices_MPIAIJ_SingleIS_Local */
150   Mat_SubMat       *submatis1;
151 } Mat_SeqAIJ;
152 
153 /*
154   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
155 */
156 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
157 {
158   PetscErrorCode ierr;
159   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
160   if (A->singlemalloc) {
161     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
162   } else {
163     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
164     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
165     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
166   }
167   return 0;
168 }
169 /*
170     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
171     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
172 */
173 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
174   if (NROW >= RMAX) { \
175     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
176     /* there is no extra room in row, therefore enlarge */ \
177     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
178     datatype *new_a; \
179  \
180     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); \
181     /* malloc new storage space */ \
182     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
183  \
184     /* copy over old data into new slots */ \
185     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
186     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
187     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
188     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
189     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
190     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
191     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
192     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
193     /* free up old matrix storage */ \
194     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
195     AA                = new_a; \
196     Ain->a            = (MatScalar*) new_a;                   \
197     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
198     Ain->singlemalloc = PETSC_TRUE; \
199  \
200     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
201     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
202     Ain->maxnz += BS2*CHUNKSIZE; \
203     Ain->reallocs++; \
204   } \
205 
206 
207 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
208 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
209 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
210 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
211 
212 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
213 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
214 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
215 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
216 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
217 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
218 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
219 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
220 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
221 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
222 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
223 
224 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
225 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
226 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
227 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
228 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
229 
230 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
231 
232 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
233 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
234 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
235 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
236 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
237 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
238 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
239 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
240 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
241 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
242 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
243 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
244 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
245 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
246 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
247 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
248 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
249 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
250 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
251 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
252 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
253 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
254 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
255 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
256 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
257 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
258 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
259 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
260 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
261 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
262 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
263 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
264 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
265 
266 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
267 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
268 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
269 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
270 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
271 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
272 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
273 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
274 
275 PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
276 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat,Mat,PetscReal,Mat*);
277 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
278 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
279 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
280 
281 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
282 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
283 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
284 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
285 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
286 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
287 
288 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
289 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
290 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
291 
292 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
293 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
294 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
295 
296 PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
297 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
298 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
299 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
300 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
301 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
302 
303 PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
304 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
305 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
306 
307 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
308 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
309 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
310 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
311 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
312 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
313 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
314 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
315 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
316 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
317 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
318 PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
319 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
320 
321 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
322 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
323 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
324 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
325 
326 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
327 
328 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
329 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
330 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
331 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
332 PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
333 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
334 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
335 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
336 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
337 
338 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
339 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
340 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
341 
342 PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
343 
344 /*
345     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
346 
347   Input Parameters:
348 +  nnz - the number of entries
349 .  r - the array of vector values
350 .  xv - the matrix values for the row
351 -  xi - the column indices of the nonzeros in the row
352 
353   Output Parameter:
354 .  sum - negative the sum of results
355 
356   PETSc compile flags:
357 +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
358 -   PETSC_KERNEL_USE_UNROLL_2 -
359 
360 .seealso: PetscSparseDensePlusDot()
361 
362 */
363 #if defined(PETSC_KERNEL_USE_UNROLL_4)
364 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
365     if (nnz > 0) { \
366       switch (nnz & 0x3) { \
367       case 3: sum -= *xv++ *r[*xi++]; \
368       case 2: sum -= *xv++ *r[*xi++]; \
369       case 1: sum -= *xv++ *r[*xi++]; \
370         nnz       -= 4;} \
371       while (nnz > 0) { \
372         sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \
373                xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \
374         xv += 4; xi += 4; nnz -= 4; }}}
375 
376 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
377 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
378     PetscInt __i,__i1,__i2; \
379     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
380                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
381     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
382 
383 #else
384 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
385     PetscInt __i; \
386     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
387 #endif
388 
389 
390 
391 /*
392     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
393 
394   Input Parameters:
395 +  nnz - the number of entries
396 .  r - the array of vector values
397 .  xv - the matrix values for the row
398 -  xi - the column indices of the nonzeros in the row
399 
400   Output Parameter:
401 .  sum - the sum of results
402 
403   PETSc compile flags:
404 +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
405 -   PETSC_KERNEL_USE_UNROLL_2 -
406 
407 .seealso: PetscSparseDenseMinusDot()
408 
409 */
410 #if defined(PETSC_KERNEL_USE_UNROLL_4)
411 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
412     if (nnz > 0) { \
413       switch (nnz & 0x3) { \
414       case 3: sum += *xv++ *r[*xi++]; \
415       case 2: sum += *xv++ *r[*xi++]; \
416       case 1: sum += *xv++ *r[*xi++]; \
417         nnz       -= 4;} \
418       while (nnz > 0) { \
419         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
420                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
421         xv += 4; xi += 4; nnz -= 4; }}}
422 
423 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
424 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
425     PetscInt __i,__i1,__i2; \
426     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
427                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
428     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
429 
430 #else
431 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
432     PetscInt __i; \
433     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
434 #endif
435 
436 
437 /*
438     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
439 
440   Input Parameters:
441 +  nnz - the number of entries
442 .  r - the array of vector values
443 .  xv - the matrix values for the row
444 -  xi - the column indices of the nonzeros in the row
445 
446   Output Parameter:
447 .  max - the max of results
448 
449 .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
450 
451 */
452 #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
453     PetscInt __i; \
454     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
455 
456 /*
457  Add column indices into table for counting the max nonzeros of merged rows
458  */
459 #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) {       \
460     PetscInt _j,_row,_nz,*_col;                     \
461     for (_row=0; _row<nrows; _row++) {\
462       _nz = mat->i[_row+1] - mat->i[_row]; \
463       for (_j=0; _j<_nz; _j++) {                \
464         _col = _j + mat->j + mat->i[_row];       \
465         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
466       }                                                                 \
467     }                                                                   \
468 }
469 
470 /*
471  Add column indices into table for counting the nonzeros of merged rows
472  */
473 #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) {    \
474   PetscInt _j,_row,_nz,*_col,_i;                      \
475     for (_i=0; _i<nrows; _i++) {\
476       _row = rows[_i]; \
477       _nz = mat->i[_row+1] - mat->i[_row]; \
478       for (_j=0; _j<_nz; _j++) {                \
479         _col = _j + mat->j + mat->i[_row];       \
480         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
481       }                                                                 \
482     }                                                                   \
483 }
484 
485 #endif
486