1b8a66259SBarry Smith 22d40f771SBarry Smith #if !defined(__AIJ_H) 32d40f771SBarry Smith #define __AIJ_H 4e6b907acSBarry Smith 5af0996ceSBarry Smith #include <petsc/private/matimpl.h> 6eca6b86aSHong Zhang #include <petscctable.h> 78f690400SShri Abhyankar 8d67d9f35SJunchao Zhang /* Operations provided by MATSEQAIJ and its subclasses */ 9d67d9f35SJunchao Zhang typedef struct { 10d67d9f35SJunchao Zhang PetscErrorCode (*getarray)(Mat,PetscScalar **); 11d67d9f35SJunchao Zhang PetscErrorCode (*restorearray)(Mat,PetscScalar **); 12d67d9f35SJunchao Zhang PetscErrorCode (*getarrayread)(Mat,const PetscScalar **); 13d67d9f35SJunchao Zhang PetscErrorCode (*restorearrayread)(Mat,const PetscScalar **); 14d67d9f35SJunchao Zhang PetscErrorCode (*getarraywrite)(Mat,PetscScalar **); 15d67d9f35SJunchao Zhang PetscErrorCode (*restorearraywrite)(Mat,PetscScalar **); 16d67d9f35SJunchao Zhang } Mat_SeqAIJOps; 17d67d9f35SJunchao Zhang 18e6b907acSBarry Smith /* 19e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 20e6b907acSBarry Smith */ 21421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 22ace3abfcSBarry Smith PetscBool roworiented; /* if true, row-oriented input, default */ \ 23e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 2428b2fa4aSMatthew Knepley PetscInt nounused; /* -1 generate error on unused space */ \ 25ace3abfcSBarry Smith PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \ 26e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */ \ 27e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */ \ 28e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */ \ 29846b4da1SFande Kong PetscInt *ipre; /* space preallocated for each row by user */ \ 30ace3abfcSBarry Smith PetscBool free_imax_ilen; \ 31e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 32e6b907acSBarry Smith as more values are set than were prealloced */\ 33e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */ \ 34ace3abfcSBarry Smith PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \ 35ace3abfcSBarry Smith PetscBool ignorezeroentries; \ 36ace3abfcSBarry Smith PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 37ace3abfcSBarry Smith PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 38e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 39e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 40e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 41e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 42e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 437b083b7cSBarry Smith PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 44ace3abfcSBarry Smith PetscBool free_diag; \ 45421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 46e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 474fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 48ace3abfcSBarry Smith PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 4917df9f7cSHong Zhang Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \ 5017df9f7cSHong Zhang means that this shares some data structures with the parent including diag, ilen, imax, i, j */\ 51d67d9f35SJunchao Zhang Mat_SubSppt *submatis1; /* used by MatCreateSubMatrices_MPIXAIJ_Local */ \ 52d67d9f35SJunchao Zhang Mat_SeqAIJOps ops[1] /* operations for SeqAIJ and its subclasses */ 53e6b907acSBarry Smith 5453565b12SHong Zhang typedef struct { 5553565b12SHong Zhang MatTransposeColoring matcoloring; 5653565b12SHong Zhang Mat Bt_den; /* dense matrix of B^T */ 5753565b12SHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 5853565b12SHong Zhang PetscBool usecoloring; 5953565b12SHong Zhang } Mat_MatMatTransMult; 6053565b12SHong Zhang 616d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */ 626d373c3eSHong Zhang Mat At; /* transpose of the first matrix */ 632cff0574SHong Zhang Mat mA; /* maij matrix of A */ 642cff0574SHong Zhang Vec bt,ct; /* vectors to hold locally transposed arrays of B and C */ 654222ddf1SHong Zhang PetscBool updateAt; /* flg to avoid recomputing At in MatProductNumeric_AtB_SeqAIJ_SeqAIJ() */ 666718818eSStefano Zampini /* used by PtAP */ 676718818eSStefano Zampini void *data; 686718818eSStefano Zampini PetscErrorCode (*destroy)(void*); 692cff0574SHong Zhang } Mat_MatTransMatMult; 702cff0574SHong Zhang 7153565b12SHong Zhang typedef struct { 7253565b12SHong Zhang PetscInt *api,*apj; /* symbolic structure of A*P */ 7353565b12SHong Zhang PetscScalar *apa; /* temporary array for storing one row of A*P */ 743cdca5ebSHong Zhang } Mat_AP; 7553565b12SHong Zhang 7653565b12SHong Zhang typedef struct { 7753565b12SHong Zhang MatTransposeColoring matcoloring; 78257c235dSHong Zhang Mat Rt; /* sparse or dense matrix of R^T */ 7953565b12SHong Zhang Mat RARt; /* dense matrix of R*A*R^T */ 803b1b9624SHong Zhang Mat ARt; /* A*R^T used for the case -matrart_color_art */ 8153565b12SHong Zhang MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 826718818eSStefano Zampini /* free intermediate products needed for PtAP */ 836718818eSStefano Zampini void *data; 846718818eSStefano Zampini PetscErrorCode (*destroy)(void*); 8553565b12SHong Zhang } Mat_RARt; 8653565b12SHong Zhang 876d0b6147SHong Zhang typedef struct { 886d0b6147SHong Zhang Mat BC; /* temp matrix for storing B*C */ 896d0b6147SHong Zhang } Mat_MatMatMatMult; 906d0b6147SHong Zhang 912d40f771SBarry Smith /* 92ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 93e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 94dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 955768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 962d40f771SBarry Smith */ 97d35516d3SLois Curfman McInnes 98e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 99b8a66259SBarry Smith typedef struct { 1004108e4d5SBarry Smith MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 101f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 102ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 103f0d39aaaSBarry Smith 104ace3abfcSBarry Smith PetscBool use; 105e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 106e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 107e6b907acSBarry Smith PetscInt limit; /* inode limit */ 108e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 109ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 110a02bda8eSBarry Smith PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 1114108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 112e6b907acSBarry Smith 1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 1175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool); 1185a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 1195a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool); 1205a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 1215a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 122f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat,PetscScalar**); 123f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat,PetscScalar**); 124e6b907acSBarry Smith 125e6b907acSBarry Smith typedef struct { 12654f21887SBarry Smith SEQAIJHEADER(MatScalar); 1274108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 12854f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 12971f1c65dSBarry Smith 13071f1c65dSBarry Smith PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 131ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 132bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 1332291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 13461ecd0c6SBarry Smith PetscBool diagonaldense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */ 13571f1c65dSBarry Smith PetscScalar fshift,omega; /* last used omega and fshift */ 136ec8511deSBarry Smith } Mat_SeqAIJ; 137b8a66259SBarry Smith 138e6b907acSBarry Smith /* 139e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 140e6b907acSBarry Smith */ 141*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 142d0f46423SBarry Smith { 143e6b907acSBarry Smith PetscErrorCode ierr; 144e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 145e6b907acSBarry Smith if (A->singlemalloc) { 146e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 147e6b907acSBarry Smith } else { 148c31cb41cSBarry Smith if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 149c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 150c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 151e6b907acSBarry Smith } 152e6b907acSBarry Smith return 0; 153e6b907acSBarry Smith } 154e6b907acSBarry Smith /* 155e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 156357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 157e6b907acSBarry Smith */ 158fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 159fef13f97SBarry Smith if (NROW >= RMAX) { \ 160fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 161fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 162f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=NULL,*new_j=NULL; \ 163fef13f97SBarry Smith datatype *new_a; \ 164fef13f97SBarry Smith \ 1659ace16cdSJacob Faibussowitsch PetscAssertFalse(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); \ 166fef13f97SBarry Smith /* malloc new storage space */ \ 167dcca6d9dSJed Brown ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \ 168fef13f97SBarry Smith \ 169fef13f97SBarry Smith /* copy over old data into new slots */ \ 170fef13f97SBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 171fef13f97SBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 172580bdb30SBarry Smith ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \ 173fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 174580bdb30SBarry Smith ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \ 175580bdb30SBarry Smith ierr = PetscArraycpy(new_a,AA,BS2*(AI[ROW]+NROW));CHKERRQ(ierr); \ 176580bdb30SBarry Smith ierr = PetscArrayzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE);CHKERRQ(ierr); \ 177580bdb30SBarry Smith ierr = PetscArraycpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len);CHKERRQ(ierr); \ 178fef13f97SBarry Smith /* free up old matrix storage */ \ 179fef13f97SBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 180fef13f97SBarry Smith AA = new_a; \ 181fef13f97SBarry Smith Ain->a = (MatScalar*) new_a; \ 182fef13f97SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 183fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 184fef13f97SBarry Smith \ 185fef13f97SBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 186fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 187fef13f97SBarry Smith Ain->maxnz += BS2*CHUNKSIZE; \ 188fef13f97SBarry Smith Ain->reallocs++; \ 189fef13f97SBarry Smith } \ 19017454e89SShri Abhyankar 191876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \ 192720833daSHong Zhang if (NROW >= RMAX) { \ 193720833daSHong Zhang Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 194720833daSHong Zhang /* there is no extra room in row, therefore enlarge */ \ 195f4259b30SLisandro Dalcin PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=NULL,*new_j=NULL; \ 196720833daSHong Zhang \ 1979ace16cdSJacob Faibussowitsch PetscAssertFalse(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); \ 198720833daSHong Zhang /* malloc new storage space */ \ 199720833daSHong Zhang ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \ 200720833daSHong Zhang ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\ 201720833daSHong Zhang \ 202720833daSHong Zhang /* copy over old data into new slots */ \ 203720833daSHong Zhang for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 204720833daSHong Zhang for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 205580bdb30SBarry Smith ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \ 206720833daSHong Zhang len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 207580bdb30SBarry Smith ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \ 208876c6284SHong Zhang \ 209720833daSHong Zhang /* free up old matrix storage */ \ 210720833daSHong Zhang ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 211876c6284SHong Zhang Ain->a = NULL; \ 212720833daSHong Zhang AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 213720833daSHong Zhang Ain->singlemalloc = PETSC_FALSE; \ 214876c6284SHong Zhang Ain->free_a = PETSC_FALSE; \ 215720833daSHong Zhang \ 216876c6284SHong Zhang RP = AJ + AI[ROW]; \ 217720833daSHong Zhang RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 218720833daSHong Zhang Ain->maxnz += BS2*CHUNKSIZE; \ 219720833daSHong Zhang Ain->reallocs++; \ 220720833daSHong Zhang } \ 221e6b907acSBarry Smith 222cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 2261df811f5SHong Zhang 2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*); 2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2375a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**); 23808480c60SBarry Smith 239b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat,Vec,Vec); 240b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqAIJ_Inode(Mat,Vec,Vec); 241b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat,Vec,Vec,Vec); 242b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ_Inode(Mat,Vec,Vec,Vec); 243b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat,Vec,Vec); 244b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 246b215bc84SStefano Zampini PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ_Inode(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 24708480c60SBarry Smith 2485a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool); 2493a7fca6bSBarry Smith 2505a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 2525a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 2535a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*); 2545008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*); 2552462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**); 2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 2705a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2715a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 2725a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 2735a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 2745a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2755a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 27779c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*); 27893dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring); 279f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring); 280a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt); 28152f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat,PetscViewer); 28252f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat,PetscViewer); 2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 2845a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 2857bab7c10SHong Zhang 286df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE) 2874222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Transpose_AIJ_AIJ(Mat); 288df97dc6dSFande Kong #endif 2894222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat); 2904222ddf1SHong Zhang 2914222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_SeqAIJ_SeqAIJ(Mat); 2924222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat); 2934222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat); 2944222ddf1SHong Zhang 2954222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat); 2964222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,PetscReal,Mat); 2974222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat); 2984222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat); 2994222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat); 3004222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat); 3014222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat); 3024222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat); 3034222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat); 3044222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE) 3054222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat); 3064222ddf1SHong Zhang #endif 3074222ddf1SHong Zhang 3085a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 309df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,Mat); 3104222ddf1SHong Zhang 3114099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat); 3125a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat); 3132b8ad9a3SHong Zhang 3144222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat); 3155a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 3165a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat); 31753565b12SHong Zhang 3184222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat); 3194222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat); 3204222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat); 3215a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 32255bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat); 32355bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat); 3245df89d91SHong Zhang 3254222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat); 3265a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 3276718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void*); 3283bf78175SHong Zhang 3294222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat); 3305a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 3315a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 3325a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 3335a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 3345df89d91SHong Zhang 3354222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat); 3365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat); 3377bab7c10SHong Zhang 338679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat,PetscInt,PetscInt,PetscRandom); 3395a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3405a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 3415a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 342db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar); 34387c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec); 34487c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode); 3455a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 3465a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3475a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3485a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3495a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3506378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 3516378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 3525a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 3535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat); 3545a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 3559af31e4aSHong Zhang 3565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 3575a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 358a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 359a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 360019b515eSShri Abhyankar 3615a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 3629f5f6813SShri Abhyankar 363f2fbf96bSVaclav Hapla #if defined(PETSC_HAVE_MATLAB_ENGINE) 364388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*); 365388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*); 366f2fbf96bSVaclav Hapla #endif 367cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*); 368cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*); 369388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*); 370388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*); 371388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*); 372d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK) 373d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*); 374d24d4204SJose E. Roman #endif 375388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*); 376cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*); 3774dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*); 3784a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*); 379388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*); 3805a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 3815a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 3835a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 3843fa6b06aSMark Adams PETSC_EXTERN PetscErrorCode MatZeroEntries_SeqAIJ(Mat); 3852f947c57SVictor Minden 386b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 3879c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 3889c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 38970f19b1fSKris Buschelman 39053dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat); 391f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*); 3920fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat); 393f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat); 39463a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]); 395feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 39653dd7562SDmitry Karpeev 397a3bb6f32SFande Kong PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat,ISLocalToGlobalMapping*); 398e4e71118SStefano Zampini PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],MatType,Mat); 399a3bb6f32SFande Kong 400003131ecSBarry Smith /* 401003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 402003131ecSBarry Smith 403003131ecSBarry Smith Input Parameters: 404003131ecSBarry Smith + nnz - the number of entries 405003131ecSBarry Smith . r - the array of vector values 406003131ecSBarry Smith . xv - the matrix values for the row 407003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 408003131ecSBarry Smith 409003131ecSBarry Smith Output Parameter: 410003131ecSBarry Smith . sum - negative the sum of results 411003131ecSBarry Smith 412003131ecSBarry Smith PETSc compile flags: 4137b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4147b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4157b42bb93SJunchao Zhang 4167b42bb93SJunchao Zhang Developer Notes: 4177b42bb93SJunchao Zhang The macro changes sum but not other parameters 418003131ecSBarry Smith 419003131ecSBarry Smith .seealso: PetscSparseDensePlusDot() 420003131ecSBarry Smith 421003131ecSBarry Smith */ 422519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 423003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 424003131ecSBarry Smith if (nnz > 0) { \ 4257b42bb93SJunchao Zhang PetscInt nnz2=nnz,rem=nnz&0x3; \ 4267b42bb93SJunchao Zhang switch (rem) { \ 427003131ecSBarry Smith case 3: sum -= *xv++ *r[*xi++]; \ 428003131ecSBarry Smith case 2: sum -= *xv++ *r[*xi++]; \ 429003131ecSBarry Smith case 1: sum -= *xv++ *r[*xi++]; \ 4307b42bb93SJunchao Zhang nnz2 -= rem;} \ 4317b42bb93SJunchao Zhang while (nnz2 > 0) { \ 4327b42bb93SJunchao Zhang sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 4337b42bb93SJunchao Zhang xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 4347b42bb93SJunchao Zhang xv += 4; xi += 4; nnz2 -= 4; \ 4357b42bb93SJunchao Zhang } \ 4367b42bb93SJunchao Zhang xv -= nnz; xi -= nnz; \ 4377b42bb93SJunchao Zhang } \ 4387b42bb93SJunchao Zhang } 439003131ecSBarry Smith 440003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 441003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 442003131ecSBarry Smith PetscInt __i,__i1,__i2; \ 443003131ecSBarry Smith for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 444003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 445003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 446003131ecSBarry Smith 447003131ecSBarry Smith #else 448003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 449003131ecSBarry Smith PetscInt __i; \ 450003131ecSBarry Smith for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];} 451003131ecSBarry Smith #endif 452003131ecSBarry Smith 453003131ecSBarry Smith /* 454003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 455003131ecSBarry Smith 456003131ecSBarry Smith Input Parameters: 457003131ecSBarry Smith + nnz - the number of entries 458003131ecSBarry Smith . r - the array of vector values 459003131ecSBarry Smith . xv - the matrix values for the row 460003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 461003131ecSBarry Smith 462003131ecSBarry Smith Output Parameter: 463003131ecSBarry Smith . sum - the sum of results 464003131ecSBarry Smith 465003131ecSBarry Smith PETSc compile flags: 4667b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4677b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4687b42bb93SJunchao Zhang 4697b42bb93SJunchao Zhang Developer Notes: 4707b42bb93SJunchao Zhang The macro changes sum but not other parameters 471003131ecSBarry Smith 472003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot() 473003131ecSBarry Smith 474003131ecSBarry Smith */ 475519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 476003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 477003131ecSBarry Smith if (nnz > 0) { \ 4787b42bb93SJunchao Zhang PetscInt nnz2=nnz,rem=nnz&0x3; \ 4797b42bb93SJunchao Zhang switch (rem) { \ 480003131ecSBarry Smith case 3: sum += *xv++ *r[*xi++]; \ 481003131ecSBarry Smith case 2: sum += *xv++ *r[*xi++]; \ 482003131ecSBarry Smith case 1: sum += *xv++ *r[*xi++]; \ 4837b42bb93SJunchao Zhang nnz2 -= rem;} \ 4847b42bb93SJunchao Zhang while (nnz2 > 0) { \ 485003131ecSBarry Smith sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 486003131ecSBarry Smith xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 4877b42bb93SJunchao Zhang xv += 4; xi += 4; nnz2 -= 4; \ 4887b42bb93SJunchao Zhang } \ 4897b42bb93SJunchao Zhang xv -= nnz; xi -= nnz; \ 4907b42bb93SJunchao Zhang } \ 4917b42bb93SJunchao Zhang } 492003131ecSBarry Smith 493003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 494003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 495003131ecSBarry Smith PetscInt __i,__i1,__i2; \ 496003131ecSBarry Smith for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 497003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 498003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 499003131ecSBarry Smith 50099acd6aaSStefano Zampini #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) 50154e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz)) 502d68ff82bSRichard Tran Mills 503003131ecSBarry Smith #else 504003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 505003131ecSBarry Smith PetscInt __i; \ 506003131ecSBarry Smith for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];} 507003131ecSBarry Smith #endif 508003131ecSBarry Smith 50999acd6aaSStefano 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) 51054e8760dSRichard Tran Mills #include <immintrin.h> 51154e8760dSRichard Tran Mills #if !defined(_MM_SCALE_8) 51254e8760dSRichard Tran Mills #define _MM_SCALE_8 8 51354e8760dSRichard Tran Mills #endif 51454e8760dSRichard Tran Mills 515*9fbee547SJacob Faibussowitsch static inline void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n) 51654e8760dSRichard Tran Mills { 51754e8760dSRichard Tran Mills __m512d vec_x,vec_y,vec_vals; 51854e8760dSRichard Tran Mills __m256i vec_idx; 51954e8760dSRichard Tran Mills __mmask8 mask; 52054e8760dSRichard Tran Mills PetscInt j; 52154e8760dSRichard Tran Mills 52254e8760dSRichard Tran Mills vec_y = _mm512_setzero_pd(); 52354e8760dSRichard Tran Mills for (j=0; j<(n>>3); j++) { 524ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const*)aj); 525ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 52654e8760dSRichard Tran Mills vec_x = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8); 52754e8760dSRichard Tran Mills vec_y = _mm512_fmadd_pd(vec_x,vec_vals,vec_y); 52854e8760dSRichard Tran Mills aj += 8; aa += 8; 52954e8760dSRichard Tran Mills } 53054e8760dSRichard Tran Mills /* masked load does not work on KNL, it requires avx512vl */ 53154e8760dSRichard Tran Mills if ((n&0x07)>2) { 53254e8760dSRichard Tran Mills mask = (__mmask8)(0xff >> (8-(n&0x07))); 533ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const*)aj); 534ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 53554e8760dSRichard Tran Mills vec_x = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8); 53654e8760dSRichard Tran Mills vec_y = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask); 53754e8760dSRichard Tran Mills } else if ((n&0x07)==2) { 53854e8760dSRichard Tran Mills *sum += aa[0]*x[aj[0]]; 53954e8760dSRichard Tran Mills *sum += aa[1]*x[aj[1]]; 54054e8760dSRichard Tran Mills } else if ((n&0x07)==1) { 54154e8760dSRichard Tran Mills *sum += aa[0]*x[aj[0]]; 54254e8760dSRichard Tran Mills } 54354e8760dSRichard Tran Mills if (n>2) *sum += _mm512_reduce_add_pd(vec_y); 54454e8760dSRichard Tran Mills /* 54554e8760dSRichard Tran Mills for (j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]]; 54654e8760dSRichard Tran Mills */ 54754e8760dSRichard Tran Mills } 54854e8760dSRichard Tran Mills #endif 549b434eb95SMatthew G. Knepley 550b434eb95SMatthew G. Knepley /* 551b434eb95SMatthew G. Knepley PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 552b434eb95SMatthew G. Knepley 553b434eb95SMatthew G. Knepley Input Parameters: 554b434eb95SMatthew G. Knepley + nnz - the number of entries 555b434eb95SMatthew G. Knepley . r - the array of vector values 556b434eb95SMatthew G. Knepley . xv - the matrix values for the row 557b434eb95SMatthew G. Knepley - xi - the column indices of the nonzeros in the row 558b434eb95SMatthew G. Knepley 559b434eb95SMatthew G. Knepley Output Parameter: 560b434eb95SMatthew G. Knepley . max - the max of results 561b434eb95SMatthew G. Knepley 562b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot() 563b434eb95SMatthew G. Knepley 564b434eb95SMatthew G. Knepley */ 565b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \ 566b434eb95SMatthew G. Knepley PetscInt __i; \ 5675e792707SMatthew G. Knepley for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));} 568b434eb95SMatthew G. Knepley 5694b38b95cSHong Zhang /* 5704b38b95cSHong Zhang Add column indices into table for counting the max nonzeros of merged rows 5714b38b95cSHong Zhang */ 5724b38b95cSHong Zhang #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) { \ 5734b38b95cSHong Zhang PetscInt _j,_row,_nz,*_col; \ 574ec07b8f8SHong Zhang if (mat) { \ 5754b38b95cSHong Zhang for (_row=0; _row<nrows; _row++) { \ 5764b38b95cSHong Zhang _nz = mat->i[_row+1] - mat->i[_row]; \ 5774b38b95cSHong Zhang for (_j=0; _j<_nz; _j++) { \ 5784b38b95cSHong Zhang _col = _j + mat->j + mat->i[_row]; \ 5794b38b95cSHong Zhang PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \ 5804b38b95cSHong Zhang } \ 5814b38b95cSHong Zhang } \ 582ec07b8f8SHong Zhang } \ 5834b38b95cSHong Zhang } 5844b38b95cSHong Zhang 5850ca7d551SHong Zhang /* 5860ca7d551SHong Zhang Add column indices into table for counting the nonzeros of merged rows 5870ca7d551SHong Zhang */ 5880ca7d551SHong Zhang #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) { \ 5890ca7d551SHong Zhang PetscInt _j,_row,_nz,*_col,_i; \ 5900ca7d551SHong Zhang for (_i=0; _i<nrows; _i++) {\ 5910ca7d551SHong Zhang _row = rows[_i]; \ 5920ca7d551SHong Zhang _nz = mat->i[_row+1] - mat->i[_row]; \ 5930ca7d551SHong Zhang for (_j=0; _j<_nz; _j++) { \ 5940ca7d551SHong Zhang _col = _j + mat->j + mat->i[_row]; \ 5950ca7d551SHong Zhang PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \ 5960ca7d551SHong Zhang } \ 5970ca7d551SHong Zhang } \ 5980ca7d551SHong Zhang } 5990ca7d551SHong Zhang 6002d40f771SBarry Smith #endif 601