xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 55bea0eb881ce1ddd4964f044db496e4529865de)
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 
4953565b12SHong Zhang typedef struct {
5053565b12SHong Zhang   PetscInt    *api,*apj;       /* symbolic structure of A*P */
5153565b12SHong Zhang   PetscScalar *apa;            /* temporary array for storing one row of A*P */
5253565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
5353565b12SHong Zhang } Mat_PtAP;
5453565b12SHong Zhang 
5553565b12SHong Zhang typedef struct {
5653565b12SHong Zhang   MatTransposeColoring matcoloring;
57257c235dSHong Zhang   Mat                  Rt;    /* sparse or dense matrix of R^T */
5853565b12SHong Zhang   Mat                  RARt;  /* dense matrix of R*A*R^T */
593b1b9624SHong Zhang   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
6053565b12SHong Zhang   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
6153565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
6253565b12SHong Zhang } Mat_RARt;
6353565b12SHong Zhang 
646d0b6147SHong Zhang typedef struct {
656d0b6147SHong Zhang   Mat BC;               /* temp matrix for storing B*C */
666d0b6147SHong Zhang   PetscErrorCode (*destroy)(Mat);
676d0b6147SHong Zhang } Mat_MatMatMatMult;
686d0b6147SHong Zhang 
692d40f771SBarry Smith /*
70ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
71e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
72dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
735768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
742d40f771SBarry Smith */
75d35516d3SLois Curfman McInnes 
76e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
77b8a66259SBarry Smith typedef struct {
784108e4d5SBarry Smith   MatScalar *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
79f0d39aaaSBarry Smith   PetscInt  bdiagsize;                         /* length of bdiag and ibdiag */
80ace3abfcSBarry Smith   PetscBool ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
81f0d39aaaSBarry Smith 
82ace3abfcSBarry Smith   PetscBool use;
83e6b907acSBarry Smith   PetscInt  node_count;                     /* number of inodes */
84e6b907acSBarry Smith   PetscInt  *size;                          /* size of each inode */
85e6b907acSBarry Smith   PetscInt  limit;                          /* inode limit */
86e6b907acSBarry Smith   PetscInt  max_limit;                      /* maximum supported inode limit */
87ace3abfcSBarry Smith   PetscBool checked;                        /* if inodes have been checked for */
884108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
89e6b907acSBarry Smith 
905a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
915a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
925a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
935a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
945a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
955a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
965a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
975a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
985a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
99e6b907acSBarry Smith 
100e6b907acSBarry Smith typedef struct {
10154f21887SBarry Smith   SEQAIJHEADER(MatScalar);
1024108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
10354f21887SBarry Smith   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
10471f1c65dSBarry Smith 
10571f1c65dSBarry Smith   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
106ace3abfcSBarry Smith   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
107bbead8a2SBarry Smith   PetscScalar *ibdiag;                        /* inverses of block diagonals */
1082291e4fdSJed Brown   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
10971f1c65dSBarry Smith   PetscScalar fshift,omega;                   /* last used omega and fshift */
11071f1c65dSBarry Smith 
1113a7fca6bSBarry Smith   ISColoring coloring;                        /* set with MatADSetColoring() used by MatADSetValues() */
1120b7e3e3dSHong Zhang 
1130b7e3e3dSHong Zhang   PetscScalar       *matmult_abdense;    /* used by MatMatMult() */
11453565b12SHong Zhang   Mat_PtAP          *ptap;               /* used by MatPtAP() */
1156d0b6147SHong Zhang   Mat_MatMatMatMult *matmatmatmult;      /* used by MatMatMatMult() */
116257c235dSHong Zhang   Mat_RARt          *rart;               /* used by MatRARt() */
117ec8511deSBarry Smith } Mat_SeqAIJ;
118b8a66259SBarry Smith 
119e6b907acSBarry Smith /*
120e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
121e6b907acSBarry Smith */
122e6b907acSBarry Smith #undef __FUNCT__
123e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ"
124d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
125d0f46423SBarry Smith {
126e6b907acSBarry Smith   PetscErrorCode ierr;
127e6b907acSBarry Smith   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
128e6b907acSBarry Smith   if (A->singlemalloc) {
129e6b907acSBarry Smith     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
130e6b907acSBarry Smith   } else {
131c31cb41cSBarry Smith     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
132c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
133c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
134e6b907acSBarry Smith   }
135e6b907acSBarry Smith   return 0;
136e6b907acSBarry Smith }
137e6b907acSBarry Smith /*
138e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
139357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
140e6b907acSBarry Smith */
141fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
142fef13f97SBarry Smith   if (NROW >= RMAX) { \
143fef13f97SBarry Smith     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
144fef13f97SBarry Smith     /* there is no extra room in row, therefore enlarge */ \
145fef13f97SBarry Smith     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
146fef13f97SBarry Smith     datatype *new_a; \
147fef13f97SBarry Smith  \
148fef13f97SBarry Smith     if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
149fef13f97SBarry Smith     /* malloc new storage space */ \
150fef13f97SBarry Smith     ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr); \
151fef13f97SBarry Smith  \
152fef13f97SBarry Smith     /* copy over old data into new slots */ \
153fef13f97SBarry Smith     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
154fef13f97SBarry Smith     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
155fef13f97SBarry Smith     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
156fef13f97SBarry Smith     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
157fef13f97SBarry Smith     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
158fef13f97SBarry Smith     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
159fef13f97SBarry Smith     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
160fef13f97SBarry Smith     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
161fef13f97SBarry Smith     /* free up old matrix storage */ \
162fef13f97SBarry Smith     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
163fef13f97SBarry Smith     AA                = new_a; \
164fef13f97SBarry Smith     Ain->a            = (MatScalar*) new_a;                   \
165fef13f97SBarry Smith     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
166fef13f97SBarry Smith     Ain->singlemalloc = PETSC_TRUE; \
167fef13f97SBarry Smith  \
168fef13f97SBarry Smith     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
169fef13f97SBarry Smith     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
170fef13f97SBarry Smith     Ain->maxnz += BS2*CHUNKSIZE; \
171fef13f97SBarry Smith     Ain->reallocs++; \
172fef13f97SBarry Smith   } \
17317454e89SShri Abhyankar 
174e6b907acSBarry Smith 
1758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
1765a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
1775a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
1785a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
1791df811f5SHong Zhang 
1805a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
1815a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
1825a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
1835a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
1845a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
1855a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
1865a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
1875a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
1885a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
1895a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
1905a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
19108480c60SBarry Smith 
1925a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
1935a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
1945a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
1955a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
1965a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
19708480c60SBarry Smith 
1985a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
1995a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
2005a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
2013a7fca6bSBarry Smith 
2025a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2035a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
2045a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2055a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
2065008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
2075a576424SJed Brown PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
2085a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2095a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2105a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2115a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2125a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
2135a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
2145a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
2155a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
2165a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
2175a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
2185a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
2195a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
2205a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2225a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg);
2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
2325a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
2337bab7c10SHong Zhang 
2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2365a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
2375a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
2385a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
2395a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
2405a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
2422b8ad9a3SHong Zhang 
2435a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2445a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
24853565b12SHong Zhang 
2495a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
250*55bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
251*55bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
2525a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
253*55bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
254*55bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
2555df89d91SHong Zhang 
2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2585a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2605a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2615a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
2635a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
2645a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
2655df89d91SHong Zhang 
2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
2697bab7c10SHong Zhang 
2705a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
2715a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
2725a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
2735a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
2745a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2755a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2785a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
2795a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
2805a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
2819af31e4aSHong Zhang 
2825a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
2845a576424SJed Brown PETSC_INTERN PetscErrorCode Mat_CheckInode(Mat,PetscBool);
2855a576424SJed Brown PETSC_INTERN PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool);
286019b515eSShri Abhyankar 
2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
2889f5f6813SShri Abhyankar 
2898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
2908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
2918cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
2925a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
2935a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2945a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
2982f947c57SVictor Minden 
29970f19b1fSKris Buschelman 
300003131ecSBarry Smith /*
301003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
302003131ecSBarry Smith 
303003131ecSBarry Smith   Input Parameters:
304003131ecSBarry Smith +  nnz - the number of entries
305003131ecSBarry Smith .  r - the array of vector values
306003131ecSBarry Smith .  xv - the matrix values for the row
307003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
308003131ecSBarry Smith 
309003131ecSBarry Smith   Output Parameter:
310003131ecSBarry Smith .  sum - negative the sum of results
311003131ecSBarry Smith 
312003131ecSBarry Smith   PETSc compile flags:
313b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
314e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
315003131ecSBarry Smith 
316003131ecSBarry Smith .seealso: PetscSparseDensePlusDot()
317003131ecSBarry Smith 
318003131ecSBarry Smith */
319519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
320003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
321003131ecSBarry Smith     if (nnz > 0) { \
322003131ecSBarry Smith       switch (nnz & 0x3) { \
323003131ecSBarry Smith       case 3: sum -= *xv++ *r[*xi++]; \
324003131ecSBarry Smith       case 2: sum -= *xv++ *r[*xi++]; \
325003131ecSBarry Smith       case 1: sum -= *xv++ *r[*xi++]; \
326003131ecSBarry Smith         nnz       -= 4;} \
327003131ecSBarry Smith       while (nnz > 0) { \
328003131ecSBarry Smith         sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \
329003131ecSBarry Smith                xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \
330003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
331003131ecSBarry Smith 
332003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
333003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
334003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
335003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
336003131ecSBarry Smith                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
337003131ecSBarry Smith     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
338003131ecSBarry Smith 
339003131ecSBarry Smith #else
340003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
341003131ecSBarry Smith     PetscInt __i; \
342003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
343003131ecSBarry Smith #endif
344003131ecSBarry Smith 
345003131ecSBarry Smith 
346003131ecSBarry Smith 
347003131ecSBarry Smith /*
348003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
349003131ecSBarry Smith 
350003131ecSBarry Smith   Input Parameters:
351003131ecSBarry Smith +  nnz - the number of entries
352003131ecSBarry Smith .  r - the array of vector values
353003131ecSBarry Smith .  xv - the matrix values for the row
354003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
355003131ecSBarry Smith 
356003131ecSBarry Smith   Output Parameter:
357003131ecSBarry Smith .  sum - the sum of results
358003131ecSBarry Smith 
359003131ecSBarry Smith   PETSc compile flags:
360b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
361e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
362003131ecSBarry Smith 
363003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot()
364003131ecSBarry Smith 
365003131ecSBarry Smith */
366519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
367003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
368003131ecSBarry Smith     if (nnz > 0) { \
369003131ecSBarry Smith       switch (nnz & 0x3) { \
370003131ecSBarry Smith       case 3: sum += *xv++ *r[*xi++]; \
371003131ecSBarry Smith       case 2: sum += *xv++ *r[*xi++]; \
372003131ecSBarry Smith       case 1: sum += *xv++ *r[*xi++]; \
373003131ecSBarry Smith         nnz       -= 4;} \
374003131ecSBarry Smith       while (nnz > 0) { \
375003131ecSBarry Smith         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
376003131ecSBarry Smith                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
377003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
378003131ecSBarry Smith 
379003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
380003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
381003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
382003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
383003131ecSBarry Smith                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
384003131ecSBarry Smith     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
385003131ecSBarry Smith 
386003131ecSBarry Smith #else
387003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
388003131ecSBarry Smith     PetscInt __i; \
389003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
390003131ecSBarry Smith #endif
391003131ecSBarry Smith 
3922d40f771SBarry Smith #endif
393