Actual source code: aij.h

petsc-3.12.5 2020-03-29
Report Typos and Errors


  5:  #include <petsc/private/matimpl.h>
  6:  #include <petscctable.h>

  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 */

 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;

 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;

 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;

 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;

 73: typedef struct {
 74:   Mat BC;               /* temp matrix for storing B*C */
 75:   PetscErrorCode (*destroy)(Mat);
 76: } Mat_MatMatMatMult;

 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: */

 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 */

 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;

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: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat,PetscScalar**);
110: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat,PetscScalar**);

112: typedef struct {
113:   SEQAIJHEADER(MatScalar);
114:   Mat_SeqAIJ_Inode inode;
115:   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */

117:   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
118:   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
119:   PetscScalar *ibdiag;                        /* inverses of block diagonals */
120:   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
121:   PetscBool   diagonaldense;                  /* all entries along the diagonal have been set; i.e. no missing diagonal terms */
122:   PetscScalar fshift,omega;                   /* last used omega and fshift */

124:   ISColoring  coloring;                       /* set with MatADSetColoring() used by MatADSetValues() */

126:   PetscScalar         *matmult_abdense;    /* used by MatMatMult() */
127:   Mat_AP              *ap;                 /* used by MatPtAP() */
128:   Mat_MatMatMatMult   *matmatmatmult;      /* used by MatMatMatMult() */
129:   Mat_RARt            *rart;               /* used by MatRARt() */
130:   Mat_MatMatTransMult *abt;                /* used by MatMatTransposeMult() */
131:   Mat_MatTransMatMult *atb;                /* used by MatTransposeMatMult() */
132: } Mat_SeqAIJ;

134: /*
135:   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
136: */
137: PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
138: {
140:   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
141:   if (A->singlemalloc) {
142:     PetscFree3(*a,*j,*i);
143:   } else {
144:     if (A->free_a)  {PetscFree(*a);}
145:     if (A->free_ij) {PetscFree(*j);}
146:     if (A->free_ij) {PetscFree(*i);}
147:   }
148:   return 0;
149: }
150: /*
151:     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
152:     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
153: */
154: #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
155:   if (NROW >= RMAX) { \
156:     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
157:     /* there is no extra room in row, therefore enlarge */ \
158:     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
159:     datatype *new_a; \
160:  \
161:     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); \
162:     /* malloc new storage space */ \
163:     PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i); \
164:  \
165:     /* copy over old data into new slots */ \
166:     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
167:     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
168:     PetscArraycpy(new_j,AJ,AI[ROW]+NROW); \
169:     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
170:     PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len); \
171:     PetscArraycpy(new_a,AA,BS2*(AI[ROW]+NROW));    \
172:     PetscArrayzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE); \
173:     PetscArraycpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len);  \
174:     /* free up old matrix storage */ \
175:     MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i); \
176:     AA                = new_a; \
177:     Ain->a            = (MatScalar*) new_a;                   \
178:     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
179:     Ain->singlemalloc = PETSC_TRUE; \
180:  \
181:     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
182:     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
183:     Ain->maxnz += BS2*CHUNKSIZE; \
184:     Ain->reallocs++; \
185:   } \

187: #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \
188:   if (NROW >= RMAX) { \
189:     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
190:     /* there is no extra room in row, therefore enlarge */ \
191:     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
192:  \
193:     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); \
194:     /* malloc new storage space */ \
195:     PetscMalloc1(new_nz,&new_j); \
196:     PetscMalloc1(AM+1,&new_i);\
197:  \
198:     /* copy over old data into new slots */ \
199:     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
200:     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
201:     PetscArraycpy(new_j,AJ,AI[ROW]+NROW); \
202:     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
203:     PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len); \
204:  \
205:     /* free up old matrix storage */ \
206:     MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i); \
207:     Ain->a            = NULL;                   \
208:     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
209:     Ain->singlemalloc = PETSC_FALSE; \
210:     Ain->free_a       = PETSC_FALSE;             \
211:  \
212:     RP          = AJ + AI[ROW];    \
213:     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
214:     Ain->maxnz += BS2*CHUNKSIZE; \
215:     Ain->reallocs++; \
216:   } \

218: PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
219: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
220: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
221: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);

