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