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