223: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
224: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
225: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
226: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
227: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
228: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
229: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
230: PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
231: PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
232: PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
233: PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);

235: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
236: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
237: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
238: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
239: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

241: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);

243: PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
244: PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
245: PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
246: PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
247: PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
248: PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**);
249: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
250: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
251: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
252: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
253: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
254: PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
255: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
256: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
257: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
258: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
259: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
260: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
261: PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
262: PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
263: PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
264: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
265: PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
266: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
267: PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
268: PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
269: PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
270: PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
271: PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
272: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
273: PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
274: PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat,PetscViewer);
275: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat,PetscViewer);
276: PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
277: PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);

279: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
280: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
281: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,PetscReal,Mat*);
282: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat*);
283: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
284: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
285: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
286: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
287: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat*);
288: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
289: #if defined(PETSC_HAVE_HYPRE)
290: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
291: #endif
292: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ_Combined(Mat,Mat,PetscReal,Mat*);
293: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
294: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,Mat);
295: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
296: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
297: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Combined(Mat,Mat,Mat);

299: PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
300: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
301: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
302: PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);

304: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
305: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
306: PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
307: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
308: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
309: PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);

311: PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
312: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
313: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
314: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat);

316: PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
317: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
318: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);

320: PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
321: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
322: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
323: PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
324: PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
325: PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);

327: PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
328: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
329: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);

331: PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat,PetscInt,PetscInt,PetscRandom);
332: PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
333: PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
334: PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
335: PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
336: PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
337: PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
338: PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
339: PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
340: PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
341: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
342: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
343: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
344: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
345: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
346: PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
347: PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);

349: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
350: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
351: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
352: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);

354: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);

356: PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
357: #if defined(PETSC_HAVE_MATLAB_ENGINE)
358: PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
359: PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
360: #endif
361: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
362: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
363: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
364: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
365: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
366: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
367: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
368: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*);
369: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
370: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
371: PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
372: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
373: PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
374: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
375: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);

377: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
378: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
379: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

381: PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
382: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
383: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
384: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
385: PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]);
386: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);

388: PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat,ISLocalToGlobalMapping*);

390: /*
391:     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage

393:   Input Parameters:
394: +  nnz - the number of entries
395: .  r - the array of vector values
396: .  xv - the matrix values for the row
397: -  xi - the column indices of the nonzeros in the row

399:   Output Parameter:
400: .  sum - negative the sum of results

402:   PETSc compile flags:
403: +   PETSC_KERNEL_USE_UNROLL_4
404: -   PETSC_KERNEL_USE_UNROLL_2

406:   Developer Notes:
407:     The macro changes sum but not other parameters

409: .seealso: PetscSparseDensePlusDot()

411: */
412: #if defined(PETSC_KERNEL_USE_UNROLL_4)
413: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
414:     if (nnz > 0) { \
415:       PetscInt nnz2=nnz,rem=nnz&0x3; \
416:       switch (rem) { \
417:       case 3: sum -= *xv++ *r[*xi++]; \
418:       case 2: sum -= *xv++ *r[*xi++]; \
419:       case 1: sum -= *xv++ *r[*xi++]; \
420:         nnz2      -= rem;} \
421:       while (nnz2 > 0) { \
422:         sum -=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
423:                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
424:         xv += 4; xi += 4; nnz2 -= 4; \
425:       } \
426:       xv -= nnz; xi -= nnz; \
427:     } \
428:   }

430: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
431: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
432:     PetscInt __i,__i1,__i2; \
433:     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
434:                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
435:     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}

437: #else
438: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
439:     PetscInt __i; \
440:     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
441: #endif



