1 2 #if !defined(__AIJ_H) 3 #define __AIJ_H 4 5 #include <petsc/private/matimpl.h> 6 #include <petscctable.h> 7 8 /* 9 Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats 10 */ 11 #define SEQAIJHEADER(datatype) \ 12 PetscBool roworiented; /* if true, row-oriented input, default */ \ 13 PetscInt nonew; /* 1 don't add new nonzeros, -1 generate error on new */ \ 14 PetscInt nounused; /* -1 generate error on unused space */ \ 15 PetscBool singlemalloc; /* if true a, i, and j have been obtained with one big malloc */ \ 16 PetscInt maxnz; /* allocated nonzeros */ \ 17 PetscInt *imax; /* maximum space allocated for each row */ \ 18 PetscInt *ilen; /* actual length of each row */ \ 19 PetscBool free_imax_ilen; \ 20 PetscInt reallocs; /* number of mallocs done during MatSetValues() \ 21 as more values are set than were prealloced */\ 22 PetscInt rmax; /* max nonzeros in any row */ \ 23 PetscBool keepnonzeropattern; /* keeps matrix structure same in calls to MatZeroRows()*/ \ 24 PetscBool ignorezeroentries; \ 25 PetscBool free_ij; /* free the column indices j and row offsets i when the matrix is destroyed */ \ 26 PetscBool free_a; /* free the numerical values when matrix is destroy */ \ 27 Mat_CompressedRow compressedrow; /* use compressed row format */ \ 28 PetscInt nz; /* nonzeros */ \ 29 PetscInt *i; /* pointer to beginning of each row */ \ 30 PetscInt *j; /* column values: j + i[k] - 1 is start of row k */ \ 31 PetscInt *diag; /* pointers to diagonal elements */ \ 32 PetscInt nonzerorowcnt; /* how many rows have nonzero entries */ \ 33 PetscBool free_diag; \ 34 datatype *a; /* nonzero elements */ \ 35 PetscScalar *solve_work; /* work space used in MatSolve */ \ 36 IS row, col, icol; /* index sets, used for reorderings */ \ 37 PetscBool pivotinblocks; /* pivot inside factorization of each diagonal block */ \ 38 Mat parent; /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....); \ 39 means that this shares some data structures with the parent including diag, ilen, imax, i, j */\ 40 Mat_SubSppt *submatis1 /* used by MatCreateSubMatrices_MPIXAIJ_Local */ 41 42 typedef struct { 43 MatTransposeColoring matcoloring; 44 Mat Bt_den; /* dense matrix of B^T */ 45 Mat ABt_den; /* dense matrix of A*B^T */ 46 PetscBool usecoloring; 47 PetscErrorCode (*destroy)(Mat); 48 } Mat_MatMatTransMult; 49 50 typedef struct { /* for MatTransposeMatMult_SeqAIJ_SeqDense() */ 51 Mat mA; /* maij matrix of A */ 52 Vec bt,ct; /* vectors to hold locally transposed arrays of B and C */ 53 PetscErrorCode (*destroy)(Mat); 54 } Mat_MatTransMatMult; 55 56 typedef struct { 57 PetscInt *api,*apj; /* symbolic structure of A*P */ 58 PetscScalar *apa; /* temporary array for storing one row of A*P */ 59 PetscErrorCode (*destroy)(Mat); 60 } Mat_PtAP; 61 62 typedef struct { 63 MatTransposeColoring matcoloring; 64 Mat Rt; /* sparse or dense matrix of R^T */ 65 Mat RARt; /* dense matrix of R*A*R^T */ 66 Mat ARt; /* A*R^T used for the case -matrart_color_art */ 67 MatScalar *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 68 PetscErrorCode (*destroy)(Mat); 69 } Mat_RARt; 70 71 typedef struct { 72 Mat BC; /* temp matrix for storing B*C */ 73 PetscErrorCode (*destroy)(Mat); 74 } Mat_MatMatMatMult; 75 76 /* 77 MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix 78 format) or compressed sparse row (CSR). The i[] and j[] arrays start at 0. For example, 79 j[i[k]+p] is the pth column in row k. Note that the diagonal 80 matrix elements are stored with the rest of the nonzeros (not separately). 81 */ 82 83 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */ 84 typedef struct { 85 MatScalar *bdiag,*ibdiag,*ssor_work; /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */ 86 PetscInt bdiagsize; /* length of bdiag and ibdiag */ 87 PetscBool ibdiagvalid; /* do ibdiag[] and bdiag[] contain the most recent values */ 88 89 PetscBool use; 90 PetscInt node_count; /* number of inodes */ 91 PetscInt *size; /* size of each inode */ 92 PetscInt limit; /* inode limit */ 93 PetscInt max_limit; /* maximum supported inode limit */ 94 PetscBool checked; /* if inodes have been checked for */ 95 PetscObjectState mat_nonzerostate; /* non-zero state when inodes were checked for */ 96 } Mat_SeqAIJ_Inode; 97 98 PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer); 99 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType); 100 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat); 101 PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat); 102 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool); 103 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*); 104 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool); 105 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*); 106 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*); 107 108 typedef struct { 109 SEQAIJHEADER(MatScalar); 110 Mat_SeqAIJ_Inode inode; 111 MatScalar *saved_values; /* location for stashing nonzero values of matrix */ 112 113 PetscScalar *idiag,*mdiag,*ssor_work; /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */ 114 PetscBool idiagvalid; /* current idiag[] and mdiag[] are valid */ 115 PetscScalar *ibdiag; /* inverses of block diagonals */ 116 PetscBool ibdiagvalid; /* inverses of block diagonals are valid. */ 117 PetscScalar fshift,omega; /* last used omega and fshift */ 118 119 ISColoring coloring; /* set with MatADSetColoring() used by MatADSetValues() */ 120 121 PetscScalar *matmult_abdense; /* used by MatMatMult() */ 122 Mat_PtAP *ptap; /* used by MatPtAP() */ 123 Mat_MatMatMatMult *matmatmatmult; /* used by MatMatMatMult() */ 124 Mat_RARt *rart; /* used by MatRARt() */ 125 Mat_MatMatTransMult *abt; /* used by MatMatTransposeMult() */ 126 } Mat_SeqAIJ; 127 128 /* 129 Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 130 */ 131 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) 132 { 133 PetscErrorCode ierr; 134 Mat_SeqAIJ *A = (Mat_SeqAIJ*) AA->data; 135 if (A->singlemalloc) { 136 ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr); 137 } else { 138 if (A->free_a) {ierr = PetscFree(*a);CHKERRQ(ierr);} 139 if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);} 140 if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);} 141 } 142 return 0; 143 } 144 /* 145 Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types 146 This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar 147 */ 148 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \ 149 if (NROW >= RMAX) { \ 150 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 151 /* there is no extra room in row, therefore enlarge */ \ 152 PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 153 datatype *new_a; \ 154 \ 155 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); \ 156 /* malloc new storage space */ \ 157 ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \ 158 \ 159 /* copy over old data into new slots */ \ 160 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 161 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 162 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 163 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 164 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 165 ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \ 166 ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \ 167 ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr); \ 168 /* free up old matrix storage */ \ 169 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 170 AA = new_a; \ 171 Ain->a = (MatScalar*) new_a; \ 172 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 173 Ain->singlemalloc = PETSC_TRUE; \ 174 \ 175 RP = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \ 176 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 177 Ain->maxnz += BS2*CHUNKSIZE; \ 178 Ain->reallocs++; \ 179 } \ 180 181 #define MatSeqXAIJReallocateAIJ_structure_only(Amat,AM,BS2,NROW,ROW,COL,RMAX,AI,AJ,RP,AIMAX,NONEW,datatype) \ 182 if (NROW >= RMAX) { \ 183 Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \ 184 /* there is no extra room in row, therefore enlarge */ \ 185 PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \ 186 \ 187 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); \ 188 /* malloc new storage space */ \ 189 ierr = PetscMalloc1(new_nz,&new_j);CHKERRQ(ierr); \ 190 ierr = PetscMalloc1(AM+1,&new_i);CHKERRQ(ierr);\ 191 \ 192 /* copy over old data into new slots */ \ 193 for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \ 194 for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \ 195 ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \ 196 len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \ 197 ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \ 198 \ 199 /* free up old matrix storage */ \ 200 ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \ 201 Ain->a = NULL; \ 202 AI = Ain->i = new_i; AJ = Ain->j = new_j; \ 203 Ain->singlemalloc = PETSC_FALSE; \ 204 Ain->free_a = PETSC_FALSE; \ 205 \ 206 RP = AJ + AI[ROW]; \ 207 RMAX = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \ 208 Ain->maxnz += BS2*CHUNKSIZE; \ 209 Ain->reallocs++; \ 210 } \ 211 212 PETSC_INTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*); 213 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 214 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 215 PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*); 216 217 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 218 PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 219 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*); 220 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*); 221 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 222 PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 223 PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*); 224 PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure); 225 PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*); 226 PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat); 227 PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**); 228 229 PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec); 230 PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 231 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec); 232 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec); 233 PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 234 235 PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool); 236 237 PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 238 PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]); 239 PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]); 240 PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*); 241 PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*); 242 PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscBool,PetscInt,PetscInt,PetscInt**,PetscInt**); 243 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*); 244 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*); 245 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*); 246 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*); 247 PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*); 248 PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*); 249 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec); 250 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec); 251 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec); 252 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec); 253 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec); 254 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec); 255 PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec); 256 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 257 PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 258 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec); 259 PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec); 260 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec); 261 PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec); 262 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat); 263 PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat); 264 PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*); 265 PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring); 266 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring); 267 PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt); 268 PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer); 269 PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat); 270 271 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 272 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 273 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat,Mat,PetscReal,Mat*); 274 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*); 275 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*); 276 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*); 277 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*); 278 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 279 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat,Mat,Mat); 280 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat); 281 282 PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 283 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat,Mat,PetscReal,Mat*); 284 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*); 285 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 286 PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat); 287 288 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 289 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*); 290 PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*); 291 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 292 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat); 293 PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat); 294 295 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 296 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 297 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 298 299 PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*); 300 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*); 301 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat); 302 303 PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 304 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*); 305 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat); 306 PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring); 307 PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat); 308 PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat); 309 310 PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 311 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*); 312 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat); 313 314 PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 315 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 316 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 317 PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure); 318 PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 319 PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 320 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 321 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*); 322 PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 323 PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*); 324 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 325 PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat); 326 PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer); 327 328 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat); 329 PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat); 330 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat); 331 PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat); 332 333 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*); 334 335 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*); 336 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*); 337 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*); 338 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS); 339 PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 340 PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 341 PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat); 342 PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType); 343 344 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 345 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 346 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 347 348 PETSC_INTERN PetscErrorCode MatSetSeqMat_SeqAIJ(Mat,IS,IS,MatStructure,Mat); 349 PETSC_INTERN PetscErrorCode MatDestroySubMatrices_Private(Mat_SubSppt*); 350 PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Submatrices(Mat); 351 PETSC_INTERN PetscErrorCode MatDestroy_Dummy_Submatrices(Mat); 352 PETSC_INTERN PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 353 354 /* 355 PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage 356 357 Input Parameters: 358 + nnz - the number of entries 359 . r - the array of vector values 360 . xv - the matrix values for the row 361 - xi - the column indices of the nonzeros in the row 362 363 Output Parameter: 364 . sum - negative the sum of results 365 366 PETSc compile flags: 367 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 368 - PETSC_KERNEL_USE_UNROLL_2 - 369 370 .seealso: PetscSparseDensePlusDot() 371 372 */ 373 #if defined(PETSC_KERNEL_USE_UNROLL_4) 374 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 375 if (nnz > 0) { \ 376 switch (nnz & 0x3) { \ 377 case 3: sum -= *xv++ *r[*xi++]; \ 378 case 2: sum -= *xv++ *r[*xi++]; \ 379 case 1: sum -= *xv++ *r[*xi++]; \ 380 nnz -= 4;} \ 381 while (nnz > 0) { \ 382 sum -= xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \ 383 xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \ 384 xv += 4; xi += 4; nnz -= 4; }}} 385 386 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 387 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 388 PetscInt __i,__i1,__i2; \ 389 for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 390 sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 391 if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];} 392 393 #else 394 #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \ 395 PetscInt __i; \ 396 for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];} 397 #endif 398 399 400 401 /* 402 PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage 403 404 Input Parameters: 405 + nnz - the number of entries 406 . r - the array of vector values 407 . xv - the matrix values for the row 408 - xi - the column indices of the nonzeros in the row 409 410 Output Parameter: 411 . sum - the sum of results 412 413 PETSc compile flags: 414 + PETSC_KERNEL_USE_UNROLL_4 - don't use this; it changes nnz and hence is WRONG 415 - PETSC_KERNEL_USE_UNROLL_2 - 416 417 .seealso: PetscSparseDenseMinusDot() 418 419 */ 420 #if defined(PETSC_KERNEL_USE_UNROLL_4) 421 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 422 if (nnz > 0) { \ 423 switch (nnz & 0x3) { \ 424 case 3: sum += *xv++ *r[*xi++]; \ 425 case 2: sum += *xv++ *r[*xi++]; \ 426 case 1: sum += *xv++ *r[*xi++]; \ 427 nnz -= 4;} \ 428 while (nnz > 0) { \ 429 sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 430 xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 431 xv += 4; xi += 4; nnz -= 4; }}} 432 433 #elif defined(PETSC_KERNEL_USE_UNROLL_2) 434 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 435 PetscInt __i,__i1,__i2; \ 436 for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 437 sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 438 if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 439 440 #else 441 #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \ 442 PetscInt __i; \ 443 for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];} 444 #endif 445 446 447 /* 448 PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage 449 450 Input Parameters: 451 + nnz - the number of entries 452 . r - the array of vector values 453 . xv - the matrix values for the row 454 - xi - the column indices of the nonzeros in the row 455 456 Output Parameter: 457 . max - the max of results 458 459 .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot() 460 461 */ 462 #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \ 463 PetscInt __i; \ 464 for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));} 465 466 /* 467 Add column indices into table for counting the max nonzeros of merged rows 468 */ 469 #define MatRowMergeMax_SeqAIJ(mat,nrows,ta) { \ 470 PetscInt _j,_row,_nz,*_col; \ 471 if (mat) { \ 472 for (_row=0; _row<nrows; _row++) { \ 473 _nz = mat->i[_row+1] - mat->i[_row]; \ 474 for (_j=0; _j<_nz; _j++) { \ 475 _col = _j + mat->j + mat->i[_row]; \ 476 PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \ 477 } \ 478 } \ 479 } \ 480 } 481 482 /* 483 Add column indices into table for counting the nonzeros of merged rows 484 */ 485 #define MatMergeRows_SeqAIJ(mat,nrows,rows,ta) { \ 486 PetscInt _j,_row,_nz,*_col,_i; \ 487 for (_i=0; _i<nrows; _i++) {\ 488 _row = rows[_i]; \ 489 _nz = mat->i[_row+1] - mat->i[_row]; \ 490 for (_j=0; _j<_nz; _j++) { \ 491 _col = _j + mat->j + mat->i[_row]; \ 492 PetscTableAdd(ta,*_col+1,1,INSERT_VALUES); \ 493 } \ 494 } \ 495 } 496 497 #endif 498