xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 8fcddce65efd55a8fe3f87d4c08c15577ce4cbef)
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   PetscInt  *ipre;                    /* space preallocated for each row by user */ \
20   PetscBool free_imax_ilen;  \
21   PetscInt  reallocs;                 /* number of mallocs done during MatSetValues() \
22                                         as more values are set than were prealloced */\
23   PetscInt          rmax;             /* max nonzeros in any row */ \
24   PetscBool         keepnonzeropattern;   /* keeps matrix structure same in calls to MatZeroRows()*/ \
25   PetscBool         ignorezeroentries; \
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   Mat_SubSppt       *submatis1         /* used by MatCreateSubMatrices_MPIXAIJ_Local */
42 
43 typedef struct {
44   MatTransposeColoring matcoloring;
45   Mat                  Bt_den;       /* dense matrix of B^T */
46   Mat                  ABt_den;      /* dense matrix of A*B^T */
47   PetscBool            usecoloring;
48   PetscErrorCode (*destroy)(Mat);
49 } Mat_MatMatTransMult;
50 
51 typedef struct { /* used by MatTransposeMatMult() */
52   Mat          At;           /* transpose of the first matrix */
53   Mat          mA;           /* maij matrix of A */
54   Vec          bt,ct;        /* vectors to hold locally transposed arrays of B and C */
55   PetscErrorCode (*destroy)(Mat);
56 } Mat_MatTransMatMult;
57 
58 typedef struct {
59   PetscInt    *api,*apj;       /* symbolic structure of A*P */
60   PetscScalar *apa;            /* temporary array for storing one row of A*P */
61   PetscErrorCode (*destroy)(Mat);
62 } Mat_AP;
63 
64 typedef struct {
65   MatTransposeColoring matcoloring;
66   Mat                  Rt;    /* sparse or dense matrix of R^T */
67   Mat                  RARt;  /* dense matrix of R*A*R^T */
68   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
69   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
70   PetscErrorCode (*destroy)(Mat);
71 } Mat_RARt;
72 
73 typedef struct {
74   Mat BC;               /* temp matrix for storing B*C */
75   PetscErrorCode (*destroy)(Mat);
76 } Mat_MatMatMatMult;
77 
78 /*
79   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
80   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
81   j[i[k]+p] is the pth column in row k.  Note that the diagonal
82   matrix elements are stored with the rest of the nonzeros (not separately).
83 */
84 
85 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
86 typedef struct {
87   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
88   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
89   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
90 
91   PetscBool        use;
92   PetscInt         node_count;                     /* number of inodes */
93   PetscInt         *size;                          /* size of each inode */
94   PetscInt         limit;                          /* inode limit */
95   PetscInt         max_limit;                      /* maximum supported inode limit */
96   PetscBool        checked;                        /* if inodes have been checked for */
97   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
98 } Mat_SeqAIJ_Inode;
99 
100 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
101 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
102 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
103 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
104 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
105 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
106 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
107 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
108 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
109 
110 typedef struct {
111   SEQAIJHEADER(MatScalar);
112   Mat_SeqAIJ_Inode inode;
113   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
114 
115   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
116   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
117   PetscScalar *ibdiag;                        /* inverses of block diagonals */
118   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
119   PetscBool   diagonaldense;                  /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
120   PetscScalar fshift,omega;                   /* last used omega and fshift */
121 
122   ISColoring  coloring;                       /* set with MatADSetColoring() used by MatADSetValues() */
123 
124   PetscScalar         *matmult_abdense;    /* used by MatMatMult() */
125   Mat_AP              *ap;                 /* used by MatPtAP() */
126   Mat_MatMatMatMult   *matmatmatmult;      /* used by MatMatMatMult() */
127   Mat_RARt            *rart;               /* used by MatRARt() */
128   Mat_MatMatTransMult *abt;                /* used by MatMatTransposeMult() */
129   Mat_MatTransMatMult *atb;                /* used by MatTransposeMatMult() */
130 } Mat_SeqAIJ;
131 
132 /*
133   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
134 */
135 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
136 {
137   PetscErrorCode ierr;
138   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
139   if (A->singlemalloc) {
140     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
141   } else {
142     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
143     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
144     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
145   }
146   return 0;
147 }
148 /*
149     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
150     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
151 */
152 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
153   if (NROW >= RMAX) { \
154     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
155     /* there is no extra room in row, therefore enlarge */ \
156     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
157     datatype *new_a; \
158  \
159     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); \
160     /* malloc new storage space */ \
161     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
162  \
163     /* copy over old data into new slots */ \
164     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
165     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
166     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
167     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
168     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
169     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
170     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
171     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
172     /* free up old matrix storage */ \
173     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
174     AA                = new_a; \
175     Ain->a            = (MatScalar*) new_a;                   \
176     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
177     Ain->singlemalloc = PETSC_TRUE; \
178  \
179     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
180     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
181     Ain->maxnz += BS2*CHUNKSIZE; \
182     Ain->reallocs++; \
183   } \
184 
185 #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \
186   if (NROW >= RMAX) { \
187     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
188     /* there is no extra room in row, therefore enlarge */ \
189     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
190  \
191     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); \
192     /* malloc new storage space */ \
193     ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \
194     ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\
195  \
196     /* copy over old data into new slots */ \
197     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
198     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
199     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
200     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
201     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
202  \
203     /* free up old matrix storage */ \
204     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
205     Ain->a            = NULL;                   \
206     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
207     Ain->singlemalloc = PETSC_FALSE; \
208     Ain->free_a       = PETSC_FALSE;             \
209  \
210     RP          = AJ + AI[ROW];    \
211     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
212     Ain->maxnz += BS2*CHUNKSIZE; \
213     Ain->reallocs++; \
214   } \
215 
216 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
217 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
218 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
219 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
220 
221 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
222 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
223 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
224 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
225 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
226 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
227 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
228 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
229 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
230 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
231 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
232 
233 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
234 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
235 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
236 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
237 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
238 
239 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
240 
241 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
242 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
243 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ_FAST(Mat,MatReuse,Mat *);
244 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
245 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
246 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
247 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
248 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
249 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
250 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
251 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
252 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
253 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
254 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
255 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
256 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
257 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
258 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
259 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
260 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
261 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
262 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
263 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
264 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
265 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
266 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
267 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
268 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
269 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
270 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
271 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
272 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
273 PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat,PetscViewer);
274 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat,PetscViewer);
275 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
276 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
277 
278 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
279 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
280 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat*);
281 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
282 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
283 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
284 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
285 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat*);
286 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ_Combined(Mat,Mat,PetscReal,Mat*);
287 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
288 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
289 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
290 
291 PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
292 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
293 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
294 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
295 
296 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
297 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
298 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
299 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
300 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
301 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
302 
303 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
304 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
305 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
306 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat);
307 
308 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
309 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
310 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
311 
312 PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
313 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
314 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
315 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
316 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
317 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
318 
319 PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
320 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
321 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
322 
323 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
324 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
325 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
326 PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
327 PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
328 PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
329 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
330 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
331 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
332 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
333 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
334 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
335 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
336 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
337 PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
338 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
339 
340 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
341 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
342 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
343 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
344 
345 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
346 
347 PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
348 #if defined(PETSC_HAVE_MATLAB_ENGINE)
349 PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
350 PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
351 #endif
352 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
353 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
354 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
355 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
356 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
357 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
358 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
359 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*);
360 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
361 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
362 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
363 PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
364 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
365 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
366 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
367 
368 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
369 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
370 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
371 
372 PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
373 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
374 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
375 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
376 PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]);
377 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
378 
379 /*
380     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
381 
382   Input Parameters:
383 +  nnz - the number of entries
384 .  r - the array of vector values
385 .  xv - the matrix values for the row
386 -  xi - the column indices of the nonzeros in the row
387 
388   Output Parameter:
389 .  sum - negative the sum of results
390 
391   PETSc compile flags:
392 +   PETSC_KERNEL_USE_UNROLL_4
393 -   PETSC_KERNEL_USE_UNROLL_2
394 
395   Developer Notes:
396     The macro changes sum but not other parameters
397 
398 .seealso: PetscSparseDensePlusDot()
399 
400 */
401 #if defined(PETSC_KERNEL_USE_UNROLL_4)
402 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
403     if (nnz > 0) { \
404       PetscInt nnz2=nnz,rem=nnz&0x3; \
405       switch (rem) { \
406       case 3: sum -= *xv++ *r[*xi++]; \
407       case 2: sum -= *xv++ *r[*xi++]; \
408       case 1: sum -= *xv++ *r[*xi++]; \
409         nnz2      -= rem;} \
410       while (nnz2 > 0) { \
411         sum -=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
412                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
413         xv += 4; xi += 4; nnz2 -= 4; \
414       } \
415       xv -= nnz; xi -= nnz; \
416     } \
417   }
418 
419 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
420 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
421     PetscInt __i,__i1,__i2; \
422     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
423                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
424     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
425 
426 #else
427 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
428     PetscInt __i; \
429     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
430 #endif
431 
432 
433 
434 /*
435     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
436 
437   Input Parameters:
438 +  nnz - the number of entries
439 .  r - the array of vector values
440 .  xv - the matrix values for the row
441 -  xi - the column indices of the nonzeros in the row
442 
443   Output Parameter:
444 .  sum - the sum of results
445 
446   PETSc compile flags:
447 +   PETSC_KERNEL_USE_UNROLL_4
448 -   PETSC_KERNEL_USE_UNROLL_2
449 
450   Developer Notes:
451     The macro changes sum but not other parameters
452 
453 .seealso: PetscSparseDenseMinusDot()
454 
455 */
456 #if defined(PETSC_KERNEL_USE_UNROLL_4)
457 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
458     if (nnz > 0) { \
459       PetscInt nnz2=nnz,rem=nnz&0x3; \
460       switch (rem) { \
461       case 3: sum += *xv++ *r[*xi++]; \
462       case 2: sum += *xv++ *r[*xi++]; \
463       case 1: sum += *xv++ *r[*xi++]; \
464         nnz2      -= rem;} \
465       while (nnz2 > 0) { \
466         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
467                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
468         xv += 4; xi += 4; nnz2 -= 4; \
469       } \
470       xv -= nnz; xi -= nnz; \
471     } \
472   }
473 
474 #elif defined(PETSC_KERNEL_USE_UNROLL_2)
475 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
476     PetscInt __i,__i1,__i2; \
477     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
478                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
479     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
480 
481 #elif defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
482 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))
483 
484 #else
485 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
486     PetscInt __i; \
487     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
488 #endif
489 
490 #if defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
491   #include <immintrin.h>
492   #if !defined(_MM_SCALE_8)
493   #define _MM_SCALE_8    8
494   #endif
495 
496 PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
497 {
498   __m512d  vec_x,vec_y,vec_vals;
499   __m256i  vec_idx;
500   __mmask8 mask;
501   PetscInt j;
502 
503   vec_y = _mm512_setzero_pd();
504   for (j=0; j<(n>>3); j++) {
505     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
506     vec_vals = _mm512_loadu_pd(aa);
507     vec_x    = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8);
508     vec_y    = _mm512_fmadd_pd(vec_x,vec_vals,vec_y);
509     aj += 8; aa += 8;
510   }
511   /* masked load does not work on KNL, it requires avx512vl */
512   if ((n&0x07)>2) {
513     mask     = (__mmask8)(0xff >> (8-(n&0x07)));
514     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
515     vec_vals = _mm512_loadu_pd(aa);
516     vec_x    = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8);
517     vec_y    = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask);
518   } else if ((n&0x07)==2) {
519     *sum += aa[0]*x[aj[0]];
520     *sum += aa[1]*x[aj[1]];
521   } else if ((n&0x07)==1) {
522     *sum += aa[0]*x[aj[0]];
523   }
524   if (n>2) *sum += _mm512_reduce_add_pd(vec_y);
525 /*
526   for(j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
527 */
528 }
529 #endif
530 
531 /*
532     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
533 
534   Input Parameters:
535 +  nnz - the number of entries
536 .  r - the array of vector values
537 .  xv - the matrix values for the row
538 -  xi - the column indices of the nonzeros in the row
539 
540   Output Parameter:
541 .  max - the max of results
542 
543 .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
544 
545 */
546 #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
547     PetscInt __i; \
548     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
549 
550 /*
551  Add column indices into table for counting the max nonzeros of merged rows
552  */
553 #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) {       \
554     PetscInt _j,_row,_nz,*_col;                     \
555     if (mat) { \
556       for (_row=0; _row<nrows; _row++) {   \
557         _nz = mat->i[_row+1] - mat->i[_row];    \
558         for (_j=0; _j<_nz; _j++) {               \
559           _col = _j + mat->j + mat->i[_row];       \
560           PetscTableAdd(ta,*_col+1,1,INSERT_VALUES);                    \
561         }                                                               \
562       }                                                                 \
563     }    \
564 }
565 
566 /*
567  Add column indices into table for counting the nonzeros of merged rows
568  */
569 #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) {    \
570   PetscInt _j,_row,_nz,*_col,_i;                      \
571     for (_i=0; _i<nrows; _i++) {\
572       _row = rows[_i]; \
573       _nz = mat->i[_row+1] - mat->i[_row]; \
574       for (_j=0; _j<_nz; _j++) {                \
575         _col = _j + mat->j + mat->i[_row];       \
576         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
577       }                                                                 \
578     }                                                                   \
579 }
580 
581 #endif
582