xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 6378f32dd8656d85e3a9cd8ab0fbf0edbabc7e18)
1b8a66259SBarry Smith 
22d40f771SBarry Smith #if !defined(__AIJ_H)
32d40f771SBarry Smith #define __AIJ_H
4e6b907acSBarry Smith 
5b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
68f690400SShri Abhyankar 
7e6b907acSBarry Smith /*
8e6b907acSBarry Smith     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
9e6b907acSBarry Smith */
10421e10b8SBarry Smith #define SEQAIJHEADER(datatype)        \
11ace3abfcSBarry Smith   PetscBool roworiented;              /* if true, row-oriented input, default */ \
12e6b907acSBarry Smith   PetscInt  nonew;                    /* 1 don't add new nonzeros, -1 generate error on new */ \
1328b2fa4aSMatthew Knepley   PetscInt  nounused;                 /* -1 generate error on unused space */ \
14ace3abfcSBarry Smith   PetscBool singlemalloc;             /* if true a, i, and j have been obtained with one big malloc */ \
15e6b907acSBarry Smith   PetscInt  maxnz;                    /* allocated nonzeros */ \
16e6b907acSBarry Smith   PetscInt  *imax;                    /* maximum space allocated for each row */ \
17e6b907acSBarry Smith   PetscInt  *ilen;                    /* actual length of each row */ \
18ace3abfcSBarry Smith   PetscBool free_imax_ilen;  \
19e6b907acSBarry Smith   PetscInt  reallocs;                 /* number of mallocs done during MatSetValues() \
20e6b907acSBarry Smith                                         as more values are set than were prealloced */\
21e6b907acSBarry Smith   PetscInt          rmax;             /* max nonzeros in any row */ \
22ace3abfcSBarry Smith   PetscBool         keepnonzeropattern;   /* keeps matrix structure same in calls to MatZeroRows()*/ \
23ace3abfcSBarry Smith   PetscBool         ignorezeroentries; \
24e6b907acSBarry Smith   PetscInt          *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */ \
25e6b907acSBarry Smith   Mat               XtoY;             /* used by MatAXPY() */ \
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 */                  \
33ace3abfcSBarry Smith   PetscBool         free_diag;         \
34421e10b8SBarry Smith   datatype          *a;               /* nonzero elements */                               \
35e6b907acSBarry Smith   PetscScalar       *solve_work;      /* work space used in MatSolve */                    \
364fd072dbSBarry Smith   IS                row, col, icol;   /* index sets, used for reorderings */ \
37ace3abfcSBarry Smith   PetscBool         pivotinblocks;    /* pivot inside factorization of each diagonal block */ \
384fd072dbSBarry Smith   Mat               parent             /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....);
394fd072dbSBarry Smith                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */
40e6b907acSBarry Smith 
4153565b12SHong Zhang typedef struct {
4253565b12SHong Zhang   MatTransposeColoring matcoloring;
4353565b12SHong Zhang   Mat                  Bt_den;       /* dense matrix of B^T */
4453565b12SHong Zhang   Mat                  ABt_den;      /* dense matrix of A*B^T */
4553565b12SHong Zhang   PetscBool            usecoloring;
4653565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
4753565b12SHong Zhang } Mat_MatMatTransMult;
4853565b12SHong Zhang 
492cff0574SHong Zhang typedef struct { /* for MatTransposeMatMult_SeqAIJ_SeqDense() */
502cff0574SHong Zhang   Mat          mA;           /* maij matrix of A */
512cff0574SHong Zhang   Vec          bt,ct;        /* vectors to hold locally transposed arrays of B and C */
522cff0574SHong Zhang   PetscErrorCode (*destroy)(Mat);
532cff0574SHong Zhang } Mat_MatTransMatMult;
542cff0574SHong Zhang 
5553565b12SHong Zhang typedef struct {
5653565b12SHong Zhang   PetscInt    *api,*apj;       /* symbolic structure of A*P */
5753565b12SHong Zhang   PetscScalar *apa;            /* temporary array for storing one row of A*P */
5853565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
5953565b12SHong Zhang } Mat_PtAP;
6053565b12SHong Zhang 
6153565b12SHong Zhang typedef struct {
6253565b12SHong Zhang   MatTransposeColoring matcoloring;
63257c235dSHong Zhang   Mat                  Rt;    /* sparse or dense matrix of R^T */
6453565b12SHong Zhang   Mat                  RARt;  /* dense matrix of R*A*R^T */
653b1b9624SHong Zhang   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
6653565b12SHong Zhang   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
6753565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
6853565b12SHong Zhang } Mat_RARt;
6953565b12SHong Zhang 
706d0b6147SHong Zhang typedef struct {
716d0b6147SHong Zhang   Mat BC;               /* temp matrix for storing B*C */
726d0b6147SHong Zhang   PetscErrorCode (*destroy)(Mat);
736d0b6147SHong Zhang } Mat_MatMatMatMult;
746d0b6147SHong Zhang 
755cc03489SHong Zhang typedef struct { /* used by MatGetRedundantMatrix() for reusing matredundant */
765cc03489SHong Zhang   PetscInt     nzlocal,nsends,nrecvs;
775cc03489SHong Zhang   PetscMPIInt  *send_rank,*recv_rank;
785cc03489SHong Zhang   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
795cc03489SHong Zhang   PetscScalar  *sbuf_a,**rbuf_a;
800b291e46SHong Zhang   PetscSubcomm psubcomm;
813c79b8e7SHong Zhang   IS           isrow,iscol;
823c79b8e7SHong Zhang   Mat          *matseq;
835cc03489SHong Zhang   PetscErrorCode (*Destroy)(Mat);
845cc03489SHong Zhang } Mat_Redundant;
855cc03489SHong Zhang 
862d40f771SBarry Smith /*
87ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
88e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
89dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
905768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
912d40f771SBarry Smith */
92d35516d3SLois Curfman McInnes 
93e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
94b8a66259SBarry Smith typedef struct {
954108e4d5SBarry Smith   MatScalar *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
96f0d39aaaSBarry Smith   PetscInt  bdiagsize;                         /* length of bdiag and ibdiag */
97ace3abfcSBarry Smith   PetscBool ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
98f0d39aaaSBarry Smith 
99ace3abfcSBarry Smith   PetscBool use;
100e6b907acSBarry Smith   PetscInt  node_count;                     /* number of inodes */
101e6b907acSBarry Smith   PetscInt  *size;                          /* size of each inode */
102e6b907acSBarry Smith   PetscInt  limit;                          /* inode limit */
103e6b907acSBarry Smith   PetscInt  max_limit;                      /* maximum supported inode limit */
104ace3abfcSBarry Smith   PetscBool checked;                        /* if inodes have been checked for */
1054108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
106e6b907acSBarry Smith 
1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
1085a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
1095a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
1115a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
1135a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
116e6b907acSBarry Smith 
117e6b907acSBarry Smith typedef struct {
11854f21887SBarry Smith   SEQAIJHEADER(MatScalar);
1194108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
12054f21887SBarry Smith   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
12171f1c65dSBarry Smith 
12271f1c65dSBarry Smith   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
123ace3abfcSBarry Smith   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
124bbead8a2SBarry Smith   PetscScalar *ibdiag;                        /* inverses of block diagonals */
1252291e4fdSJed Brown   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
12671f1c65dSBarry Smith   PetscScalar fshift,omega;                   /* last used omega and fshift */
12771f1c65dSBarry Smith 
1283a7fca6bSBarry Smith   ISColoring coloring;                        /* set with MatADSetColoring() used by MatADSetValues() */
1290b7e3e3dSHong Zhang 
1300b7e3e3dSHong Zhang   PetscScalar       *matmult_abdense;    /* used by MatMatMult() */
13153565b12SHong Zhang   Mat_PtAP          *ptap;               /* used by MatPtAP() */
1326d0b6147SHong Zhang   Mat_MatMatMatMult *matmatmatmult;      /* used by MatMatMatMult() */
133257c235dSHong Zhang   Mat_RARt          *rart;               /* used by MatRARt() */
13440192850SHong Zhang   Mat_MatMatTransMult *abt;              /* used by MatMatTransposeMult() */
1355cc03489SHong Zhang   Mat_Redundant       *redundant;        /* used by MatGetRedundantMatrix() */
136ec8511deSBarry Smith } Mat_SeqAIJ;
137b8a66259SBarry Smith 
138e6b907acSBarry Smith /*
139e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
140e6b907acSBarry Smith */
141e6b907acSBarry Smith #undef __FUNCT__
142e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ"
143d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
144d0f46423SBarry Smith {
145e6b907acSBarry Smith   PetscErrorCode ierr;
146e6b907acSBarry Smith   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
147e6b907acSBarry Smith   if (A->singlemalloc) {
148e6b907acSBarry Smith     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
149e6b907acSBarry Smith   } else {
150c31cb41cSBarry Smith     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
151c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
152c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
153e6b907acSBarry Smith   }
154e6b907acSBarry Smith   return 0;
155e6b907acSBarry Smith }
156e6b907acSBarry Smith /*
157e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
158357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
159e6b907acSBarry Smith */
160fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
161fef13f97SBarry Smith   if (NROW >= RMAX) { \
162fef13f97SBarry Smith     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
163fef13f97SBarry Smith     /* there is no extra room in row, therefore enlarge */ \
164fef13f97SBarry Smith     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
165fef13f97SBarry Smith     datatype *new_a; \
166fef13f97SBarry Smith  \
167fef13f97SBarry Smith     if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
168fef13f97SBarry Smith     /* malloc new storage space */ \
169fef13f97SBarry Smith     ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr); \
170fef13f97SBarry Smith  \
171fef13f97SBarry Smith     /* copy over old data into new slots */ \
172fef13f97SBarry Smith     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
173fef13f97SBarry Smith     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
174fef13f97SBarry Smith     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
175fef13f97SBarry Smith     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
176fef13f97SBarry Smith     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
177fef13f97SBarry Smith     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
178fef13f97SBarry Smith     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
179fef13f97SBarry Smith     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
180fef13f97SBarry Smith     /* free up old matrix storage */ \
181fef13f97SBarry Smith     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
182fef13f97SBarry Smith     AA                = new_a; \
183fef13f97SBarry Smith     Ain->a            = (MatScalar*) new_a;                   \
184fef13f97SBarry Smith     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
185fef13f97SBarry Smith     Ain->singlemalloc = PETSC_TRUE; \
186fef13f97SBarry Smith  \
187fef13f97SBarry Smith     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
188fef13f97SBarry Smith     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
189fef13f97SBarry Smith     Ain->maxnz += BS2*CHUNKSIZE; \
190fef13f97SBarry Smith     Ain->reallocs++; \
191fef13f97SBarry Smith   } \
19217454e89SShri Abhyankar 
193e6b907acSBarry Smith 
1948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
1955a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
1965a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
1975a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
1981df811f5SHong Zhang 
1995a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
2005a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
2015a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
2025a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
2035a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2045a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2055a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
2065a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
2075a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
2085a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2095a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
21008480c60SBarry Smith 
2115a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
2125a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2135a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
2145a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2155a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
21608480c60SBarry Smith 
2175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
2185a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
2195a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
2203a7fca6bSBarry Smith 
2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2225a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
2255008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
2375a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
2385a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
2395a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
2405a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2425a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
2435a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
2445a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
24879c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
2495a576424SJed Brown PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
25079c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringApply_SeqAIJ(Mat,MatFDColoring,Vec,MatStructure*,void*);
2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
2525a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
2537bab7c10SHong Zhang 
2545a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2555a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
2622b8ad9a3SHong Zhang 
2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2644a1b09b7SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat,Mat,PetscReal,Mat*);
2655a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
26853565b12SHong Zhang 
2695a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
27055bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
27155bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
2725a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
27355bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
27455bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
2755df89d91SHong Zhang 
2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2785a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2793bf78175SHong Zhang 
2803bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
2813bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
2823bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
2833bf78175SHong Zhang 
2845a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2855a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
2885a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
2895a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
2905df89d91SHong Zhang 
2915a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
2925a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
2935a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
2947bab7c10SHong Zhang 
2955a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
2985a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
2995a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3005a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3015a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
3025a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
303*6378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
304*6378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
3055a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
3065a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
3075a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
3089af31e4aSHong Zhang 
3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
3115a576424SJed Brown PETSC_INTERN PetscErrorCode Mat_CheckInode(Mat,PetscBool);
3125a576424SJed Brown PETSC_INTERN PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool);
313019b515eSShri Abhyankar 
3145a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
3159f5f6813SShri Abhyankar 
3168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
3178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
3188cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
3195a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
3205a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3215a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3228cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
3235a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
3245a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
3252f947c57SVictor Minden 
32670f19b1fSKris Buschelman 
327003131ecSBarry Smith /*
328003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
329003131ecSBarry Smith 
330003131ecSBarry Smith   Input Parameters:
331003131ecSBarry Smith +  nnz - the number of entries
332003131ecSBarry Smith .  r - the array of vector values
333003131ecSBarry Smith .  xv - the matrix values for the row
334003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
335003131ecSBarry Smith 
336003131ecSBarry Smith   Output Parameter:
337003131ecSBarry Smith .  sum - negative the sum of results
338003131ecSBarry Smith 
339003131ecSBarry Smith   PETSc compile flags:
340b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
341e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
342003131ecSBarry Smith 
343003131ecSBarry Smith .seealso: PetscSparseDensePlusDot()
344003131ecSBarry Smith 
345003131ecSBarry Smith */
346519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
347003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
348003131ecSBarry Smith     if (nnz > 0) { \
349003131ecSBarry Smith       switch (nnz & 0x3) { \
350003131ecSBarry Smith       case 3: sum -= *xv++ *r[*xi++]; \
351003131ecSBarry Smith       case 2: sum -= *xv++ *r[*xi++]; \
352003131ecSBarry Smith       case 1: sum -= *xv++ *r[*xi++]; \
353003131ecSBarry Smith         nnz       -= 4;} \
354003131ecSBarry Smith       while (nnz > 0) { \
355003131ecSBarry Smith         sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \
356003131ecSBarry Smith                xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \
357003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
358003131ecSBarry Smith 
359003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
360003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
361003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
362003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
363003131ecSBarry Smith                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
364003131ecSBarry Smith     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
365003131ecSBarry Smith 
366003131ecSBarry Smith #else
367003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
368003131ecSBarry Smith     PetscInt __i; \
369003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
370003131ecSBarry Smith #endif
371003131ecSBarry Smith 
372003131ecSBarry Smith 
373003131ecSBarry Smith 
374003131ecSBarry Smith /*
375003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
376003131ecSBarry Smith 
377003131ecSBarry Smith   Input Parameters:
378003131ecSBarry Smith +  nnz - the number of entries
379003131ecSBarry Smith .  r - the array of vector values
380003131ecSBarry Smith .  xv - the matrix values for the row
381003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
382003131ecSBarry Smith 
383003131ecSBarry Smith   Output Parameter:
384003131ecSBarry Smith .  sum - the sum of results
385003131ecSBarry Smith 
386003131ecSBarry Smith   PETSc compile flags:
387b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
388e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
389003131ecSBarry Smith 
390003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot()
391003131ecSBarry Smith 
392003131ecSBarry Smith */
393519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
394003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
395003131ecSBarry Smith     if (nnz > 0) { \
396003131ecSBarry Smith       switch (nnz & 0x3) { \
397003131ecSBarry Smith       case 3: sum += *xv++ *r[*xi++]; \
398003131ecSBarry Smith       case 2: sum += *xv++ *r[*xi++]; \
399003131ecSBarry Smith       case 1: sum += *xv++ *r[*xi++]; \
400003131ecSBarry Smith         nnz       -= 4;} \
401003131ecSBarry Smith       while (nnz > 0) { \
402003131ecSBarry Smith         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
403003131ecSBarry Smith                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
404003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
405003131ecSBarry Smith 
406003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
407003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
408003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
409003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
410003131ecSBarry Smith                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
411003131ecSBarry Smith     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
412003131ecSBarry Smith 
413003131ecSBarry Smith #else
414003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
415003131ecSBarry Smith     PetscInt __i; \
416003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
417003131ecSBarry Smith #endif
418003131ecSBarry Smith 
419b434eb95SMatthew G. Knepley 
420b434eb95SMatthew G. Knepley /*
421b434eb95SMatthew G. Knepley     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
422b434eb95SMatthew G. Knepley 
423b434eb95SMatthew G. Knepley   Input Parameters:
424b434eb95SMatthew G. Knepley +  nnz - the number of entries
425b434eb95SMatthew G. Knepley .  r - the array of vector values
426b434eb95SMatthew G. Knepley .  xv - the matrix values for the row
427b434eb95SMatthew G. Knepley -  xi - the column indices of the nonzeros in the row
428b434eb95SMatthew G. Knepley 
429b434eb95SMatthew G. Knepley   Output Parameter:
430b434eb95SMatthew G. Knepley .  max - the max of results
431b434eb95SMatthew G. Knepley 
432b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
433b434eb95SMatthew G. Knepley 
434b434eb95SMatthew G. Knepley */
435b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
436b434eb95SMatthew G. Knepley     PetscInt __i; \
4375e792707SMatthew G. Knepley     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
438b434eb95SMatthew G. Knepley 
4392d40f771SBarry Smith #endif
440