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