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