Actual source code: aij.h

petsc-3.14.0 2020-09-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: } Mat_MatMatTransMult;

 50: typedef struct { /* used by MatTransposeMatMult() */
 51:   Mat          At;           /* transpose of the first matrix */
 52:   Mat          mA;           /* maij matrix of A */
 53:   Vec          bt,ct;        /* vectors to hold locally transposed arrays of B and C */
 54:   PetscBool    updateAt;     /* flg to avoid recomputing At in MatProductNumeric_AtB_SeqAIJ_SeqAIJ() */
 55:   /* used by PtAP */
 56:   void           *data;
 57:   PetscErrorCode (*destroy)(void*);
 58: } Mat_MatTransMatMult;

 60: typedef struct {
 61:   PetscInt    *api,*apj;       /* symbolic structure of A*P */
 62:   PetscScalar *apa;            /* temporary array for storing one row of A*P */
 63: } Mat_AP;

 65: typedef struct {
 66:   MatTransposeColoring matcoloring;
 67:   Mat                  Rt;    /* sparse or dense matrix of R^T */
 68:   Mat                  RARt;  /* dense matrix of R*A*R^T */
 69:   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
 70:   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
 71:   /* free intermediate products needed for PtAP */
 72:   void                 *data;
 73:   PetscErrorCode       (*destroy)(void*);
 74: } Mat_RARt;

 76: typedef struct {
 77:   Mat BC;               /* temp matrix for storing B*C */
 78: } Mat_MatMatMatMult;

 80: /*
 81:   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
 82:   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
 83:   j[i[k]+p] is the pth column in row k.  Note that the diagonal
 84:   matrix elements are stored with the rest of the nonzeros (not separately).
 85: */

 87: /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
 88: typedef struct {
 89:   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
 90:   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
 91:   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */

 93:   PetscBool        use;
 94:   PetscInt         node_count;                     /* number of inodes */
 95:   PetscInt         *size;                          /* size of each inode */
 96:   PetscInt         limit;                          /* inode limit */
 97:   PetscInt         max_limit;                      /* maximum supported inode limit */
 98:   PetscBool        checked;                        /* if inodes have been checked for */
 99:   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
100: } Mat_SeqAIJ_Inode;

102: PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
103: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
104: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
105: PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
106: PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
107: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
108: PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
109: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
110: PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
111: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat,PetscScalar**);
112: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat,PetscScalar**);

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

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

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

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

211: PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
212: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
213: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
214: PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);

216: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
217: PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
218: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
219: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
220: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
221: PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
222: PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
223: PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
224: PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
225: PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
226: PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);

228: PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
229: PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
230: PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
231: PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
232: PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

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

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

272: #if defined(PETSC_HAVE_HYPRE)
273: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat);
274: #endif
275: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat);

277: PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat);
278: PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat);
279: PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat);

281: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
282: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,PetscReal,Mat);
283: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat);
284: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat);
285: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat);
286: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat);
287: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat);
288: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat);
289: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat);
290: #if defined(PETSC_HAVE_HYPRE)
291: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
292: #endif

294: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
295: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,Mat);

297: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat);
298: PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,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 MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
312: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
313: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void*);

315: PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat);
316: PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
317: PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
318: PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
319: PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);

321: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat);
322: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);

324: PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat,PetscInt,PetscInt,PetscRandom);
325: PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
326: PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
327: PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
328: PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar);
329: PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec);
330: PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode);
331: PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
332: PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
333: PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
334: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
335: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
336: PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
337: PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
338: PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
339: PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
340: PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);

342: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
343: PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
344: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
345: PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);

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

349: #if defined(PETSC_HAVE_MATLAB_ENGINE)
350: PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*);
351: PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*);
352: #endif
353: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
354: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
355: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);
356: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
357: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
358: #if defined(PETSC_HAVE_SCALAPACK)
359: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
360: #endif
361: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
362: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
363: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*);
364: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*);
365: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*);
366: PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
367: PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
368: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
369: PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
370: PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat);

372: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
373: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
374: PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

376: PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat);
377: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*);
378: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat);
379: PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat);
380: PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]);
381: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);

383: PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat,ISLocalToGlobalMapping*);
384: PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],MatType,Mat);

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

389:   Input Parameters:
390: +  nnz - the number of entries
391: .  r - the array of vector values
392: .  xv - the matrix values for the row
393: -  xi - the column indices of the nonzeros in the row

395:   Output Parameter:
396: .  sum - negative the sum of results

398:   PETSc compile flags:
399: +   PETSC_KERNEL_USE_UNROLL_4
400: -   PETSC_KERNEL_USE_UNROLL_2

402:   Developer Notes:
403:     The macro changes sum but not other parameters

405: .seealso: PetscSparseDensePlusDot()

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

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

433: #else
434: #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
435:     PetscInt __i; \
436:     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
437: #endif



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

444:   Input Parameters:
445: +  nnz - the number of entries
446: .  r - the array of vector values
447: .  xv - the matrix values for the row
448: -  xi - the column indices of the nonzeros in the row

450:   Output Parameter:
451: .  sum - the sum of results

453:   PETSc compile flags:
454: +   PETSC_KERNEL_USE_UNROLL_4
455: -   PETSC_KERNEL_USE_UNROLL_2

457:   Developer Notes:
458:     The macro changes sum but not other parameters

460: .seealso: PetscSparseDenseMinusDot()

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

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

488: #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)
489: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz))

491: #else
492: #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
493:     PetscInt __i; \
494:     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
495: #endif

497: #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)
498:   #include <immintrin.h>
499:   #if !defined(_MM_SCALE_8)
500:   #define _MM_SCALE_8    8
501:   #endif

503: PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n)
504: {
505:   __m512d  vec_x,vec_y,vec_vals;
506:   __m256i  vec_idx;
507:   __mmask8 mask;
508:   PetscInt j;

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

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

541:   Input Parameters:
542: +  nnz - the number of entries
543: .  r - the array of vector values
544: .  xv - the matrix values for the row
545: -  xi - the column indices of the nonzeros in the row

547:   Output Parameter:
548: .  max - the max of results

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

552: */
553: #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
554:     PetscInt __i; \
555:     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}

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

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

588: #endif