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 8e6b907acSBarry Smith /* 9e6b907acSBarry Smith Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 10e6b907acSBarry Smith */ 11421e10b8SBarry Smith #define SEQAIJHEADER(datatype) \ 12ace3abfcSBarry Smith PetscBool roworiented; /* if true, row-oriented input, default */ \ 13e6b907acSBarry Smith PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 1428b2fa4aSMatthew Knepley PetscInt nounused; /* -1 generate error on unused space */ \ 15ace3abfcSBarry Smith PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \ 16e6b907acSBarry Smith PetscInt maxnz; /* allocated nonzeros */ \ 17e6b907acSBarry Smith PetscInt *imax; /* maximum space allocated for each row */ \ 18e6b907acSBarry Smith PetscInt *ilen; /* actual length of each row */ \ 19846b4da1SFande Kong PetscInt *ipre; /* space preallocated for each row by user */ \ 20ace3abfcSBarry Smith PetscBool free_imax_ilen; \ 21e6b907acSBarry Smith PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 22e6b907acSBarry Smith as more values are set than were prealloced */\ 23e6b907acSBarry Smith PetscInt rmax; /* max nonzeros in any row */ \ 24ace3abfcSBarry Smith PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \ 25ace3abfcSBarry Smith PetscBool ignorezeroentries; \ 26ace3abfcSBarry Smith PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 27ace3abfcSBarry Smith PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 28e6b907acSBarry Smith Mat_CompressedRow compressedrow; /* use compressed row format */ \ 29e6b907acSBarry Smith PetscInt nz; /* nonzeros */ \ 30e6b907acSBarry Smith PetscInt *i; /* pointer to beginning of each row */ \ 31e6b907acSBarry Smith PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 32e6b907acSBarry Smith PetscInt *diag; /* pointers to diagonal elements */ \ 337b083b7cSBarry Smith PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 34ace3abfcSBarry Smith PetscBool free_diag; \ 35421e10b8SBarry Smith datatype *a; /* nonzero elements */ \ 36e6b907acSBarry Smith PetscScalar *solve_work; /* work space used in MatSolve */ \ 374fd072dbSBarry Smith IS row, col, icol; /* index sets, used for reorderings */ \ 38ace3abfcSBarry Smith PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 3917df9f7cSHong Zhang Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \ 4017df9f7cSHong Zhang means that this shares some data structures with the parent including diag, ilen, imax, i, j */\ 415c39f6d9SHong Zhang Mat_SubSppt *submatis1 /* used by MatCreateSubMatrices_MPIXAIJ_Local */ 42e6b907acSBarry Smith 4353565b12SHong Zhang typedef struct { 4453565b12SHong Zhang MatTransposeColoring matcoloring; 4553565b12SHong Zhang Mat Bt_den; /* dense matrix of B^T */ 4653565b12SHong Zhang Mat ABt_den; /* dense matrix of A*B^T */ 4753565b12SHong Zhang PetscBool usecoloring; 4853565b12SHong Zhang PetscErrorCode (*destroy)(Mat); 4953565b12SHong Zhang } Mat_MatMatTransMult; 5053565b12SHong Zhang 516d373c3eSHong Zhang typedef struct { /* used by MatTransposeMatMult() */ 526d373c3eSHong Zhang Mat At; /* transpose of the first matrix */ 532cff0574SHong Zhang Mat mA; /* maij matrix of A */ 542cff0574SHong Zhang Vec bt,ct; /* vectors to hold locally transposed arrays of B and C */ 552cff0574SHong Zhang PetscErrorCode (*destroy)(Mat); 562cff0574SHong Zhang } Mat_MatTransMatMult; 572cff0574SHong Zhang 5853565b12SHong Zhang typedef struct { 5953565b12SHong Zhang PetscInt *api,*apj; /* symbolic structure of A*P */ 6053565b12SHong Zhang PetscScalar *apa; /* temporary array for storing one row of A*P */ 6153565b12SHong Zhang PetscErrorCode (*destroy)(Mat); 623cdca5ebSHong Zhang } Mat_AP; 6353565b12SHong Zhang 6453565b12SHong Zhang typedef struct { 6553565b12SHong Zhang MatTransposeColoring matcoloring; 66257c235dSHong Zhang Mat Rt; /* sparse or dense matrix of R^T */ 6753565b12SHong Zhang Mat RARt; /* dense matrix of R*A*R^T */ 683b1b9624SHong Zhang Mat ARt; /* A*R^T used for the case -matrart_color_art */ 6953565b12SHong Zhang MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 7053565b12SHong Zhang PetscErrorCode (*destroy)(Mat); 7153565b12SHong Zhang } Mat_RARt; 7253565b12SHong Zhang 736d0b6147SHong Zhang typedef struct { 746d0b6147SHong Zhang Mat BC; /* temp matrix for storing B*C */ 756d0b6147SHong Zhang PetscErrorCode (*destroy)(Mat); 766d0b6147SHong Zhang } Mat_MatMatMatMult; 776d0b6147SHong Zhang 782d40f771SBarry Smith /* 79ec8511deSBarry Smith MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 80e6b907acSBarry Smith format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 81dfbc5765Svictorle j[i[k]+p] is the pth column in row k. Note that the diagonal 825768c4f9SLois Curfman McInnes matrix elements are stored with the rest of the nonzeros (not separately). 832d40f771SBarry Smith */ 84d35516d3SLois Curfman McInnes 85e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 86b8a66259SBarry Smith typedef struct { 874108e4d5SBarry Smith MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 88f0d39aaaSBarry Smith PetscInt bdiagsize; /* length of bdiag and ibdiag */ 89ace3abfcSBarry Smith PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 90f0d39aaaSBarry Smith 91ace3abfcSBarry Smith PetscBool use; 92e6b907acSBarry Smith PetscInt node_count; /* number of inodes */ 93e6b907acSBarry Smith PetscInt *size; /* size of each inode */ 94e6b907acSBarry Smith PetscInt limit; /* inode limit */ 95e6b907acSBarry Smith PetscInt max_limit; /* maximum supported inode limit */ 96ace3abfcSBarry Smith PetscBool checked; /* if inodes have been checked for */ 97a02bda8eSBarry Smith PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 984108e4d5SBarry Smith } Mat_SeqAIJ_Inode; 99e6b907acSBarry Smith 1005a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool); 1055a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 1065a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool); 1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 109*f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat,PetscScalar**); 110*f38c1e66SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat,PetscScalar**); 111e6b907acSBarry Smith 112e6b907acSBarry Smith typedef struct { 11354f21887SBarry Smith SEQAIJHEADER(MatScalar); 1144108e4d5SBarry Smith Mat_SeqAIJ_Inode inode; 11554f21887SBarry Smith MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 11671f1c65dSBarry Smith 11771f1c65dSBarry Smith PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 118ace3abfcSBarry Smith PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 119bbead8a2SBarry Smith PetscScalar *ibdiag; /* inverses of block diagonals */ 1202291e4fdSJed Brown PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 12161ecd0c6SBarry Smith PetscBool diagonaldense; /* all entries along the diagonal have been set; i.e. no missing diagonal terms */ 12271f1c65dSBarry Smith PetscScalar fshift,omega; /* last used omega and fshift */ 12371f1c65dSBarry Smith 1243a7fca6bSBarry Smith ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 1250b7e3e3dSHong Zhang 1260b7e3e3dSHong Zhang PetscScalar *matmult_abdense; /* used by MatMatMult() */ 1273cdca5ebSHong Zhang Mat_AP *ap; /* used by MatPtAP() */ 1286d0b6147SHong Zhang Mat_MatMatMatMult *matmatmatmult; /* used by MatMatMatMult() */ 129257c235dSHong Zhang Mat_RARt *rart; /* used by MatRARt() */ 13040192850SHong Zhang Mat_MatMatTransMult *abt; /* used by MatMatTransposeMult() */ 1316d373c3eSHong Zhang Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult() */ 132ec8511deSBarry Smith } Mat_SeqAIJ; 133b8a66259SBarry Smith 134e6b907acSBarry Smith /* 135e6b907acSBarry Smith Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 136e6b907acSBarry Smith */ 137d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 138d0f46423SBarry Smith { 139e6b907acSBarry Smith PetscErrorCode ierr; 140e6b907acSBarry Smith Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 141e6b907acSBarry Smith if (A->singlemalloc) { 142e6b907acSBarry Smith ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 143e6b907acSBarry Smith } else { 144c31cb41cSBarry Smith if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 145c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 146c31cb41cSBarry Smith if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 147e6b907acSBarry Smith } 148e6b907acSBarry Smith return 0; 149e6b907acSBarry Smith } 150e6b907acSBarry Smith /* 151e6b907acSBarry Smith Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 152357fed3dSBarry Smith This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 153e6b907acSBarry Smith */ 154fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 155fef13f97SBarry Smith if (NROW >= RMAX) { \ 156fef13f97SBarry Smith Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 157fef13f97SBarry Smith /* there is no extra room in row, therefore enlarge */ \ 158fef13f97SBarry Smith PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 159fef13f97SBarry Smith datatype *new_a; \ 160fef13f97SBarry Smith \ 1615aae435aSMatthew G. Knepley if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \ 162fef13f97SBarry Smith /* malloc new storage space */ \ 163dcca6d9dSJed Brown ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \ 164fef13f97SBarry Smith \ 165fef13f97SBarry Smith /* copy over old data into new slots */ \ 166fef13f97SBarry Smith for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 167fef13f97SBarry Smith for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 168580bdb30SBarry Smith ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \ 169fef13f97SBarry Smith len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 170580bdb30SBarry Smith ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \ 171580bdb30SBarry Smith ierr = PetscArraycpy(new_a,AA,BS2*(AI[ROW]+NROW));CHKERRQ(ierr); \ 172580bdb30SBarry Smith ierr = PetscArrayzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE);CHKERRQ(ierr); \ 173580bdb30SBarry Smith ierr = PetscArraycpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len);CHKERRQ(ierr); \ 174fef13f97SBarry Smith /* free up old matrix storage */ \ 175fef13f97SBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 176fef13f97SBarry Smith AA = new_a; \ 177fef13f97SBarry Smith Ain->a = (MatScalar*) new_a; \ 178fef13f97SBarry Smith AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 179fef13f97SBarry Smith Ain->singlemalloc = PETSC_TRUE; \ 180fef13f97SBarry Smith \ 181fef13f97SBarry Smith RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 182fef13f97SBarry Smith RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 183fef13f97SBarry Smith Ain->maxnz += BS2*CHUNKSIZE; \ 184fef13f97SBarry Smith Ain->reallocs++; \ 185fef13f97SBarry Smith } \ 18617454e89SShri Abhyankar 187876c6284SHong Zhang #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \ 188720833daSHong Zhang if (NROW >= RMAX) { \ 189720833daSHong Zhang Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 190720833daSHong Zhang /* there is no extra room in row, therefore enlarge */ \ 191720833daSHong Zhang PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 192720833daSHong Zhang \ 193720833daSHong Zhang if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc\nUse MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE) to turn off this check",ROW,COL); \ 194720833daSHong Zhang /* malloc new storage space */ \ 195720833daSHong Zhang ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \ 196720833daSHong Zhang ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\ 197720833daSHong Zhang \ 198720833daSHong Zhang /* copy over old data into new slots */ \ 199720833daSHong Zhang for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 200720833daSHong Zhang for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 201580bdb30SBarry Smith ierr = PetscArraycpy(new_j,AJ,AI[ROW]+NROW);CHKERRQ(ierr); \ 202720833daSHong Zhang len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 203580bdb30SBarry Smith ierr = PetscArraycpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len);CHKERRQ(ierr); \ 204876c6284SHong Zhang \ 205720833daSHong Zhang /* free up old matrix storage */ \ 206720833daSHong Zhang ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 207876c6284SHong Zhang Ain->a = NULL; \ 208720833daSHong Zhang AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 209720833daSHong Zhang Ain->singlemalloc = PETSC_FALSE; \ 210876c6284SHong Zhang Ain->free_a = PETSC_FALSE; \ 211720833daSHong Zhang \ 212876c6284SHong Zhang RP = AJ + AI[ROW]; \ 213720833daSHong Zhang RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 214720833daSHong Zhang Ain->maxnz += BS2*CHUNKSIZE; \ 215720833daSHong Zhang Ain->reallocs++; \ 216720833daSHong Zhang } \ 217e6b907acSBarry Smith 218cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 2195a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 2205a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 2221df811f5SHong Zhang 2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*); 2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**); 23408480c60SBarry Smith 2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 2375a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 2385a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 2395a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 24008480c60SBarry Smith 2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool); 2423a7fca6bSBarry Smith 2435a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 2445a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*); 2475008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*); 2482462f5fdSStefano Zampini PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**); 2495a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 2505a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 2525a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 2535a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 2545a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 2555a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 27079c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*); 27193dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring); 272f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring); 273a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt); 27452f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_AIJ_HDF5(Mat,PetscViewer); 27552f91c60SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ_Binary(Mat,PetscViewer); 2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 2775a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 2787bab7c10SHong Zhang 2795a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 2805a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 281df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,PetscReal,Mat*); 2824099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat*); 2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*); 2845a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*); 2855a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*); 2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*); 287d7ed1a4dSandi selinger PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat,Mat,PetscReal,Mat*); 288df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*); 289df97dc6dSFande Kong #if defined(PETSC_HAVE_HYPRE) 290df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 291df97dc6dSFande Kong #endif 292d013fc79Sandi selinger PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ_Combined(Mat,Mat,PetscReal,Mat*); 2935a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 294df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat,Mat,Mat); 2954099cc6bSBarry Smith PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat); 2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat); 297df97dc6dSFande Kong PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Combined(Mat,Mat,Mat); 2982b8ad9a3SHong Zhang 2995a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3005a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*); 3015a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 3025a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat); 30353565b12SHong Zhang 3045a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 30555bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*); 30655bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*); 3075a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 30855bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat); 30955bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat); 3105df89d91SHong Zhang 3115a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3125a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 3135a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 3146d373c3eSHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat); 3153bf78175SHong Zhang 3163bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*); 3173bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*); 3183bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat); 3193bf78175SHong Zhang 3205a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3215a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 3225a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 3235a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 3245a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 3255a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 3265df89d91SHong Zhang 3275a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 3285a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*); 3295a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat); 3307bab7c10SHong Zhang 331679944adSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat,PetscInt,PetscInt,PetscRandom); 3325a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3335a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 3345a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 335db63039fSRichard Tran Mills PETSC_INTERN PetscErrorCode MatScale_SeqAIJ(Mat,PetscScalar); 33687c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalScale_SeqAIJ(Mat,Vec,Vec); 33787c2a1d7SRichard Tran Mills PETSC_INTERN PetscErrorCode MatDiagonalSet_SeqAIJ(Mat,Vec,InsertMode); 3385a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 3395a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3405a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3415a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3425a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 3436378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 3446378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 3455a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 3465a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat); 3475a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 3489af31e4aSHong Zhang 3495a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 3505a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 351a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 352a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 353019b515eSShri Abhyankar 3545a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 3559f5f6813SShri Abhyankar 356388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 357f2fbf96bSVaclav Hapla #if defined(PETSC_HAVE_MATLAB_ENGINE) 358388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEnginePut_SeqAIJ(PetscObject,void*); 359388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatlabEngineGet_SeqAIJ(PetscObject,void*); 360f2fbf96bSVaclav Hapla #endif 361cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*); 362cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*); 363388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*); 364388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*); 365388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*); 366388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*); 367cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*); 3684dfdc2d9SRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJSELL(Mat,MatType,MatReuse,Mat*); 3694a2a386eSRichard Tran Mills PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJMKL(Mat,MatType,MatReuse,Mat*); 370388d47a6SSatish Balay PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJViennaCL(Mat,MatType,MatReuse,Mat*); 3715a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 3725a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3735a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 3748cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 3755a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 3762f947c57SVictor Minden 377b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 3789c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 3799c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 38070f19b1fSKris Buschelman 38153dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat); 382f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt*); 3830fb991dcSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat); 384f68bb481SHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrix_Dummy(Mat); 38563a75b2aSHong Zhang PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Dummy(PetscInt, Mat*[]); 386feb78a15SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 38753dd7562SDmitry Karpeev 388a3bb6f32SFande Kong PETSC_INTERN PetscErrorCode MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat,ISLocalToGlobalMapping*); 389a3bb6f32SFande Kong 390003131ecSBarry Smith /* 391003131ecSBarry Smith PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 392003131ecSBarry Smith 393003131ecSBarry Smith Input Parameters: 394003131ecSBarry Smith + nnz - the number of entries 395003131ecSBarry Smith . r - the array of vector values 396003131ecSBarry Smith . xv - the matrix values for the row 397003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 398003131ecSBarry Smith 399003131ecSBarry Smith Output Parameter: 400003131ecSBarry Smith . sum - negative the sum of results 401003131ecSBarry Smith 402003131ecSBarry Smith PETSc compile flags: 4037b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4047b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4057b42bb93SJunchao Zhang 4067b42bb93SJunchao Zhang Developer Notes: 4077b42bb93SJunchao Zhang The macro changes sum but not other parameters 408003131ecSBarry Smith 409003131ecSBarry Smith .seealso: PetscSparseDensePlusDot() 410003131ecSBarry Smith 411003131ecSBarry Smith */ 412519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 413003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 414003131ecSBarry Smith if (nnz > 0) { \ 4157b42bb93SJunchao Zhang PetscInt nnz2=nnz,rem=nnz&0x3; \ 4167b42bb93SJunchao Zhang switch (rem) { \ 417003131ecSBarry Smith case 3: sum -= *xv++ *r[*xi++]; \ 418003131ecSBarry Smith case 2: sum -= *xv++ *r[*xi++]; \ 419003131ecSBarry Smith case 1: sum -= *xv++ *r[*xi++]; \ 4207b42bb93SJunchao Zhang nnz2 -= rem;} \ 4217b42bb93SJunchao Zhang while (nnz2 > 0) { \ 4227b42bb93SJunchao Zhang sum -= xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 4237b42bb93SJunchao Zhang xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 4247b42bb93SJunchao Zhang xv += 4; xi += 4; nnz2 -= 4; \ 4257b42bb93SJunchao Zhang } \ 4267b42bb93SJunchao Zhang xv -= nnz; xi -= nnz; \ 4277b42bb93SJunchao Zhang } \ 4287b42bb93SJunchao Zhang } 429003131ecSBarry Smith 430003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 431003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 432003131ecSBarry Smith PetscInt __i,__i1,__i2; \ 433003131ecSBarry Smith for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 434003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 435003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 436003131ecSBarry Smith 437003131ecSBarry Smith #else 438003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 439003131ecSBarry Smith PetscInt __i; \ 440003131ecSBarry Smith for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];} 441003131ecSBarry Smith #endif 442003131ecSBarry Smith 443003131ecSBarry Smith 444003131ecSBarry Smith 445003131ecSBarry Smith /* 446003131ecSBarry Smith PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 447003131ecSBarry Smith 448003131ecSBarry Smith Input Parameters: 449003131ecSBarry Smith + nnz - the number of entries 450003131ecSBarry Smith . r - the array of vector values 451003131ecSBarry Smith . xv - the matrix values for the row 452003131ecSBarry Smith - xi - the column indices of the nonzeros in the row 453003131ecSBarry Smith 454003131ecSBarry Smith Output Parameter: 455003131ecSBarry Smith . sum - the sum of results 456003131ecSBarry Smith 457003131ecSBarry Smith PETSc compile flags: 4587b42bb93SJunchao Zhang + PETSC_KERNEL_USE_UNROLL_4 4597b42bb93SJunchao Zhang - PETSC_KERNEL_USE_UNROLL_2 4607b42bb93SJunchao Zhang 4617b42bb93SJunchao Zhang Developer Notes: 4627b42bb93SJunchao Zhang The macro changes sum but not other parameters 463003131ecSBarry Smith 464003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot() 465003131ecSBarry Smith 466003131ecSBarry Smith */ 467519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4) 468003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 469003131ecSBarry Smith if (nnz > 0) { \ 4707b42bb93SJunchao Zhang PetscInt nnz2=nnz,rem=nnz&0x3; \ 4717b42bb93SJunchao Zhang switch (rem) { \ 472003131ecSBarry Smith case 3: sum += *xv++ *r[*xi++]; \ 473003131ecSBarry Smith case 2: sum += *xv++ *r[*xi++]; \ 474003131ecSBarry Smith case 1: sum += *xv++ *r[*xi++]; \ 4757b42bb93SJunchao Zhang nnz2 -= rem;} \ 4767b42bb93SJunchao Zhang while (nnz2 > 0) { \ 477003131ecSBarry Smith sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 478003131ecSBarry Smith xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 4797b42bb93SJunchao Zhang xv += 4; xi += 4; nnz2 -= 4; \ 4807b42bb93SJunchao Zhang } \ 4817b42bb93SJunchao Zhang xv -= nnz; xi -= nnz; \ 4827b42bb93SJunchao Zhang } \ 4837b42bb93SJunchao Zhang } 484003131ecSBarry Smith 485003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2) 486003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 487003131ecSBarry Smith PetscInt __i,__i1,__i2; \ 488003131ecSBarry Smith for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 489003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 490003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 491003131ecSBarry Smith 492ef4c6e81SRichard Tran Mills #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) 49354e8760dSRichard Tran Mills #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) PetscSparseDensePlusDot_AVX512_Private(&(sum),(r),(xv),(xi),(nnz)) 494d68ff82bSRichard Tran Mills 495003131ecSBarry Smith #else 496003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 497003131ecSBarry Smith PetscInt __i; \ 498003131ecSBarry Smith for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];} 499003131ecSBarry Smith #endif 500003131ecSBarry Smith 501ef4c6e81SRichard Tran Mills #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) 50254e8760dSRichard Tran Mills #include <immintrin.h> 50354e8760dSRichard Tran Mills #if !defined(_MM_SCALE_8) 50454e8760dSRichard Tran Mills #define _MM_SCALE_8 8 50554e8760dSRichard Tran Mills #endif 50654e8760dSRichard Tran Mills 50754e8760dSRichard Tran Mills PETSC_STATIC_INLINE void PetscSparseDensePlusDot_AVX512_Private(PetscScalar *sum,const PetscScalar *x,const MatScalar *aa,const PetscInt *aj,PetscInt n) 50854e8760dSRichard Tran Mills { 50954e8760dSRichard Tran Mills __m512d vec_x,vec_y,vec_vals; 51054e8760dSRichard Tran Mills __m256i vec_idx; 51154e8760dSRichard Tran Mills __mmask8 mask; 51254e8760dSRichard Tran Mills PetscInt j; 51354e8760dSRichard Tran Mills 51454e8760dSRichard Tran Mills vec_y = _mm512_setzero_pd(); 51554e8760dSRichard Tran Mills for (j=0; j<(n>>3); j++) { 516ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const*)aj); 517ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 51854e8760dSRichard Tran Mills vec_x = _mm512_i32gather_pd(vec_idx,x,_MM_SCALE_8); 51954e8760dSRichard Tran Mills vec_y = _mm512_fmadd_pd(vec_x,vec_vals,vec_y); 52054e8760dSRichard Tran Mills aj += 8; aa += 8; 52154e8760dSRichard Tran Mills } 52254e8760dSRichard Tran Mills /* masked load does not work on KNL, it requires avx512vl */ 52354e8760dSRichard Tran Mills if ((n&0x07)>2) { 52454e8760dSRichard Tran Mills mask = (__mmask8)(0xff >> (8-(n&0x07))); 525ef588d5cSRichard Tran Mills vec_idx = _mm256_loadu_si256((__m256i const*)aj); 526ef588d5cSRichard Tran Mills vec_vals = _mm512_loadu_pd(aa); 52754e8760dSRichard Tran Mills vec_x = _mm512_mask_i32gather_pd(vec_x,mask,vec_idx,x,_MM_SCALE_8); 52854e8760dSRichard Tran Mills vec_y = _mm512_mask3_fmadd_pd(vec_x,vec_vals,vec_y,mask); 52954e8760dSRichard Tran Mills } else if ((n&0x07)==2) { 53054e8760dSRichard Tran Mills *sum += aa[0]*x[aj[0]]; 53154e8760dSRichard Tran Mills *sum += aa[1]*x[aj[1]]; 53254e8760dSRichard Tran Mills } else if ((n&0x07)==1) { 53354e8760dSRichard Tran Mills *sum += aa[0]*x[aj[0]]; 53454e8760dSRichard Tran Mills } 53554e8760dSRichard Tran Mills if (n>2) *sum += _mm512_reduce_add_pd(vec_y); 53654e8760dSRichard Tran Mills /* 53754e8760dSRichard Tran Mills for(j=0;j<(n&0x07);j++) *sum += aa[j]*x[aj[j]]; 53854e8760dSRichard Tran Mills */ 53954e8760dSRichard Tran Mills } 54054e8760dSRichard Tran Mills #endif 541b434eb95SMatthew G. Knepley 542b434eb95SMatthew G. Knepley /* 543b434eb95SMatthew G. Knepley PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 544b434eb95SMatthew G. Knepley 545b434eb95SMatthew G. Knepley Input Parameters: 546b434eb95SMatthew G. Knepley + nnz - the number of entries 547b434eb95SMatthew G. Knepley . r - the array of vector values 548b434eb95SMatthew G. Knepley . xv - the matrix values for the row 549b434eb95SMatthew G. Knepley - xi - the column indices of the nonzeros in the row 550b434eb95SMatthew G. Knepley 551b434eb95SMatthew G. Knepley Output Parameter: 552b434eb95SMatthew G. Knepley . max - the max of results 553b434eb95SMatthew G. Knepley 554b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot() 555b434eb95SMatthew G. Knepley 556b434eb95SMatthew G. Knepley */ 557b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \ 558b434eb95SMatthew G. Knepley PetscInt __i; \ 5595e792707SMatthew G. Knepley for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));} 560b434eb95SMatthew G. Knepley 5614b38b95cSHong Zhang /* 5624b38b95cSHong Zhang Add column indices into table for counting the max nonzeros of merged rows 5634b38b95cSHong Zhang */ 5644b38b95cSHong Zhang #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) { \ 5654b38b95cSHong Zhang PetscInt _j,_row,_nz,*_col; \ 566ec07b8f8SHong Zhang if (mat) { \ 5674b38b95cSHong Zhang for (_row=0; _row<nrows; _row++) { \ 5684b38b95cSHong Zhang _nz = mat->i[_row+1] - mat->i[_row]; \ 5694b38b95cSHong Zhang for (_j=0; _j<_nz; _j++) { \ 5704b38b95cSHong Zhang _col = _j + mat->j + mat->i[_row]; \ 5714b38b95cSHong Zhang PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \ 5724b38b95cSHong Zhang } \ 5734b38b95cSHong Zhang } \ 574ec07b8f8SHong Zhang } \ 5754b38b95cSHong Zhang } 5764b38b95cSHong Zhang 5770ca7d551SHong Zhang /* 5780ca7d551SHong Zhang Add column indices into table for counting the nonzeros of merged rows 5790ca7d551SHong Zhang */ 5800ca7d551SHong Zhang #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) { \ 5810ca7d551SHong Zhang PetscInt _j,_row,_nz,*_col,_i; \ 5820ca7d551SHong Zhang for (_i=0; _i<nrows; _i++) {\ 5830ca7d551SHong Zhang _row = rows[_i]; \ 5840ca7d551SHong Zhang _nz = mat->i[_row+1] - mat->i[_row]; \ 5850ca7d551SHong Zhang for (_j=0; _j<_nz; _j++) { \ 5860ca7d551SHong Zhang _col = _j + mat->j + mat->i[_row]; \ 5870ca7d551SHong Zhang PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \ 5880ca7d551SHong Zhang } \ 5890ca7d551SHong Zhang } \ 5900ca7d551SHong Zhang } 5910ca7d551SHong Zhang 5922d40f771SBarry Smith #endif 593