1a4963045SJacob Faibussowitsch #pragma once 2e6b907acSBarry Smith 3af0996ceSBarry Smith #include <petsc/private/matimpl.h> 4eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h> 526cec326SBarry Smith #include <petsc/private/hashmapijv.h> 68f690400SShri Abhyankar 75b17b1ebSStefano Zampini /* 85b17b1ebSStefano Zampini Used by MatCreateSubMatrices_MPIXAIJ_Local() 95b17b1ebSStefano Zampini */ 105b17b1ebSStefano Zampini typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */ 115b17b1ebSStefano Zampini PetscInt id; /* index of submats, only submats[0] is responsible for deleting some arrays below */ 126497c311SBarry Smith PetscMPIInt nrqs, nrqr; 135b17b1ebSStefano Zampini PetscInt **rbuf1, **rbuf2, **rbuf3, **sbuf1, **sbuf2; 145b17b1ebSStefano Zampini PetscInt **ptr; 155b17b1ebSStefano Zampini PetscInt *tmp; 165b17b1ebSStefano Zampini PetscInt *ctr; 176497c311SBarry Smith PetscMPIInt *pa; /* process array */ 186497c311SBarry Smith PetscInt *req_size; 196497c311SBarry Smith PetscMPIInt *req_source1, *req_source2; 205b17b1ebSStefano Zampini PetscBool allcolumns, allrows; 215b17b1ebSStefano Zampini PetscBool singleis; 226497c311SBarry Smith PetscMPIInt *row2proc; /* row to process (MPI rank) map */ 235b17b1ebSStefano Zampini PetscInt nstages; 245b17b1ebSStefano Zampini #if defined(PETSC_USE_CTABLE) 255b17b1ebSStefano Zampini PetscHMapI cmap, rmap; 265b17b1ebSStefano Zampini PetscInt *cmap_loc, *rmap_loc; 275b17b1ebSStefano Zampini #else 285b17b1ebSStefano Zampini PetscInt *cmap, *rmap; 295b17b1ebSStefano Zampini #endif 305b17b1ebSStefano Zampini PetscErrorCode (*destroy)(Mat); 315b17b1ebSStefano Zampini } Mat_SubSppt; 325b17b1ebSStefano Zampini 33d67d9f35SJunchao Zhang /* Operations provided by MATSEQAIJ and its subclasses */ 34d67d9f35SJunchao Zhang typedef struct { 35d67d9f35SJunchao Zhang PetscErrorCode (*getarray)(Mat, PetscScalar **); 36d67d9f35SJunchao Zhang PetscErrorCode (*restorearray)(Mat, PetscScalar **); 37d67d9f35SJunchao Zhang PetscErrorCode (*getarrayread)(Mat, const PetscScalar **); 38d67d9f35SJunchao Zhang PetscErrorCode (*restorearrayread)(Mat, const PetscScalar **); 39d67d9f35SJunchao Zhang PetscErrorCode (*getarraywrite)(Mat, PetscScalar **); 40d67d9f35SJunchao Zhang PetscErrorCode (*restorearraywrite)(Mat, PetscScalar **); 417ee59b9bSJunchao Zhang PetscErrorCode (*getcsrandmemtype)(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *); 42d67d9f35SJunchao Zhang } Mat_SeqAIJOps; 43d67d9f35SJunchao Zhang 44e6b907acSBarry Smith /* 45e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 46e6b907acSBarry Smith */ 47421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 48ace3abfcSBarry Smith PetscBool roworiented; /* if true, row-oriented input, default */ \ 49e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 5028b2fa4aSMatthew Knepley PetscInt nounused; /* -1 generate error on unused space */ \ 51e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */ \ 52e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */ \ 53e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */ \ 54846b4da1SFande Kong PetscInt *ipre; /* space preallocated for each row by user */ \ 55ace3abfcSBarry Smith PetscBool free_imax_ilen; \ 56e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 57e6b907acSBarry Smith as more values are set than were prealloced */ \ 58e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */ \ 59*0b4b7b1cSBarry Smith PetscBool keepnonzeropattern; /* keeps matrix nonzero structure same in calls to MatZeroRows()*/ \ 60ace3abfcSBarry Smith PetscBool ignorezeroentries; \ 61ace3abfcSBarry Smith PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 62ace3abfcSBarry Smith PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 63e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 64e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 65e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 66e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 67e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 687b083b7cSBarry Smith PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 69ace3abfcSBarry Smith PetscBool free_diag; \ 70421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 71e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 724fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 73ace3abfcSBarry Smith PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 7417df9f7cSHong Zhang Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \ 7517df9f7cSHong Zhang means that this shares some data structures with the parent including diag, ilen, imax, i, j */ \ 76d67d9f35SJunchao Zhang Mat_SubSppt *submatis1; /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \ 77d67d9f35SJunchao Zhang Mat_SeqAIJOps ops[1] /* operations for SeqAIJ and its subclasses */ 78e6b907acSBarry Smith 7953565b12SHong Zhang typedef struct { 8053565b12SHong Zhang MatTransposeColoring matcoloring; 8153565b12SHong Zhang Mat Bt_den; /* dense matrix of B^T */ 8253565b12SHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 8353565b12SHong Zhang PetscBool usecoloring; 8453565b12SHong Zhang } Mat_MatMatTransMult; 8553565b12SHong Zhang 866d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */ 876d373c3eSHong Zhang Mat At; /* transpose of the first matrix */ 882cff0574SHong Zhang Mat mA; /* maij matrix of A */ 892cff0574SHong Zhang Vec bt, ct; /* vectors to hold locally transposed arrays of B and C */ 906718818eSStefano Zampini /* used by PtAP */ 916718818eSStefano Zampini void *data; 926718818eSStefano Zampini PetscErrorCode (*destroy)(void *); 932cff0574SHong Zhang } Mat_MatTransMatMult; 942cff0574SHong Zhang 9553565b12SHong Zhang typedef struct { 9653565b12SHong Zhang PetscInt *api, *apj; /* symbolic structure of A*P */ 9753565b12SHong Zhang PetscScalar *apa; /* temporary array for storing one row of A*P */ 983cdca5ebSHong Zhang } Mat_AP; 9953565b12SHong Zhang 10053565b12SHong Zhang typedef struct { 10153565b12SHong Zhang MatTransposeColoring matcoloring; 102257c235dSHong Zhang Mat Rt; /* sparse or dense matrix of R^T */ 10353565b12SHong Zhang Mat RARt; /* dense matrix of R*A*R^T */ 1043b1b9624SHong Zhang Mat ARt; /* A*R^T used for the case -matrart_color_art */ 10553565b12SHong Zhang MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 1066718818eSStefano Zampini /* free intermediate products needed for PtAP */ 1076718818eSStefano Zampini void *data; 1086718818eSStefano Zampini PetscErrorCode (*destroy)(void *); 10953565b12SHong Zhang } Mat_RARt; 11053565b12SHong Zhang 1116d0b6147SHong Zhang typedef struct { 1126d0b6147SHong Zhang Mat BC; /* temp matrix for storing B*C */ 1136d0b6147SHong Zhang } Mat_MatMatMatMult; 1146d0b6147SHong Zhang 1152d40f771SBarry Smith /* 116ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 117e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 118dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 1195768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 1202d40f771SBarry Smith */ 121d35516d3SLois Curfman McInnes 122e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 123b8a66259SBarry Smith typedef struct { 1244108e4d5SBarry Smith MatScalar *bdiag, *ibdiag, *ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 125f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 126ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 127f0d39aaaSBarry Smith 128ace3abfcSBarry Smith PetscBool use; 129e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 130e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 131e6b907acSBarry Smith PetscInt limit; /* inode limit */ 132e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 133ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 134a02bda8eSBarry Smith PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 1354108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 136e6b907acSBarry Smith 1375a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat, PetscViewer); 1385a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat, MatAssemblyType); 1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 1405a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 1415a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat, MatOption, PetscBool); 1425a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat, MatDuplicateOption, Mat *); 1435a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat, Mat, MatDuplicateOption, PetscBool); 1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat, Mat, const MatFactorInfo *); 145f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat, PetscScalar **); 146f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat, PetscScalar **); 147e6b907acSBarry Smith 148e6b907acSBarry Smith typedef struct { 14954f21887SBarry Smith SEQAIJHEADER(MatScalar); 1504108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 15154f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 15271f1c65dSBarry Smith 15371f1c65dSBarry Smith PetscScalar *idiag, *mdiag, *ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 154ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 155bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 1562291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 15761ecd0c6SBarry Smith PetscBool diagonaldense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */ 15871f1c65dSBarry Smith PetscScalar fshift, omega; /* last used omega and fshift */ 159394ed5ebSJunchao Zhang 16026cec326SBarry Smith /* MatSetValues() via hash related fields */ 16126cec326SBarry Smith PetscHMapIJV ht; 16226cec326SBarry Smith PetscInt *dnz; 16326cec326SBarry Smith struct _MatOps cops; 164ec8511deSBarry Smith } Mat_SeqAIJ; 165b8a66259SBarry Smith 1662c4ab24aSJunchao Zhang typedef struct { 1672c4ab24aSJunchao Zhang PetscInt nz; /* nz of the matrix after assembly */ 1682c4ab24aSJunchao Zhang PetscCount n; /* Number of entries in MatSetPreallocationCOO() */ 1692c4ab24aSJunchao Zhang PetscCount Atot; /* Total number of valid (i.e., w/ non-negative indices) entries in the COO array */ 1702c4ab24aSJunchao Zhang PetscCount *jmap; /* perm[jmap[i]..jmap[i+1]) give indices of entries in v[] associated with i-th nonzero of the matrix */ 1712c4ab24aSJunchao Zhang PetscCount *perm; /* The permutation array in sorting (i,j) by row and then by col */ 1722c4ab24aSJunchao Zhang } MatCOOStruct_SeqAIJ; 1732c4ab24aSJunchao Zhang 174c508b908SBarry Smith #define MatSeqXAIJGetOptions_Private(A) \ 175c508b908SBarry Smith { \ 176c508b908SBarry Smith const PetscBool oldvalues = (PetscBool)(A != PETSC_NULLPTR); \ 177c508b908SBarry Smith PetscInt nonew = 0, nounused = 0; \ 178c508b908SBarry Smith PetscBool roworiented = PETSC_FALSE; \ 179c508b908SBarry Smith if (oldvalues) { \ 180c508b908SBarry Smith nonew = ((Mat_SeqAIJ *)A->data)->nonew; \ 181c508b908SBarry Smith nounused = ((Mat_SeqAIJ *)A->data)->nounused; \ 182c508b908SBarry Smith roworiented = ((Mat_SeqAIJ *)A->data)->roworiented; \ 183f5729728SPierre Jolivet } \ 184f5729728SPierre Jolivet (void)0 185c508b908SBarry Smith 186c508b908SBarry Smith #define MatSeqXAIJRestoreOptions_Private(A) \ 187c508b908SBarry Smith if (oldvalues) { \ 188c508b908SBarry Smith ((Mat_SeqAIJ *)A->data)->nonew = nonew; \ 189c508b908SBarry Smith ((Mat_SeqAIJ *)A->data)->nounused = nounused; \ 190c508b908SBarry Smith ((Mat_SeqAIJ *)A->data)->roworiented = roworiented; \ 191c508b908SBarry Smith } \ 192f5729728SPierre Jolivet } \ 193f5729728SPierre Jolivet (void)0 194c508b908SBarry Smith 1959f0612e4SBarry Smith static inline PetscErrorCode MatXAIJAllocatea(Mat A, PetscInt nz, PetscScalar **array) 1969f0612e4SBarry Smith { 1979f0612e4SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1989f0612e4SBarry Smith 1999f0612e4SBarry Smith PetscFunctionBegin; 2009f0612e4SBarry Smith PetscCall(PetscShmgetAllocateArray(nz, sizeof(PetscScalar), (void **)array)); 2019f0612e4SBarry Smith a->free_a = PETSC_TRUE; 2029f0612e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2039f0612e4SBarry Smith } 2049f0612e4SBarry Smith 2059f0612e4SBarry Smith static inline PetscErrorCode MatXAIJDeallocatea(Mat A, PetscScalar **array) 2069f0612e4SBarry Smith { 2079f0612e4SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 2089f0612e4SBarry Smith 2099f0612e4SBarry Smith PetscFunctionBegin; 2109f0612e4SBarry Smith if (a->free_a) PetscCall(PetscShmgetDeallocateArray((void **)array)); 2119f0612e4SBarry Smith a->free_a = PETSC_FALSE; 2129f0612e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 2139f0612e4SBarry Smith } 2149f0612e4SBarry Smith 215e6b907acSBarry Smith /* 216e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 217e6b907acSBarry Smith */ 218d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA, MatScalar **a, PetscInt **j, PetscInt **i) 219d71ae5a4SJacob Faibussowitsch { 220e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ *)AA->data; 2213ba16761SJacob Faibussowitsch 2223ba16761SJacob Faibussowitsch PetscFunctionBegin; 2239f0612e4SBarry Smith if (A->free_a) PetscCall(PetscShmgetDeallocateArray((void **)a)); 2249f0612e4SBarry Smith if (A->free_ij) PetscCall(PetscShmgetDeallocateArray((void **)j)); 2259f0612e4SBarry Smith if (A->free_ij) PetscCall(PetscShmgetDeallocateArray((void **)i)); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227e6b907acSBarry Smith } 228e6b907acSBarry Smith /* 229e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 230357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 231e6b907acSBarry Smith */ 232fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat, AM, BS2, NROW, ROW, COL, RMAX, AA, AI, AJ, RP, AP, AIMAX, NONEW, datatype) \ 233a8f51744SPierre Jolivet do { \ 234fef13f97SBarry Smith if (NROW >= RMAX) { \ 235fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 236f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 237fef13f97SBarry Smith datatype *new_a; \ 238fef13f97SBarry Smith \ 23900045ab3SPierre Jolivet 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); \ 240fef13f97SBarry Smith /* malloc new storage space */ \ 2419f0612e4SBarry Smith PetscCall(PetscShmgetAllocateArray(BS2 *new_nz, sizeof(PetscScalar), (void **)&new_a)); \ 2429f0612e4SBarry Smith PetscCall(PetscShmgetAllocateArray(new_nz, sizeof(PetscInt), (void **)&new_j)); \ 2439f0612e4SBarry Smith PetscCall(PetscShmgetAllocateArray(AM + 1, sizeof(PetscInt), (void **)&new_i)); \ 2449f0612e4SBarry Smith Ain->free_a = PETSC_TRUE; \ 2459f0612e4SBarry Smith Ain->free_ij = PETSC_TRUE; \ 246fef13f97SBarry Smith /* copy over old data into new slots */ \ 247ad540459SPierre Jolivet for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 248ad540459SPierre Jolivet for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 2499566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 250fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 2518e3a54c0SPierre Jolivet PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, PetscSafePointerPlusOffset(AJ, AI[ROW] + NROW), len)); \ 2529566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_a, AA, BS2 *(AI[ROW] + NROW))); \ 2539566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(new_a + BS2 * (AI[ROW] + NROW), BS2 * CHUNKSIZE)); \ 2548e3a54c0SPierre Jolivet PetscCall(PetscArraycpy(new_a + BS2 * (AI[ROW] + NROW + CHUNKSIZE), PetscSafePointerPlusOffset(AA, BS2 * (AI[ROW] + NROW)), BS2 * len)); \ 255fef13f97SBarry Smith /* free up old matrix storage */ \ 2569566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 257fef13f97SBarry Smith AA = new_a; \ 258fef13f97SBarry Smith Ain->a = (MatScalar *)new_a; \ 2599371c9d4SSatish Balay AI = Ain->i = new_i; \ 2609371c9d4SSatish Balay AJ = Ain->j = new_j; \ 261fef13f97SBarry Smith \ 2629371c9d4SSatish Balay RP = AJ + AI[ROW]; \ 2639371c9d4SSatish Balay AP = AA + BS2 * AI[ROW]; \ 264fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 265fef13f97SBarry Smith Ain->maxnz += BS2 * CHUNKSIZE; \ 266fef13f97SBarry Smith Ain->reallocs++; \ 267b758ae8cSBarry Smith Amat->nonzerostate++; \ 268a8f51744SPierre Jolivet } \ 269a8f51744SPierre Jolivet } while (0) 27017454e89SShri Abhyankar 271876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat, AM, BS2, NROW, ROW, COL, RMAX, AI, AJ, RP, AIMAX, NONEW, datatype) \ 272a8f51744SPierre Jolivet do { \ 273720833daSHong Zhang if (NROW >= RMAX) { \ 274720833daSHong Zhang Mat_SeqAIJ *Ain = (Mat_SeqAIJ *)Amat->data; \ 275720833daSHong Zhang /* there is no extra room in row, therefore enlarge */ \ 276f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15, new_nz = AI[AM] + CHUNKSIZE, len, *new_i = NULL, *new_j = NULL; \ 277720833daSHong Zhang \ 27800045ab3SPierre Jolivet 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); \ 279720833daSHong Zhang /* malloc new storage space */ \ 2809f0612e4SBarry Smith PetscCall(PetscShmgetAllocateArray(new_nz, sizeof(PetscInt), (void **)&new_j)); \ 2819f0612e4SBarry Smith PetscCall(PetscShmgetAllocateArray(AM + 1, sizeof(PetscInt), (void **)&new_i)); \ 2829f0612e4SBarry Smith Ain->free_a = PETSC_FALSE; \ 2839f0612e4SBarry Smith Ain->free_ij = PETSC_TRUE; \ 284720833daSHong Zhang \ 285720833daSHong Zhang /* copy over old data into new slots */ \ 286ad540459SPierre Jolivet for (ii = 0; ii < ROW + 1; ii++) new_i[ii] = AI[ii]; \ 287ad540459SPierre Jolivet for (ii = ROW + 1; ii < AM + 1; ii++) new_i[ii] = AI[ii] + CHUNKSIZE; \ 2889566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j, AJ, AI[ROW] + NROW)); \ 289720833daSHong Zhang len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 2909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_j + AI[ROW] + NROW + CHUNKSIZE, AJ + AI[ROW] + NROW, len)); \ 291876c6284SHong Zhang \ 292720833daSHong Zhang /* free up old matrix storage */ \ 2939566063dSJacob Faibussowitsch PetscCall(MatSeqXAIJFreeAIJ(A, &Ain->a, &Ain->j, &Ain->i)); \ 294876c6284SHong Zhang Ain->a = NULL; \ 2959371c9d4SSatish Balay AI = Ain->i = new_i; \ 2969371c9d4SSatish Balay AJ = Ain->j = new_j; \ 297720833daSHong Zhang \ 298876c6284SHong Zhang RP = AJ + AI[ROW]; \ 299720833daSHong Zhang RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 300720833daSHong Zhang Ain->maxnz += BS2 * CHUNKSIZE; \ 301720833daSHong Zhang Ain->reallocs++; \ 302b758ae8cSBarry Smith Amat->nonzerostate++; \ 303a8f51744SPierre Jolivet } \ 304a8f51744SPierre Jolivet } while (0) 305e6b907acSBarry Smith 306cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat, PetscInt, const PetscInt *); 307e8729f6fSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJ(Mat, PetscCount, PetscInt[], PetscInt[]); 308cbc6b225SStefano Zampini 3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat, Mat, IS, IS, const MatFactorInfo *); 3111df811f5SHong Zhang 3125a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 3135a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat, Mat, IS, const MatFactorInfo *); 3145a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 3155a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 3165a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat, MatDuplicateOption, Mat *); 3175a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat, Mat, MatStructure); 3185a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat, PetscBool *, PetscInt *); 3195a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 3205a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat, PetscInt *, PetscInt **); 32108480c60SBarry Smith 322b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat, Vec, Vec); 323b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat, Vec, Vec); 324b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat, Vec, Vec, Vec); 325b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat, Vec, Vec, Vec); 326b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat, Vec, Vec); 327b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 3285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 329b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 33008480c60SBarry Smith 3315a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat, MatOption, PetscBool); 3323a7fca6bSBarry Smith 3335a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 3345a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat, PetscInt *[], PetscInt *[]); 3357fb60732SBarry Smith PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat, PetscInt, PetscInt, PetscInt *[], PetscInt *[]); 3365a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat, Mat *); 3375008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat, MatReuse, Mat *); 3387fb60732SBarry Smith 3392462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt, PetscInt *, PetscInt *, PetscBool, PetscInt, PetscInt, PetscInt **, PetscInt **); 3405a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat, Mat, IS, IS, const MatFactorInfo *); 3415a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat, Mat, const MatFactorInfo *); 3425a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat, Mat, const MatFactorInfo *); 3435a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat, Mat, const MatFactorInfo *); 3445a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat, IS, IS, const MatFactorInfo *); 3455a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat, Vec, Vec); 3465a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat, Vec, Vec); 3475a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat, Vec, Vec); 3485a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat, Vec, Vec); 3495a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat, Vec, Vec, Vec); 3505a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat, Vec, Vec); 3515a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat, Vec, Vec); 3525a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat, Vec, Vec, Vec); 3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat, Vec, Vec, Vec); 3545a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat, Mat, Mat); 355a3d9026eSPierre Jolivet PETSC_INTERN PetscErrorCode MatMatSolveTranspose_SeqAIJ(Mat, Mat, Mat); 35679c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat, Mat, PetscBool *); 35793dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat, ISColoring, MatFDColoring); 358f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat, ISColoring, MatFDColoring); 359a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat, MatFDColoring, PetscInt); 36052f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat, PetscViewer); 36152f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat, PetscViewer); 3625a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat, PetscViewer); 3635a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 3647bab7c10SHong Zhang 365df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE) 3664222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat); 367df97dc6dSFande Kong #endif 3684222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat); 3694222ddf1SHong Zhang 3704222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat); 3714222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat); 3724222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat); 3734222ddf1SHong Zhang 3744222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3754222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, PetscReal, Mat); 3764222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat, Mat, PetscReal, Mat); 3774222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, PetscReal, Mat); 3784222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat, Mat, PetscReal, Mat); 3794222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat, Mat, PetscReal, Mat); 3804222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat, Mat, PetscReal, Mat); 3814222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat, Mat, PetscReal, Mat); 3824222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat, Mat, PetscReal, Mat); 3834222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE) 3844222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat); 3854222ddf1SHong Zhang #endif 3864222ddf1SHong Zhang 3875a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 388df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat, Mat, Mat); 3894222ddf1SHong Zhang 3904099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat, Mat, Mat); 3915a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat, Mat, Mat); 3922b8ad9a3SHong Zhang 3934222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, PetscReal, Mat); 3945a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 3955a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat, Mat, Mat); 39653565b12SHong Zhang 3974222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 3984222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, PetscReal, Mat); 3994222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, PetscReal, Mat); 4005a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 40155bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat, Mat, Mat); 40255bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat, Mat, Mat); 4035df89d91SHong Zhang 4044222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 4055a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 4066718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *); 4073bf78175SHong Zhang 4084222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat, Mat, PetscReal, Mat); 4095a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat, Mat, Mat); 4105a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat, ISColoring, MatTransposeColoring); 4115a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring, Mat, Mat); 4125a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring, Mat, Mat); 4135df89d91SHong Zhang 4144222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, PetscReal, Mat); 4155a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat, Mat, Mat, Mat); 4167bab7c10SHong Zhang 417679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat, PetscInt, PetscInt, PetscRandom); 4185a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 4195a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 4205a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 421db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat, PetscScalar); 42287c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat, Vec, Vec); 42387c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat, Vec, InsertMode); 4245a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat, PetscScalar, Mat, MatStructure); 4255a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 4265a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 4275a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 4285a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 4296378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 4306378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscInt *[], PetscBool *); 4315a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 4325a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat, PetscViewer); 4339af31e4aSHong Zhang 4345a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 4355a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 436a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 437a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 438019b515eSShri Abhyankar 4395a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat, Mat, PetscInt *); 4409f5f6813SShri Abhyankar 441d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB) 442388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject, void *); 443388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject, void *); 444f2fbf96bSVaclav Hapla #endif 445cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType, MatReuse, Mat *); 446cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType, MatReuse, Mat *); 447388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat, MatType, MatReuse, Mat *); 448388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat, MatType, MatReuse, Mat *); 449388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat, MatType, MatReuse, Mat *); 450d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK) 451d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat, MatType, MatReuse, Mat *); 452d24d4204SJose E. Roman #endif 453388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat, MatType, MatReuse, Mat *); 454cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat, MatType, MatReuse, Mat *); 4554dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat, MatType, MatReuse, Mat *); 4564a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat, MatType, MatReuse, Mat *); 457388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat, MatType, MatReuse, Mat *); 4585a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat, PetscReal, IS, IS); 4595a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat, Mat, MatReuse, PetscReal, Mat *); 4608cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 4615a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat, MatAssemblyType); 4623fa6b06aSMark Adams PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat); 4632f947c57SVictor Minden 464b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt, const PetscInt *, const PetscInt *, const PetscInt *, const PetscInt *, PetscInt *); 4659c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 4669c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 46770f19b1fSKris Buschelman 46853dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat, IS, IS, MatStructure, Mat); 46958c11ad4SPierre Jolivet PETSC_INTERN PetscErrorCode MatEliminateZeros_SeqAIJ(Mat, PetscBool); 470f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *); 4710fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat); 472f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat); 47363a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat *[]); 474feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat, IS, IS, PetscInt, MatReuse, Mat *); 47553dd7562SDmitry Karpeev 476e4e71118SStefano Zampini PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], MatType, Mat); 477a3bb6f32SFande Kong 4789de2952eSStefano Zampini PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat, ISLocalToGlobalMapping *); 4799de2952eSStefano Zampini 480003131ecSBarry Smith /* 481003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 482003131ecSBarry Smith 483003131ecSBarry Smith Input Parameters: 484003131ecSBarry Smith + nnz - the number of entries 485003131ecSBarry Smith . r - the array of vector values 486003131ecSBarry Smith . xv - the matrix values for the row 487003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 488003131ecSBarry Smith 489003131ecSBarry Smith Output Parameter: 490003131ecSBarry Smith . sum - negative the sum of results 491003131ecSBarry Smith 492003131ecSBarry Smith PETSc compile flags: 4937b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4947b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4957b42bb93SJunchao Zhang 49611a5261eSBarry Smith Developer Note: 4977b42bb93SJunchao Zhang The macro changes sum but not other parameters 498003131ecSBarry Smith 499db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()` 500003131ecSBarry Smith */ 501519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 5029371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 503a8f51744SPierre Jolivet do { \ 504003131ecSBarry Smith if (nnz > 0) { \ 5057b42bb93SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 5067b42bb93SJunchao Zhang switch (rem) { \ 507d71ae5a4SJacob Faibussowitsch case 3: \ 508d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 509d71ae5a4SJacob Faibussowitsch case 2: \ 510d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 511d71ae5a4SJacob Faibussowitsch case 1: \ 512d71ae5a4SJacob Faibussowitsch sum -= *xv++ * r[*xi++]; \ 513d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 5147b42bb93SJunchao Zhang } \ 5159371c9d4SSatish Balay while (nnz2 > 0) { \ 5169371c9d4SSatish Balay sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 5179371c9d4SSatish Balay xv += 4; \ 5189371c9d4SSatish Balay xi += 4; \ 5199371c9d4SSatish Balay nnz2 -= 4; \ 5209371c9d4SSatish Balay } \ 5219371c9d4SSatish Balay xv -= nnz; \ 5229371c9d4SSatish Balay xi -= nnz; \ 5237b42bb93SJunchao Zhang } \ 524a8f51744SPierre Jolivet } while (0) 525003131ecSBarry Smith 526003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 5279371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 528a8f51744SPierre Jolivet do { \ 529003131ecSBarry Smith PetscInt __i, __i1, __i2; \ 5309371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 5319371c9d4SSatish Balay __i1 = xi[__i]; \ 5329371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 5339371c9d4SSatish Balay sum -= (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 5349371c9d4SSatish Balay } \ 5359371c9d4SSatish Balay if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]]; \ 536a8f51744SPierre Jolivet } while (0) 537003131ecSBarry Smith 538003131ecSBarry Smith #else 5399371c9d4SSatish Balay #define PetscSparseDenseMinusDot(sum, r, xv, xi, nnz) \ 540a8f51744SPierre Jolivet do { \ 541003131ecSBarry Smith PetscInt __i; \ 5429371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum -= xv[__i] * r[xi[__i]]; \ 543a8f51744SPierre Jolivet } while (0) 544003131ecSBarry Smith #endif 545003131ecSBarry Smith 546003131ecSBarry Smith /* 547003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 548003131ecSBarry Smith 549003131ecSBarry Smith Input Parameters: 550003131ecSBarry Smith + nnz - the number of entries 551003131ecSBarry Smith . r - the array of vector values 552003131ecSBarry Smith . xv - the matrix values for the row 553003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 554003131ecSBarry Smith 555003131ecSBarry Smith Output Parameter: 556003131ecSBarry Smith . sum - the sum of results 557003131ecSBarry Smith 558003131ecSBarry Smith PETSc compile flags: 5597b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 5607b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 5617b42bb93SJunchao Zhang 56211a5261eSBarry Smith Developer Note: 5637b42bb93SJunchao Zhang The macro changes sum but not other parameters 564003131ecSBarry Smith 565db781477SPatrick Sanan .seealso: `PetscSparseDenseMinusDot()` 566003131ecSBarry Smith */ 567519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 5689371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 569a8f51744SPierre Jolivet do { \ 570003131ecSBarry Smith if (nnz > 0) { \ 5717b42bb93SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 5727b42bb93SJunchao Zhang switch (rem) { \ 573d71ae5a4SJacob Faibussowitsch case 3: \ 574d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 575d71ae5a4SJacob Faibussowitsch case 2: \ 576d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 577d71ae5a4SJacob Faibussowitsch case 1: \ 578d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 579d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 5807b42bb93SJunchao Zhang } \ 5819371c9d4SSatish Balay while (nnz2 > 0) { \ 5829371c9d4SSatish Balay sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 5839371c9d4SSatish Balay xv += 4; \ 5849371c9d4SSatish Balay xi += 4; \ 5859371c9d4SSatish Balay nnz2 -= 4; \ 5869371c9d4SSatish Balay } \ 5879371c9d4SSatish Balay xv -= nnz; \ 5889371c9d4SSatish Balay xi -= nnz; \ 5897b42bb93SJunchao Zhang } \ 590a8f51744SPierre Jolivet } while (0) 591003131ecSBarry Smith 592003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 5939371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 594a8f51744SPierre Jolivet do { \ 595003131ecSBarry Smith PetscInt __i, __i1, __i2; \ 5969371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 5979371c9d4SSatish Balay __i1 = xi[__i]; \ 5989371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 5999371c9d4SSatish Balay sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 6009371c9d4SSatish Balay } \ 6019371c9d4SSatish Balay if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \ 602a8f51744SPierre Jolivet } while (0) 603003131ecSBarry Smith 6048d031ccaSJunchao Zhang #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) 60554e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum), (r), (xv), (xi), (nnz)) 606d68ff82bSRichard Tran Mills 607003131ecSBarry Smith #else 6089371c9d4SSatish Balay #define PetscSparseDensePlusDot(sum, r, xv, xi, nnz) \ 609a8f51744SPierre Jolivet do { \ 610003131ecSBarry Smith PetscInt __i; \ 6119371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \ 612a8f51744SPierre Jolivet } while (0) 613003131ecSBarry Smith #endif 614003131ecSBarry Smith 61599acd6aaSStefano Zampini #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) 61654e8760dSRichard Tran Mills #include <immintrin.h> 61754e8760dSRichard Tran Mills #if !defined(_MM_SCALE_8) 61854e8760dSRichard Tran Mills #define _MM_SCALE_8 8 61954e8760dSRichard Tran Mills #endif 62054e8760dSRichard Tran Mills 621d71ae5a4SJacob Faibussowitsch static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum, const PetscScalar *x, const MatScalar *aa, const PetscInt *aj, PetscInt n) 622d71ae5a4SJacob Faibussowitsch { 62354e8760dSRichard Tran Mills __m512d vec_x, vec_y, vec_vals; 62454e8760dSRichard Tran Mills __m256i vec_idx; 62554e8760dSRichard Tran Mills PetscInt j; 62654e8760dSRichard Tran Mills 62754e8760dSRichard Tran Mills vec_y = _mm512_setzero_pd(); 62854e8760dSRichard Tran Mills for (j = 0; j < (n >> 3); j++) { 629ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const *)aj); 630ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 63154e8760dSRichard Tran Mills vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8); 63254e8760dSRichard Tran Mills vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y); 6339371c9d4SSatish Balay aj += 8; 6349371c9d4SSatish Balay aa += 8; 63554e8760dSRichard Tran Mills } 636851da9e1SHong Zhang #if defined(__AVX512VL__) 637851da9e1SHong Zhang /* masked load requires avx512vl, which is not supported by KNL */ 638851da9e1SHong Zhang if (n & 0x07) { 63947a5f4f7SHong Zhang __mmask8 mask; 64054e8760dSRichard Tran Mills mask = (__mmask8)(0xff >> (8 - (n & 0x07))); 641851da9e1SHong Zhang vec_idx = _mm256_mask_loadu_epi32(vec_idx, mask, aj); 642851da9e1SHong Zhang vec_vals = _mm512_mask_loadu_pd(vec_vals, mask, aa); 64354e8760dSRichard Tran Mills vec_x = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8); 64454e8760dSRichard Tran Mills vec_y = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask); 64554e8760dSRichard Tran Mills } 646851da9e1SHong Zhang *sum += _mm512_reduce_add_pd(vec_y); 647851da9e1SHong Zhang #else 648851da9e1SHong Zhang *sum += _mm512_reduce_add_pd(vec_y); 64954e8760dSRichard Tran Mills for (j = 0; j < (n & 0x07); j++) *sum += aa[j] * x[aj[j]]; 650851da9e1SHong Zhang #endif 65154e8760dSRichard Tran Mills } 65254e8760dSRichard Tran Mills #endif 653b434eb95SMatthew G. Knepley 654b434eb95SMatthew G. Knepley /* 655b434eb95SMatthew G. Knepley PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 656b434eb95SMatthew G. Knepley 657b434eb95SMatthew G. Knepley Input Parameters: 658b434eb95SMatthew G. Knepley + nnz - the number of entries 659b434eb95SMatthew G. Knepley . r - the array of vector values 660b434eb95SMatthew G. Knepley . xv - the matrix values for the row 661b434eb95SMatthew G. Knepley - xi - the column indices of the nonzeros in the row 662b434eb95SMatthew G. Knepley 663b434eb95SMatthew G. Knepley Output Parameter: 664b434eb95SMatthew G. Knepley . max - the max of results 665b434eb95SMatthew G. Knepley 666db781477SPatrick Sanan .seealso: `PetscSparseDensePlusDot()`, `PetscSparseDenseMinusDot()` 667b434eb95SMatthew G. Knepley */ 6689371c9d4SSatish Balay #define PetscSparseDenseMaxDot(max, r, xv, xi, nnz) \ 669eec179cfSJacob Faibussowitsch do { \ 670eec179cfSJacob Faibussowitsch for (PetscInt __i = 0; __i < (nnz); __i++) { max = PetscMax(PetscRealPart(max), PetscRealPart((xv)[__i] * (r)[(xi)[__i]])); } \ 671eec179cfSJacob Faibussowitsch } while (0) 672b434eb95SMatthew G. Knepley 6734b38b95cSHong Zhang /* 6744b38b95cSHong Zhang Add column indices into table for counting the max nonzeros of merged rows 6754b38b95cSHong Zhang */ 6769371c9d4SSatish Balay #define MatRowMergeMax_SeqAIJ(mat, nrows, ta) \ 677eec179cfSJacob Faibussowitsch do { \ 678f4f49eeaSPierre Jolivet if (mat) { \ 679eec179cfSJacob Faibussowitsch for (PetscInt _row = 0; _row < (nrows); _row++) { \ 680eec179cfSJacob Faibussowitsch const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \ 681eec179cfSJacob Faibussowitsch for (PetscInt _j = 0; _j < _nz; _j++) { \ 682eec179cfSJacob Faibussowitsch PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \ 683c76ffc5fSJacob Faibussowitsch PetscCall(PetscHMapISet((ta), *_col + 1, 1)); \ 6844b38b95cSHong Zhang } \ 6854b38b95cSHong Zhang } \ 686ec07b8f8SHong Zhang } \ 687eec179cfSJacob Faibussowitsch } while (0) 6884b38b95cSHong Zhang 6890ca7d551SHong Zhang /* 6900ca7d551SHong Zhang Add column indices into table for counting the nonzeros of merged rows 6910ca7d551SHong Zhang */ 6929371c9d4SSatish Balay #define MatMergeRows_SeqAIJ(mat, nrows, rows, ta) \ 693eec179cfSJacob Faibussowitsch do { \ 694eec179cfSJacob Faibussowitsch for (PetscInt _i = 0; _i < (nrows); _i++) { \ 695eec179cfSJacob Faibussowitsch const PetscInt _row = (rows)[_i]; \ 696eec179cfSJacob Faibussowitsch const PetscInt _nz = (mat)->i[_row + 1] - (mat)->i[_row]; \ 697eec179cfSJacob Faibussowitsch for (PetscInt _j = 0; _j < _nz; _j++) { \ 698eec179cfSJacob Faibussowitsch PetscInt *_col = _j + (mat)->j + (mat)->i[_row]; \ 699eec179cfSJacob Faibussowitsch PetscCall(PetscHMapISetWithMode((ta), *_col + 1, 1, INSERT_VALUES)); \ 7000ca7d551SHong Zhang } \ 7010ca7d551SHong Zhang } \ 702eec179cfSJacob Faibussowitsch } while (0) 703