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