445: /*
446:     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage

448:   Input Parameters:
449: +  nnz - the number of entries
450: .  r - the array of vector values
451: .  xv - the matrix values for the row
452: -  xi - the column indices of the nonzeros in the row

454:   Output Parameter:
455: .  sum - the sum of results

457:   PETSc compile flags:
458: +   PETSC_KERNEL_USE_UNROLL_4
459: -   PETSC_KERNEL_USE_UNROLL_2

461:   Developer Notes:
462:     The macro changes sum but not other parameters

464: .seealso: PetscSparseDenseMinusDot()

466: */
467: #if defined(PETSC_KERNEL_USE_UNROLL_4)
468: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
469:     if (nnz > 0) { \
470:       PetscInt nnz2=nnz,rem=nnz&0x3; \
471:       switch (rem) { \
472:       case 3: sum += *xv++ *r[*xi++]; \
473:       case 2: sum += *xv++ *r[*xi++]; \
474:       case 1: sum += *xv++ *r[*xi++]; \
475:         nnz2      -= rem;} \
476:       while (nnz2 > 0) { \
477:         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
478:                 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
479:         xv += 4; xi += 4; nnz2 -= 4; \
480:       } \
481:       xv -= nnz; xi -= nnz; \
482:     } \
483:   }

485: #elif defined(PETSC_KERNEL_USE_UNROLL_2)
486: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
487:     PetscInt __i,__i1,__i2; \
488:     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
489:                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
490:     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}

492: #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) && !defined(PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND)
493: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))

495: #else
496: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
497:     PetscInt __i; \
498:     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
499: #endif

501: #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) && !defined(PETSC_SKIP_IMMINTRIN_H_CUDAWORKAROUND)
502:   #include <immintrin.h>
503:   #if !defined(_MM_SCALE_8)
504:   #define _MM_SCALE_8    8
505:   #endif

507: PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
508: {
509:   __m512d  vec_x,vec_y,vec_vals;
510:   __m256i  vec_idx;
511:   __mmask8 mask;
512:   PetscInt j;

514:   vec_y = _mm512_setzero_pd();
515:   for (j=0; j<(n>>3); j++) {
516:     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
517:     vec_vals = _mm512_loadu_pd(aa);
518:     vec_x    = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8);
519:     vec_y    = _mm512_fmadd_pd(vec_x,vec_vals,vec_y);
520:     aj += 8; aa += 8;
521:   }
522:   /* masked load does not work on KNL, it requires avx512vl */
523:   if ((n&0x07)>2) {
524:     mask     = (__mmask8)(0xff >> (8-(n&0x07)));
525:     vec_idx  = _mm256_loadu_si256((__m256i const*)aj);
526:     vec_vals = _mm512_loadu_pd(aa);
527:     vec_x    = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8);
528:     vec_y    = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask);
529:   } else if ((n&0x07)==2) {
530:     *sum += aa[0]*x[aj[0]];
531:     *sum += aa[1]*x[aj[1]];
532:   } else if ((n&0x07)==1) {
533:     *sum += aa[0]*x[aj[0]];
534:   }
535:   if (n>2) *sum += _mm512_reduce_add_pd(vec_y);
536: /*
537:   for(j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]];
538: */
539: }
540: #endif

542: /*
543:     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage

545:   Input Parameters:
546: +  nnz - the number of entries
547: .  r - the array of vector values
548: .  xv - the matrix values for the row
549: -  xi - the column indices of the nonzeros in the row

551:   Output Parameter:
552: .  max - the max of results

554: .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()

556: */
557: #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
558:     PetscInt __i; \
559:     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}

561: /*
562:  Add column indices into table for counting the max nonzeros of merged rows
563:  */
564: #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) {       \
565:     PetscInt _j,_row,_nz,*_col;                     \
566:     if (mat) { \
567:       for (_row=0; _row<nrows; _row++) {   \
568:         _nz = mat->i[_row+1] - mat->i[_row];    \
569:         for (_j=0; _j<_nz; _j++) {               \
570:           _col = _j + mat->j + mat->i[_row];       \
571:           PetscTableAdd(ta,*_col+1,1,INSERT_VALUES);                    \
572:         }                                                               \
573:       }                                                                 \
574:     }    \
575: }

577: /*
578:  Add column indices into table for counting the nonzeros of merged rows
579:  */
580: #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) {    \
581:   PetscInt _j,_row,_nz,*_col,_i;                      \
582:     for (_i=0; _i<nrows; _i++) {\
583:       _row = rows[_i]; \
584:       _nz = mat->i[_row+1] - mat->i[_row]; \
585:       for (_j=0; _j<_nz; _j++) {                \
586:         _col = _j + mat->j + mat->i[_row];       \
587:         PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \
588:       }                                                                 \
589:     }                                                                   \
590: }

592: #endif