1 #pragma once 2 3 #include <petsc/private/matimpl.h> 4 #include <petsc/private/hashmapi.h> 5 #include <petsc/private/hashmapijv.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 PetscMPIInt nrqs, nrqr; 13 PetscInt **rbuf1, **rbuf2, **rbuf3, **sbuf1, **sbuf2; 14 PetscInt **ptr; 15 PetscInt *tmp; 16 PetscInt *ctr; 17 PetscMPIInt *pa; /* process array */ 18 PetscInt *req_size; 19 PetscMPIInt *req_source1, *req_source2; 20 PetscBool allcolumns, allrows; 21 PetscBool singleis; 22 PetscMPIInt *row2proc; /* row to process (MPI rank) 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 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 nonzero 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 PetscObjectState diagNonzeroState; /* nonzero state of the matrix when diag was obtained */ \ 69 PetscBool diagDense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */ \ 70 PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 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 } MatProductCtx_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 PetscCtxDestroyFn *destroy; 94 } MatProductCtx_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 } MatProductCtx_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 PetscCtxDestroyFn *destroy; 110 } MatProductCtx_RARt; 111 112 typedef struct { 113 Mat BC; /* temp matrix for storing B*C */ 114 } MatProductCtx_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 /* data for MatSOR_SeqAIJ_Inode() */ 126 MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrices */ 127 PetscInt bdiagsize; /* length of bdiag and ibdiag */ 128 PetscObjectState ibdiagState; /* state of the matrix when ibdiag[] and bdiag[] were constructed */ 129 130 PetscBool use; 131 PetscInt node_count; /* number of inodes */ 132 PetscInt *size_csr; /* inode sizes in csr with size_csr[0] = 0 and i-th node size = size_csr[i+1] - size_csr[i], to facilitate parallel computation */ 133 PetscInt limit; /* inode limit */ 134 PetscInt max_limit; /* maximum supported inode limit */ 135 PetscBool checked; /* if inodes have been checked for */ 136 PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 137 } Mat_SeqAIJ_Inode; 138 139 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer); 140 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType); 141 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 142 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 143 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool); 144 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *); 145 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool); 146 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *); 147 PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **); 148 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **); 149 150 typedef struct { 151 SEQAIJHEADER(MatScalar); 152 Mat_SeqAIJ_Inode inode; 153 MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 154 155 /* data needed for MatSOR_SeqAIJ() */ 156 PetscScalar *mdiag, *idiag; /* diagonal values, inverse of diagonal entries */ 157 PetscScalar *ssor_work; /* workspace for Eisenstat trick */ 158 PetscObjectState idiagState; /* state of the matrix when mdiag and idiag was obtained */ 159 PetscScalar fshift, omega; /* last used omega and fshift */ 160 161 PetscScalar *ibdiag; /* inverses of block diagonals */ 162 PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 163 164 /* MatSetValues() via hash related fields */ 165 PetscHMapIJV ht; 166 PetscInt *dnz; 167 struct _MatOps cops; 168 } Mat_SeqAIJ; 169 170 typedef struct { 171 PetscInt nz; /* nz of the matrix after assembly */ 172 PetscCount n; /* Number of entries in MatSetPreallocationCOO() */ 173 PetscCount Atot; /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */ 174 PetscCount *jmap; /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */ 175 PetscCount *perm; /* The permutation array in sorting (i,j) by row and then by col */ 176 } MatCOOStruct_SeqAIJ; 177 178 #define MatSeqXAIJGetOptions_Private(A) \ 179 { \ 180 const PetscBool oldvalues = (PetscBool)(A != PETSC_NULLPTR); \ 181 PetscInt nonew = 0, nounused = 0; \ 182 PetscBool roworiented = PETSC_FALSE; \ 183 if (oldvalues) { \ 184 nonew = ((Mat_SeqAIJ *)A->data)->nonew; \ 185 nounused = ((Mat_SeqAIJ *)A->data)->nounused; \ 186 roworiented = ((Mat_SeqAIJ *)A->data)->roworiented; \ 187 } \ 188 (void)0 189 190 #define MatSeqXAIJRestoreOptions_Private(A) \ 191 if (oldvalues) { \ 192 ((Mat_SeqAIJ *)A->data)->nonew = nonew; \ 193 ((Mat_SeqAIJ *)A->data)->nounused = nounused; \ 194 ((Mat_SeqAIJ *)A->data)->roworiented = roworiented; \ 195 } \ 196 } \ 197 (void)0 198 199 static inline PetscErrorCode MatXAIJAllocatea(Mat A, PetscInt nz, PetscScalar **array) 200 { 201 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 202 203 PetscFunctionBegin; 204 PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscScalar), (void **)array)); 205 a->free_a = PETSC_TRUE; 206 PetscFunctionReturn(PETSC_SUCCESS); 207 } 208 209 static inline PetscErrorCode MatXAIJDeallocatea(Mat A, PetscScalar **array) 210 { 211 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 212 213 PetscFunctionBegin; 214 if (a->free_a) PetscCall(PetscShmgetDeallocateArray((void **)array)); 215 a->free_a = PETSC_FALSE; 216 PetscFunctionReturn(PETSC_SUCCESS); 217 } 218 219 /* 220 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 221 */ 222 static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i) 223 { 224 Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data; 225 226 PetscFunctionBegin; 227 if (A->free_a) PetscCall(PetscShmgetDeallocateArray((void **)a)); 228 if (A->free_ij) PetscCall(PetscShmgetDeallocateArray((void **)j)); 229 if (A->free_ij) PetscCall(PetscShmgetDeallocateArray((void **)i)); 230 PetscFunctionReturn(PETSC_SUCCESS); 231 } 232 /* 233 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 234 This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 235 */ 236 #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \ 237 do { \ 238 if (NROW >= RMAX) { \ 239 Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 240 PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 241 datatype *new_a; \ 242 \ 243 PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc. Use MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \ 244 /* malloc new storage space */ \ 245 PetscCall(PetscShmgetAllocateArray(BS2 * new_nz, sizeof(PetscScalar), (void **)&new_a)); \ 246 PetscCall(PetscShmgetAllocateArray(new_nz, sizeof(PetscInt), (void **)&new_j)); \ 247 PetscCall(PetscShmgetAllocateArray(AM + 1, sizeof(PetscInt), (void **)&new_i)); \ 248 Ain->free_a = PETSC_TRUE; \ 249 Ain->free_ij = PETSC_TRUE; \ 250 /* copy over old data into new slots */ \ 251 for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 252 for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 253 PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 254 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 255 PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, PetscSafePointerPlusOffset(AJ, AI[ROW] + NROW), len)); \ 256 PetscCall(PetscArraycpy(new_a, AA, BS2 * (AI[ROW] + NROW))); \ 257 PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \ 258 PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), PetscSafePointerPlusOffset(AA, BS2 * (AI[ROW] + NROW)), BS2 * len)); \ 259 /* free up old matrix storage */ \ 260 PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 261 AA = new_a; \ 262 Ain->a = new_a; \ 263 AI = Ain->i = new_i; \ 264 AJ = Ain->j = new_j; \ 265 \ 266 RP = AJ + AI[ROW]; \ 267 AP = AA + BS2 * AI[ROW]; \ 268 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 269 Ain->maxnz += BS2 * CHUNKSIZE; \ 270 Ain->reallocs++; \ 271 Amat->nonzerostate++; \ 272 } \ 273 } while (0) 274 275 #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \ 276 do { \ 277 if (NROW >= RMAX) { \ 278 Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 279 /* there is no extra room in row, therefore enlarge */ \ 280 PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 281 \ 282 PetscCheck(NONEW != -2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "New nonzero at (%" PetscInt_FMT ",%" PetscInt_FMT ") caused a malloc. Use MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check", ROW, COL); \ 283 /* malloc new storage space */ \ 284 PetscCall(PetscShmgetAllocateArray(new_nz, sizeof(PetscInt), (void **)&new_j)); \ 285 PetscCall(PetscShmgetAllocateArray(AM + 1, sizeof(PetscInt), (void **)&new_i)); \ 286 Ain->free_a = PETSC_FALSE; \ 287 Ain->free_ij = PETSC_TRUE; \ 288 \ 289 /* copy over old data into new slots */ \ 290 for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 291 for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 292 PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 293 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 294 PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \ 295 \ 296 /* free up old matrix storage */ \ 297 PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 298 Ain->a = NULL; \ 299 AI = Ain->i = new_i; \ 300 AJ = Ain->j = new_j; \ 301 \ 302 RP = AJ + AI[ROW]; \ 303 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 304 Ain->maxnz += BS2 * CHUNKSIZE; \ 305 Ain->reallocs++; \ 306 Amat->nonzerostate++; \ 307 } \ 308 } while (0) 309 310 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *); 311 PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]); 312 313 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 314 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *); 315 316 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 317 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 318 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 319 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 320 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *); 321 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure); 322 PETSC_EXTERN PetscErrorCode MatGetDiagonalMarkers_SeqAIJ(Mat, const PetscInt **, PetscBool *); 323 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **); 324 325 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec); 326 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec); 327 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec); 328 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec); 329 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec); 330 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 331 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 332 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 333 334 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool); 335 336 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 337 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 338 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]); 339 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *); 340 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *); 341 342 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **); 343 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 344 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 345 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 346 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *); 347 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *); 348 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec); 349 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec); 350 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec); 351 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec); 352 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec); 353 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec); 354 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec); 355 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec); 356 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 357 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat); 358 PETSC_INTERN PetscErrorCode MatMatSolveTranspose_SeqAIJ(Mat, Mat, Mat); 359 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *); 360 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring); 361 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring); 362 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt); 363 PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer); 364 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer); 365 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer); 366 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 367 368 #if defined(PETSC_HAVE_HYPRE) 369 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat); 370 #endif 371 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat); 372 373 PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat); 374 PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat); 375 PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat); 376 377 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 378 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat); 379 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat); 380 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat); 381 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat); 382 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat); 383 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat); 384 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat); 385 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat); 386 #if defined(PETSC_HAVE_HYPRE) 387 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat); 388 #endif 389 390 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 391 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat); 392 393 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat); 394 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat); 395 396 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat); 397 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 398 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat); 399 400 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 401 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat); 402 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat); 403 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 404 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat); 405 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat); 406 407 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 408 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 409 PETSC_INTERN PetscErrorCode MatProductCtxDestroy_SeqAIJ_MatTransMatMult(void **); 410 411 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 412 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 413 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring); 414 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat); 415 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat); 416 417 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat); 418 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat); 419 420 PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom); 421 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 422 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 423 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 424 PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar); 425 PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec); 426 PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode); 427 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure); 428 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 429 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 430 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 431 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 432 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 433 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 434 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 435 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer); 436 437 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 438 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 439 440 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *); 441 442 #if defined(PETSC_HAVE_MATLAB) 443 PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *); 444 PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *); 445 #endif 446 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *); 447 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *); 448 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *); 449 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *); 450 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *); 451 #if defined(PETSC_HAVE_SCALAPACK) 452 PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *); 453 #endif 454 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *); 455 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *); 456 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *); 457 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *); 458 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *); 459 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS); 460 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *); 461 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 462 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType); 463 PETSC_INTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat); 464 465 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *); 466 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 467 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 468 469 PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat); 470 PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqAIJ(Mat, PetscBool); 471 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *); 472 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat); 473 PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat); 474 PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]); 475 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *); 476 477 PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat); 478 479 PETSC_INTERN PetscErrorCode MatResetPreallocation_SeqAIJ_Private(Mat A, PetscBool *memoryreset); 480 481 PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *); 482 483 /* 484 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 485 486 Input Parameters: 487 + nnz - the number of entries 488 . r - the array of vector values 489 . xv - the matrix values for the row 490 - xi - the column indices of the nonzeros in the row 491 492 Output Parameter: 493 . sum - negative the sum of results 494 495 PETSc compile flags: 496 + PETSC_KERNEL_USE_UNROLL_4 497 - PETSC_KERNEL_USE_UNROLL_2 498 499 Developer Note: 500 The macro changes sum but not other parameters 501 502 .seealso: `PetscSparseDensePlusDot()` 503 */ 504 #if defined(PETSC_KERNEL_USE_UNROLL_4) 505 #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 506 do { \ 507 if (nnz > 0) { \ 508 PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 509 switch (rem) { \ 510 case 3: \ 511 sum -= *xv++ * r[*xi++]; \ 512 case 2: \ 513 sum -= *xv++ * r[*xi++]; \ 514 case 1: \ 515 sum -= *xv++ * r[*xi++]; \ 516 nnz2 -= rem; \ 517 } \ 518 while (nnz2 > 0) { \ 519 sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 520 xv += 4; \ 521 xi += 4; \ 522 nnz2 -= 4; \ 523 } \ 524 xv -= nnz; \ 525 xi -= nnz; \ 526 } \ 527 } while (0) 528 529 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 530 #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 531 do { \ 532 PetscInt __i, __i1, __i2; \ 533 for (__i = 0; __i < nnz - 1; __i += 2) { \ 534 __i1 = xi[__i]; \ 535 __i2 = xi[__i + 1]; \ 536 sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 537 } \ 538 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \ 539 } while (0) 540 541 #else 542 #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 543 do { \ 544 PetscInt __i; \ 545 for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \ 546 } while (0) 547 #endif 548 549 /* 550 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 551 552 Input Parameters: 553 + nnz - the number of entries 554 . r - the array of vector values 555 . xv - the matrix values for the row 556 - xi - the column indices of the nonzeros in the row 557 558 Output Parameter: 559 . sum - the sum of results 560 561 PETSc compile flags: 562 + PETSC_KERNEL_USE_UNROLL_4 563 - PETSC_KERNEL_USE_UNROLL_2 564 565 Developer Note: 566 The macro changes sum but not other parameters 567 568 .seealso: `PetscSparseDenseMinusDot()` 569 */ 570 #if defined(PETSC_KERNEL_USE_UNROLL_4) 571 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 572 do { \ 573 if (nnz > 0) { \ 574 PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 575 switch (rem) { \ 576 case 3: \ 577 sum += *xv++ * r[*xi++]; \ 578 case 2: \ 579 sum += *xv++ * r[*xi++]; \ 580 case 1: \ 581 sum += *xv++ * r[*xi++]; \ 582 nnz2 -= rem; \ 583 } \ 584 while (nnz2 > 0) { \ 585 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 586 xv += 4; \ 587 xi += 4; \ 588 nnz2 -= 4; \ 589 } \ 590 xv -= nnz; \ 591 xi -= nnz; \ 592 } \ 593 } while (0) 594 595 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 596 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 597 do { \ 598 PetscInt __i, __i1, __i2; \ 599 for (__i = 0; __i < nnz - 1; __i += 2) { \ 600 __i1 = xi[__i]; \ 601 __i2 = xi[__i + 1]; \ 602 sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 603 } \ 604 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \ 605 } while (0) 606 607 #elif !(defined(__GNUC__) && defined(_OPENMP)) && 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) 608 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz)) 609 610 #else 611 #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 612 do { \ 613 PetscInt __i; \ 614 for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \ 615 } while (0) 616 #endif 617 618 #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) 619 #include <immintrin.h> 620 #if !defined(_MM_SCALE_8) 621 #define _MM_SCALE_8 8 622 #endif 623 624 static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n) 625 { 626 __m512d vec_x, vec_y, vec_vals; 627 __m256i vec_idx; 628 PetscInt j; 629 630 vec_y = _mm512_setzero_pd(); 631 for (j = 0; j < (n >> 3); j++) { 632 vec_idx = _mm256_loadu_si256((__m256i const *)aj); 633 vec_vals = _mm512_loadu_pd(aa); 634 vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8); 635 vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y); 636 aj += 8; 637 aa += 8; 638 } 639 #if defined(__AVX512VL__) 640 /* masked load requires avx512vl, which is not supported by KNL */ 641 if (n & 0x07) { 642 __mmask8 mask; 643 mask = (__mmask8)(0xff >> (8 - (n & 0x07))); 644 vec_idx = _mm256_mask_loadu_epi32(vec_idx, mask, aj); 645 vec_vals = _mm512_mask_loadu_pd(vec_vals, mask, aa); 646 vec_x = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8); 647 vec_y = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask); 648 } 649 *sum += _mm512_reduce_add_pd(vec_y); 650 #else 651 *sum += _mm512_reduce_add_pd(vec_y); 652 for (j = 0; j < (n & 0x07); j++) *sum += aa[j] * x[aj[j]]; 653 #endif 654 } 655 #endif 656 657 /* 658 PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 659 660 Input Parameters: 661 + nnz - the number of entries 662 . r - the array of vector values 663 . xv - the matrix values for the row 664 - xi - the column indices of the nonzeros in the row 665 666 Output Parameter: 667 . max - the max of results 668 669 .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()` 670 */ 671 #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \ 672 do { \ 673 for (PetscInt __i = 0; __i < (nnz); __i++) max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); \ 674 } while (0) 675 676 /* 677 Add column indices into table for counting the max nonzeros of merged rows 678 */ 679 #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \ 680 do { \ 681 if (mat) { \ 682 for (PetscInt _row = 0; _row < (nrows); _row++) { \ 683 const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \ 684 for (PetscInt _j = 0; _j < _nz; _j++) { \ 685 PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \ 686 PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \ 687 } \ 688 } \ 689 } \ 690 } while (0) 691 692 /* 693 Add column indices into table for counting the nonzeros of merged rows 694 */ 695 #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \ 696 do { \ 697 for (PetscInt _i = 0; _i < (nrows); _i++) { \ 698 const PetscInt _row = (rows)[_i]; \ 699 const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \ 700 for (PetscInt _j = 0; _j < _nz; _j++) { \ 701 PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \ 702 PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \ 703 } \ 704 } \ 705 } while (0) 706