xref: /petsc/src/mat/impls/aij/seq/aij.h (revision a02bda8e9177763e3a7baa8b1d4a3ca70a8d6433)
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; \
24ace3abfcSBarry Smith   PetscBool         free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
25ace3abfcSBarry Smith   PetscBool         free_a;           /* free the numerical values when matrix is destroy */ \
26e6b907acSBarry Smith   Mat_CompressedRow compressedrow;    /* use compressed row format */                      \
27e6b907acSBarry Smith   PetscInt          nz;               /* nonzeros */                                       \
28e6b907acSBarry Smith   PetscInt          *i;               /* pointer to beginning of each row */               \
29e6b907acSBarry Smith   PetscInt          *j;               /* column values: j + i[k] - 1 is start of row k */  \
30e6b907acSBarry Smith   PetscInt          *diag;            /* pointers to diagonal elements */                  \
317b083b7cSBarry Smith   PetscInt          nonzerorowcnt;    /* how many rows have nonzero entries */             \
32ace3abfcSBarry Smith   PetscBool         free_diag;         \
33421e10b8SBarry Smith   datatype          *a;               /* nonzero elements */                               \
34e6b907acSBarry Smith   PetscScalar       *solve_work;      /* work space used in MatSolve */                    \
354fd072dbSBarry Smith   IS                row, col, icol;   /* index sets, used for reorderings */ \
36ace3abfcSBarry Smith   PetscBool         pivotinblocks;    /* pivot inside factorization of each diagonal block */ \
374fd072dbSBarry Smith   Mat               parent             /* set if this matrix was formed with MatDuplicate(...,MAT_SHARE_NONZERO_PATTERN,....);
384fd072dbSBarry Smith                                          means that this shares some data structures with the parent including diag, ilen, imax, i, j */
39e6b907acSBarry Smith 
4053565b12SHong Zhang typedef struct {
4153565b12SHong Zhang   MatTransposeColoring matcoloring;
4253565b12SHong Zhang   Mat                  Bt_den;       /* dense matrix of B^T */
4353565b12SHong Zhang   Mat                  ABt_den;      /* dense matrix of A*B^T */
4453565b12SHong Zhang   PetscBool            usecoloring;
4553565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
4653565b12SHong Zhang } Mat_MatMatTransMult;
4753565b12SHong Zhang 
482cff0574SHong Zhang typedef struct { /* for MatTransposeMatMult_SeqAIJ_SeqDense() */
492cff0574SHong Zhang   Mat          mA;           /* maij matrix of A */
502cff0574SHong Zhang   Vec          bt,ct;        /* vectors to hold locally transposed arrays of B and C */
512cff0574SHong Zhang   PetscErrorCode (*destroy)(Mat);
522cff0574SHong Zhang } Mat_MatTransMatMult;
532cff0574SHong Zhang 
5453565b12SHong Zhang typedef struct {
5553565b12SHong Zhang   PetscInt    *api,*apj;       /* symbolic structure of A*P */
5653565b12SHong Zhang   PetscScalar *apa;            /* temporary array for storing one row of A*P */
5753565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
5853565b12SHong Zhang } Mat_PtAP;
5953565b12SHong Zhang 
6053565b12SHong Zhang typedef struct {
6153565b12SHong Zhang   MatTransposeColoring matcoloring;
62257c235dSHong Zhang   Mat                  Rt;    /* sparse or dense matrix of R^T */
6353565b12SHong Zhang   Mat                  RARt;  /* dense matrix of R*A*R^T */
643b1b9624SHong Zhang   Mat                  ARt;   /* A*R^T used for the case -matrart_color_art */
6553565b12SHong Zhang   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
6653565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
6753565b12SHong Zhang } Mat_RARt;
6853565b12SHong Zhang 
696d0b6147SHong Zhang typedef struct {
706d0b6147SHong Zhang   Mat BC;               /* temp matrix for storing B*C */
716d0b6147SHong Zhang   PetscErrorCode (*destroy)(Mat);
726d0b6147SHong Zhang } Mat_MatMatMatMult;
736d0b6147SHong Zhang 
742d40f771SBarry Smith /*
75ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
76e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
77dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
785768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
792d40f771SBarry Smith */
80d35516d3SLois Curfman McInnes 
81e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
82b8a66259SBarry Smith typedef struct {
834108e4d5SBarry Smith   MatScalar        *bdiag,*ibdiag,*ssor_work;        /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
84f0d39aaaSBarry Smith   PetscInt         bdiagsize;                         /* length of bdiag and ibdiag */
85ace3abfcSBarry Smith   PetscBool        ibdiagvalid;                       /* do ibdiag[] and bdiag[] contain the most recent values */
86f0d39aaaSBarry Smith 
87ace3abfcSBarry Smith   PetscBool        use;
88e6b907acSBarry Smith   PetscInt         node_count;                     /* number of inodes */
89e6b907acSBarry Smith   PetscInt         *size;                          /* size of each inode */
90e6b907acSBarry Smith   PetscInt         limit;                          /* inode limit */
91e6b907acSBarry Smith   PetscInt         max_limit;                      /* maximum supported inode limit */
92ace3abfcSBarry Smith   PetscBool        checked;                        /* if inodes have been checked for */
93*a02bda8eSBarry Smith   PetscObjectState mat_nonzerostate;               /* non-zero state when inodes were checked for */
944108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
95e6b907acSBarry Smith 
965a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
975a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
985a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
995a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
1005a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool);
1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
1025a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
105e6b907acSBarry Smith 
106e6b907acSBarry Smith typedef struct {
10754f21887SBarry Smith   SEQAIJHEADER(MatScalar);
1084108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
10954f21887SBarry Smith   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
11071f1c65dSBarry Smith 
11171f1c65dSBarry Smith   PetscScalar *idiag,*mdiag,*ssor_work;       /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
112ace3abfcSBarry Smith   PetscBool   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
113bbead8a2SBarry Smith   PetscScalar *ibdiag;                        /* inverses of block diagonals */
1142291e4fdSJed Brown   PetscBool   ibdiagvalid;                    /* inverses of block diagonals are valid. */
11571f1c65dSBarry Smith   PetscScalar fshift,omega;                   /* last used omega and fshift */
11671f1c65dSBarry Smith 
1173a7fca6bSBarry Smith   ISColoring coloring;                        /* set with MatADSetColoring() used by MatADSetValues() */
1180b7e3e3dSHong Zhang 
1190b7e3e3dSHong Zhang   PetscScalar       *matmult_abdense;    /* used by MatMatMult() */
12053565b12SHong Zhang   Mat_PtAP          *ptap;               /* used by MatPtAP() */
1216d0b6147SHong Zhang   Mat_MatMatMatMult *matmatmatmult;      /* used by MatMatMatMult() */
122257c235dSHong Zhang   Mat_RARt          *rart;               /* used by MatRARt() */
12340192850SHong Zhang   Mat_MatMatTransMult *abt;              /* used by MatMatTransposeMult() */
124ec8511deSBarry Smith } Mat_SeqAIJ;
125b8a66259SBarry Smith 
126e6b907acSBarry Smith /*
127e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
128e6b907acSBarry Smith */
129e6b907acSBarry Smith #undef __FUNCT__
130e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ"
131d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
132d0f46423SBarry Smith {
133e6b907acSBarry Smith   PetscErrorCode ierr;
134e6b907acSBarry Smith   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
135e6b907acSBarry Smith   if (A->singlemalloc) {
136e6b907acSBarry Smith     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
137e6b907acSBarry Smith   } else {
138c31cb41cSBarry Smith     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
139c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
140c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
141e6b907acSBarry Smith   }
142e6b907acSBarry Smith   return 0;
143e6b907acSBarry Smith }
144e6b907acSBarry Smith /*
145e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
146357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
147e6b907acSBarry Smith */
148fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
149fef13f97SBarry Smith   if (NROW >= RMAX) { \
150fef13f97SBarry Smith     Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data; \
151fef13f97SBarry Smith     /* there is no extra room in row, therefore enlarge */ \
152fef13f97SBarry Smith     PetscInt CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
153fef13f97SBarry Smith     datatype *new_a; \
154fef13f97SBarry Smith  \
1555aae435aSMatthew G. Knepley     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); \
156fef13f97SBarry Smith     /* malloc new storage space */ \
157dcca6d9dSJed Brown     ierr = PetscMalloc3(BS2*new_nz,&new_a,new_nz,&new_j,AM+1,&new_i);CHKERRQ(ierr); \
158fef13f97SBarry Smith  \
159fef13f97SBarry Smith     /* copy over old data into new slots */ \
160fef13f97SBarry Smith     for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
161fef13f97SBarry Smith     for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
162fef13f97SBarry Smith     ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
163fef13f97SBarry Smith     len  = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
164fef13f97SBarry Smith     ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
165fef13f97SBarry Smith     ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
166fef13f97SBarry Smith     ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr); \
167fef13f97SBarry Smith     ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
168fef13f97SBarry Smith     /* free up old matrix storage */ \
169fef13f97SBarry Smith     ierr              = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr); \
170fef13f97SBarry Smith     AA                = new_a; \
171fef13f97SBarry Smith     Ain->a            = (MatScalar*) new_a;                   \
172fef13f97SBarry Smith     AI                = Ain->i = new_i; AJ = Ain->j = new_j;  \
173fef13f97SBarry Smith     Ain->singlemalloc = PETSC_TRUE; \
174fef13f97SBarry Smith  \
175fef13f97SBarry Smith     RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
176fef13f97SBarry Smith     RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
177fef13f97SBarry Smith     Ain->maxnz += BS2*CHUNKSIZE; \
178fef13f97SBarry Smith     Ain->reallocs++; \
179fef13f97SBarry Smith   } \
18017454e89SShri Abhyankar 
181e6b907acSBarry Smith 
1828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
1835a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
1845a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
1855a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
1861df811f5SHong Zhang 
1875a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
1885a576424SJed Brown PETSC_INTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
1895a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
1905a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
1915a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
1925a576424SJed Brown PETSC_INTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
1935a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
1945a576424SJed Brown PETSC_INTERN PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
1955a576424SJed Brown PETSC_INTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool*,PetscInt*);
1965a576424SJed Brown PETSC_INTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
1975a576424SJed Brown PETSC_INTERN PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat,PetscInt*,PetscInt**);
19808480c60SBarry Smith 
1995a576424SJed Brown PETSC_INTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
2005a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2015a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
2025a576424SJed Brown PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
2035a576424SJed Brown PETSC_INTERN PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
20408480c60SBarry Smith 
2055a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_SeqAIJ(Mat,MatOption,PetscBool);
2065a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
2075a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
2083a7fca6bSBarry Smith 
2095a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2105a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
2115a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
2125a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
2135008f5a7SHong Zhang PETSC_INTERN PetscErrorCode MatTranspose_SeqAIJ(Mat,MatReuse,Mat*);
2145a576424SJed Brown PETSC_INTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
2155a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
2165a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
2175a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
2185a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
2195a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
2205a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
2215a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
2225a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
2235a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
2245a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
2255a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
2265a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
2275a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
2285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2295a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2305a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
2315a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
2325a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
2335a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
2345a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
2355a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
23679c1e64dSHong Zhang PETSC_INTERN PetscErrorCode MatEqual_SeqAIJ(Mat,Mat,PetscBool*);
23793dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat,ISColoring,MatFDColoring);
238f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat,ISColoring,MatFDColoring);
239a8971b87SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat,MatFDColoring,PetscInt);
240d1e9a80fSBarry Smith PETSC_INTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*);
2415a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
2425a576424SJed Brown PETSC_INTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
2437bab7c10SHong Zhang 
2445a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2455a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
2475a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
2485a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat,Mat,PetscReal,Mat*);
2495a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat,Mat,PetscReal,Mat*);
2505a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2515a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
2522b8ad9a3SHong Zhang 
2535a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2544a1b09b7SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat,Mat,PetscReal,Mat*);
2555a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,PetscReal,Mat*);
2565a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2575a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
25853565b12SHong Zhang 
2595a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
26055bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,PetscReal,Mat*);
26155bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,PetscReal,Mat*);
2625a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
26355bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat,Mat,Mat);
26455bea0ebSHong Zhang PETSC_INTERN PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat,Mat,Mat);
2655df89d91SHong Zhang 
2665a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2675a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2685a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2693bf78175SHong Zhang 
2703bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
2713bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
2723bf78175SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
2733bf78175SHong Zhang 
2745a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2755a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2765a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2775a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
2785a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
2795a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
2805df89d91SHong Zhang 
2815a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
2825a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,PetscReal,Mat*);
2835a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(Mat,Mat,Mat,Mat);
2847bab7c10SHong Zhang 
2855a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
2865a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
2875a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
2885a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
2895a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2905a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2915a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2925a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool*);
2936378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
2946378f32dSHong Zhang PETSC_INTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat,PetscInt,PetscBool,PetscBool,PetscInt*,const PetscInt *[],const PetscInt *[],PetscInt *[],PetscBool*);
2955a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
2965a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUp_SeqAIJ(Mat);
2975a576424SJed Brown PETSC_INTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
2989af31e4aSHong Zhang 
2995a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat);
3005a576424SJed Brown PETSC_INTERN PetscErrorCode MatSeqAIJInvalidateDiagonal_Inode(Mat);
301*a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode(Mat);
302*a02bda8eSBarry Smith PETSC_INTERN PetscErrorCode MatSeqAIJCheckInode_FactorLU(Mat);
303019b515eSShri Abhyankar 
3045a576424SJed Brown PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
3059f5f6813SShri Abhyankar 
3068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat,MatType,MatReuse,Mat*);
3078cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat,MatType,MatReuse,Mat*);
3088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat,MatType,MatReuse,Mat*);
3095a576424SJed Brown PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
3105a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3115a576424SJed Brown PETSC_INTERN PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
3128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat);
3135a576424SJed Brown PETSC_INTERN PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
3145a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
3152f947c57SVictor Minden 
316b264fe52SHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
3179c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
3189c8f2541SHong Zhang PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
31970f19b1fSKris Buschelman 
320003131ecSBarry Smith /*
321003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
322003131ecSBarry Smith 
323003131ecSBarry Smith   Input Parameters:
324003131ecSBarry Smith +  nnz - the number of entries
325003131ecSBarry Smith .  r - the array of vector values
326003131ecSBarry Smith .  xv - the matrix values for the row
327003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
328003131ecSBarry Smith 
329003131ecSBarry Smith   Output Parameter:
330003131ecSBarry Smith .  sum - negative the sum of results
331003131ecSBarry Smith 
332003131ecSBarry Smith   PETSc compile flags:
333b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
334e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
335003131ecSBarry Smith 
336003131ecSBarry Smith .seealso: PetscSparseDensePlusDot()
337003131ecSBarry Smith 
338003131ecSBarry Smith */
339519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
340003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
341003131ecSBarry Smith     if (nnz > 0) { \
342003131ecSBarry Smith       switch (nnz & 0x3) { \
343003131ecSBarry Smith       case 3: sum -= *xv++ *r[*xi++]; \
344003131ecSBarry Smith       case 2: sum -= *xv++ *r[*xi++]; \
345003131ecSBarry Smith       case 1: sum -= *xv++ *r[*xi++]; \
346003131ecSBarry Smith         nnz       -= 4;} \
347003131ecSBarry Smith       while (nnz > 0) { \
348003131ecSBarry Smith         sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] - \
349003131ecSBarry Smith                xv[2] * r[xi[2]] - xv[3] * r[xi[3]]; \
350003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
351003131ecSBarry Smith 
352003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
353003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
354003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
355003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
356003131ecSBarry Smith                                     sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
357003131ecSBarry Smith     if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
358003131ecSBarry Smith 
359003131ecSBarry Smith #else
360003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) { \
361003131ecSBarry Smith     PetscInt __i; \
362003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum -= xv[__i] * r[xi[__i]];}
363003131ecSBarry Smith #endif
364003131ecSBarry Smith 
365003131ecSBarry Smith 
366003131ecSBarry Smith 
367003131ecSBarry Smith /*
368003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
369003131ecSBarry Smith 
370003131ecSBarry Smith   Input Parameters:
371003131ecSBarry Smith +  nnz - the number of entries
372003131ecSBarry Smith .  r - the array of vector values
373003131ecSBarry Smith .  xv - the matrix values for the row
374003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
375003131ecSBarry Smith 
376003131ecSBarry Smith   Output Parameter:
377003131ecSBarry Smith .  sum - the sum of results
378003131ecSBarry Smith 
379003131ecSBarry Smith   PETSc compile flags:
380b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
381e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
382003131ecSBarry Smith 
383003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot()
384003131ecSBarry Smith 
385003131ecSBarry Smith */
386519f805aSKarl Rupp #if defined(PETSC_KERNEL_USE_UNROLL_4)
387003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
388003131ecSBarry Smith     if (nnz > 0) { \
389003131ecSBarry Smith       switch (nnz & 0x3) { \
390003131ecSBarry Smith       case 3: sum += *xv++ *r[*xi++]; \
391003131ecSBarry Smith       case 2: sum += *xv++ *r[*xi++]; \
392003131ecSBarry Smith       case 1: sum += *xv++ *r[*xi++]; \
393003131ecSBarry Smith         nnz       -= 4;} \
394003131ecSBarry Smith       while (nnz > 0) { \
395003131ecSBarry Smith         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
396003131ecSBarry Smith                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
397003131ecSBarry Smith         xv += 4; xi += 4; nnz -= 4; }}}
398003131ecSBarry Smith 
399003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
400003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
401003131ecSBarry Smith     PetscInt __i,__i1,__i2; \
402003131ecSBarry Smith     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
403003131ecSBarry Smith                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
404003131ecSBarry Smith     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
405003131ecSBarry Smith 
406003131ecSBarry Smith #else
407003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) { \
408003131ecSBarry Smith     PetscInt __i; \
409003131ecSBarry Smith     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
410003131ecSBarry Smith #endif
411003131ecSBarry Smith 
412b434eb95SMatthew G. Knepley 
413b434eb95SMatthew G. Knepley /*
414b434eb95SMatthew G. Knepley     PetscSparseDenseMaxDot - The inner kernel of a modified matrix-vector product \max_i xv[i] * r[xi[i]] for CSR storage
415b434eb95SMatthew G. Knepley 
416b434eb95SMatthew G. Knepley   Input Parameters:
417b434eb95SMatthew G. Knepley +  nnz - the number of entries
418b434eb95SMatthew G. Knepley .  r - the array of vector values
419b434eb95SMatthew G. Knepley .  xv - the matrix values for the row
420b434eb95SMatthew G. Knepley -  xi - the column indices of the nonzeros in the row
421b434eb95SMatthew G. Knepley 
422b434eb95SMatthew G. Knepley   Output Parameter:
423b434eb95SMatthew G. Knepley .  max - the max of results
424b434eb95SMatthew G. Knepley 
425b434eb95SMatthew G. Knepley .seealso: PetscSparseDensePlusDot(), PetscSparseDenseMinusDot()
426b434eb95SMatthew G. Knepley 
427b434eb95SMatthew G. Knepley */
428b434eb95SMatthew G. Knepley #define PetscSparseDenseMaxDot(max,r,xv,xi,nnz) { \
429b434eb95SMatthew G. Knepley     PetscInt __i; \
4305e792707SMatthew G. Knepley     for (__i=0; __i<nnz; __i++) max = PetscMax(PetscRealPart(max), PetscRealPart(xv[__i] * r[xi[__i]]));}
431b434eb95SMatthew G. Knepley 
4322d40f771SBarry Smith #endif
433