xref: /petsc/src/mat/impls/aij/seq/aij.h (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1b8a66259SBarry Smith 
22d40f771SBarry Smith #if !defined(__AIJ_H)
32d40f771SBarry Smith #define __AIJ_H
4e6b907acSBarry Smith 
5*b45d2f2cSJed 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;
5753565b12SHong Zhang   Mat                  Rt;    /* dense matrix of R^T */
5853565b12SHong Zhang   Mat                  RARt;  /* dense matrix of R*A*R^T */
5953565b12SHong Zhang   MatScalar            *work; /* work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
6053565b12SHong Zhang   PetscErrorCode (*destroy)(Mat);
6153565b12SHong Zhang } Mat_RARt;
6253565b12SHong Zhang 
632d40f771SBarry Smith /*
64ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
65e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
66dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
675768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
682d40f771SBarry Smith */
69d35516d3SLois Curfman McInnes 
70e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
71b8a66259SBarry Smith typedef struct {
724108e4d5SBarry Smith   MatScalar   *bdiag,*ibdiag,*ssor_work;      /* diagonal blocks of matrix used for MatSOR_SeqAIJ_Inode() */
73f0d39aaaSBarry Smith   PetscInt    bdiagsize;                       /* length of bdiag and ibdiag */
74ace3abfcSBarry Smith   PetscBool   ibdiagvalid;                     /* do ibdiag[] and bdiag[] contain the most recent values */
75f0d39aaaSBarry Smith 
76ace3abfcSBarry Smith   PetscBool  use;
77e6b907acSBarry Smith   PetscInt   node_count;                    /* number of inodes */
78e6b907acSBarry Smith   PetscInt   *size;                         /* size of each inode */
79e6b907acSBarry Smith   PetscInt   limit;                         /* inode limit */
80e6b907acSBarry Smith   PetscInt   max_limit;                     /* maximum supported inode limit */
81ace3abfcSBarry Smith   PetscBool  checked;                       /* if inodes have been checked for */
824108e4d5SBarry Smith } Mat_SeqAIJ_Inode;
83e6b907acSBarry Smith 
8409573ac7SBarry Smith extern PetscErrorCode MatView_SeqAIJ_Inode(Mat,PetscViewer);
8509573ac7SBarry Smith extern PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat,MatAssemblyType);
8609573ac7SBarry Smith extern PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat);
8709573ac7SBarry Smith extern PetscErrorCode MatCreate_SeqAIJ_Inode(Mat);
8809573ac7SBarry Smith extern PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat,MatOption,PetscBool );
8909573ac7SBarry Smith extern PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat,MatDuplicateOption,Mat*);
9009573ac7SBarry Smith extern PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscBool );
9109573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode_inplace(Mat,Mat,const MatFactorInfo*);
9209573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_Inode(Mat,Mat,const MatFactorInfo*);
93e6b907acSBarry Smith 
94e6b907acSBarry Smith typedef struct {
9554f21887SBarry Smith   SEQAIJHEADER(MatScalar);
964108e4d5SBarry Smith   Mat_SeqAIJ_Inode inode;
9754f21887SBarry Smith   MatScalar        *saved_values;             /* location for stashing nonzero values of matrix */
9871f1c65dSBarry Smith 
9971f1c65dSBarry Smith   PetscScalar      *idiag,*mdiag,*ssor_work;  /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
100ace3abfcSBarry Smith   PetscBool        idiagvalid;                     /* current idiag[] and mdiag[] are valid */
101bbead8a2SBarry Smith   PetscScalar      *ibdiag;                   /* inverses of block diagonals */
1022291e4fdSJed Brown   PetscBool        ibdiagvalid;               /* inverses of block diagonals are valid. */
10371f1c65dSBarry Smith   PetscScalar      fshift,omega;                   /* last used omega and fshift */
10471f1c65dSBarry Smith 
1053a7fca6bSBarry Smith   ISColoring       coloring;                  /* set with MatADSetColoring() used by MatADSetValues() */
1060b7e3e3dSHong Zhang 
1070b7e3e3dSHong Zhang   PetscScalar      *matmult_abdense;     /* used by MatMatMult() */
10853565b12SHong Zhang   Mat_PtAP         *ptap;                /* used by MatPtAP() */
109ec8511deSBarry Smith } Mat_SeqAIJ;
110b8a66259SBarry Smith 
111e6b907acSBarry Smith /*
112e6b907acSBarry Smith   Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
113e6b907acSBarry Smith */
114e6b907acSBarry Smith #undef __FUNCT__
115e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ"
116d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
117d0f46423SBarry Smith {
118e6b907acSBarry Smith   PetscErrorCode ierr;
119e6b907acSBarry Smith   Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
120e6b907acSBarry Smith   if (A->singlemalloc) {
121e6b907acSBarry Smith     ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
122e6b907acSBarry Smith   } else {
123c31cb41cSBarry Smith     if (A->free_a)  {ierr = PetscFree(*a);CHKERRQ(ierr);}
124c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*j);CHKERRQ(ierr);}
125c31cb41cSBarry Smith     if (A->free_ij) {ierr = PetscFree(*i);CHKERRQ(ierr);}
126e6b907acSBarry Smith   }
127e6b907acSBarry Smith   return 0;
128e6b907acSBarry Smith }
129e6b907acSBarry Smith /*
130e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
131357fed3dSBarry Smith     This is a macro because it takes the datatype as an argument which can be either a Mat or a MatScalar
132e6b907acSBarry Smith */
133fef13f97SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
134fef13f97SBarry Smith   if (NROW >= RMAX) {\
135fef13f97SBarry Smith 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
136fef13f97SBarry Smith         /* there is no extra room in row, therefore enlarge */ \
137fef13f97SBarry Smith         PetscInt   CHUNKSIZE = 15,new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
138fef13f97SBarry Smith         datatype   *new_a; \
139fef13f97SBarry Smith  \
140fef13f97SBarry Smith         if (NONEW == -2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
141fef13f97SBarry Smith         /* malloc new storage space */ \
142fef13f97SBarry Smith         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
143fef13f97SBarry Smith  \
144fef13f97SBarry Smith         /* copy over old data into new slots */ \
145fef13f97SBarry Smith         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
146fef13f97SBarry Smith         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
147fef13f97SBarry Smith         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
148fef13f97SBarry Smith         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
149fef13f97SBarry Smith         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
150fef13f97SBarry Smith         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
151fef13f97SBarry Smith         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
152fef13f97SBarry Smith         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
153fef13f97SBarry Smith         /* free up old matrix storage */ \
154fef13f97SBarry Smith         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
155fef13f97SBarry Smith         AA = new_a; \
156fef13f97SBarry Smith         Ain->a = (MatScalar*) new_a;		   \
157fef13f97SBarry Smith         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
158fef13f97SBarry Smith         Ain->singlemalloc = PETSC_TRUE; \
159fef13f97SBarry Smith  \
160fef13f97SBarry Smith         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
161fef13f97SBarry Smith         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
162fef13f97SBarry Smith         Ain->maxnz += BS2*CHUNKSIZE; \
163fef13f97SBarry Smith         Ain->reallocs++; \
164fef13f97SBarry Smith       } \
16517454e89SShri Abhyankar 
166e6b907acSBarry Smith 
167f349cc81SSatish Balay EXTERN_C_BEGIN
16809573ac7SBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,const PetscInt*);
169f349cc81SSatish Balay EXTERN_C_END
17009573ac7SBarry Smith extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
17109573ac7SBarry Smith extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
17209573ac7SBarry Smith extern PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat,Mat,IS,IS,const MatFactorInfo*);
1731df811f5SHong Zhang 
17409573ac7SBarry Smith extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
17509573ac7SBarry Smith extern PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
17609573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,const MatFactorInfo*);
17709573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
17809573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
17909573ac7SBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
18009573ac7SBarry Smith extern PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
18109573ac7SBarry Smith extern PetscErrorCode MatCopy_SeqAIJ(Mat,Mat,MatStructure);
18209573ac7SBarry Smith extern PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscBool *,PetscInt*);
18309573ac7SBarry Smith extern PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
18408480c60SBarry Smith 
18509573ac7SBarry Smith extern PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
18609573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
18709573ac7SBarry Smith extern PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
18809573ac7SBarry Smith extern PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
18909573ac7SBarry Smith extern PetscErrorCode MatSOR_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
19008480c60SBarry Smith 
19109573ac7SBarry Smith extern PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
19209573ac7SBarry Smith extern PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
19309573ac7SBarry Smith extern PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
1943a7fca6bSBarry Smith 
19509573ac7SBarry Smith extern PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
19609573ac7SBarry Smith extern PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
19709573ac7SBarry Smith extern PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
1984e938277SHong Zhang extern PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat,Mat*);
19909573ac7SBarry Smith extern PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
20009573ac7SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat,Mat,IS,IS,const MatFactorInfo*);
20109573ac7SBarry Smith extern PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
20209573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat,Mat,const MatFactorInfo*);
20309573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
20409573ac7SBarry Smith extern PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
20509573ac7SBarry Smith extern PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
20609573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_inplace(Mat,Vec,Vec);
20709573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
20809573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_Inode_inplace(Mat,Vec,Vec);
20909573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_Inode(Mat,Vec,Vec);
21009573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat,Vec,Vec);
21109573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
21209573ac7SBarry Smith extern PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
21309573ac7SBarry Smith extern PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
21409573ac7SBarry Smith extern PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
21509573ac7SBarry Smith extern PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat,Vec,Vec);
21609573ac7SBarry Smith extern PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
21709573ac7SBarry Smith extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat,Vec,Vec,Vec);
21809573ac7SBarry Smith extern PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
21909573ac7SBarry Smith extern PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat,Mat,Mat);
22009573ac7SBarry Smith extern PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
22109573ac7SBarry Smith extern PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg);
22209573ac7SBarry Smith extern PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
22309573ac7SBarry Smith extern PetscErrorCode MatLoad_SeqAIJ(Mat,PetscViewer);
22409573ac7SBarry Smith extern PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
22509573ac7SBarry Smith extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
22625023028SHong Zhang extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,PetscReal,Mat*);
2273c50cad2SHong Zhang extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat,Mat,PetscReal,Mat*);
22809573ac7SBarry Smith extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
22925023028SHong Zhang extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat,Mat,Mat);
2302b8ad9a3SHong Zhang 
23109573ac7SBarry Smith extern PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
23209573ac7SBarry Smith extern PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
23309573ac7SBarry Smith extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
23453565b12SHong Zhang extern PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(Mat,Mat,PetscReal,Mat*);
23509573ac7SBarry Smith extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
23653565b12SHong Zhang extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat,Mat,Mat);
23753565b12SHong Zhang extern PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2(Mat,Mat,Mat);
23853565b12SHong Zhang 
2392b8ad9a3SHong Zhang extern PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2402b8ad9a3SHong Zhang extern PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2415df89d91SHong Zhang 
24275648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
24375648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
24475648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
2456fc122caSHong Zhang extern PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2466fc122caSHong Zhang extern PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
2476fc122caSHong Zhang extern PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
248b9af6bddSHong Zhang extern PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat,ISColoring,MatTransposeColoring);
249b9af6bddSHong Zhang extern PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring,Mat,Mat);
250b9af6bddSHong Zhang extern PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring,Mat,Mat);
2515df89d91SHong Zhang 
25209573ac7SBarry Smith extern PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
25309573ac7SBarry Smith extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
25409573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
25509573ac7SBarry Smith extern PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
25609573ac7SBarry Smith extern PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
25709573ac7SBarry Smith extern PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
25809573ac7SBarry Smith extern PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
25909573ac7SBarry Smith extern PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
26009573ac7SBarry Smith extern PetscErrorCode MatDestroy_SeqAIJ(Mat);
26162fccb0cSShri Abhyankar extern PetscErrorCode MatSetUp_SeqAIJ(Mat);
26209573ac7SBarry Smith extern PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
2639af31e4aSHong Zhang 
26409573ac7SBarry Smith extern PetscErrorCode Mat_CheckInode(Mat,PetscBool );
26509573ac7SBarry Smith extern PetscErrorCode Mat_CheckInode_FactorLU(Mat,PetscBool );
266019b515eSShri Abhyankar 
26709573ac7SBarry Smith extern PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat,Mat,PetscInt*);
2689f5f6813SShri Abhyankar 
26997304618SKris Buschelman EXTERN_C_BEGIN
2707087cfbeSBarry Smith extern PetscErrorCode  MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*);
2717087cfbeSBarry Smith extern PetscErrorCode  MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*);
2727087cfbeSBarry Smith extern PetscErrorCode  MatConvert_SeqAIJ_SeqAIJPERM(Mat,const MatType,MatReuse,Mat*);
2737087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
2747087cfbeSBarry Smith extern PetscErrorCode  MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2758cdbd757SHong Zhang extern PetscErrorCode  MatMatMult_SeqDense_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2762b8ad9a3SHong Zhang extern PetscErrorCode  MatRARt_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
2777087cfbeSBarry Smith extern PetscErrorCode  MatCreate_SeqAIJ(Mat);
27897304618SKris Buschelman EXTERN_C_END
2797087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd_SeqAIJ(Mat,MatAssemblyType);
2807087cfbeSBarry Smith extern PetscErrorCode  MatDestroy_SeqAIJ(Mat);
2812f947c57SVictor Minden 
28270f19b1fSKris Buschelman 
283003131ecSBarry Smith /*
284003131ecSBarry Smith     PetscSparseDenseMinusDot - The inner kernel of triangular solves and Gauss-Siedel smoothing. \sum_i xv[i] * r[xi[i]] for CSR storage
285003131ecSBarry Smith 
286003131ecSBarry Smith   Input Parameters:
287003131ecSBarry Smith +  nnz - the number of entries
288003131ecSBarry Smith .  r - the array of vector values
289003131ecSBarry Smith .  xv - the matrix values for the row
290003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
291003131ecSBarry Smith 
292003131ecSBarry Smith   Output Parameter:
293003131ecSBarry Smith .  sum - negative the sum of results
294003131ecSBarry Smith 
295003131ecSBarry Smith   PETSc compile flags:
296b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -   don't use this; it changes nnz and hence is WRONG
297e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
298003131ecSBarry Smith 
299003131ecSBarry Smith .seealso: PetscSparseDensePlusDot()
300003131ecSBarry Smith 
301003131ecSBarry Smith */
302003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4
303003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
304003131ecSBarry Smith if (nnz > 0) {\
305003131ecSBarry Smith switch (nnz & 0x3) {\
306003131ecSBarry Smith case 3: sum -= *xv++ * r[*xi++];\
307003131ecSBarry Smith case 2: sum -= *xv++ * r[*xi++];\
308003131ecSBarry Smith case 1: sum -= *xv++ * r[*xi++];\
309003131ecSBarry Smith nnz -= 4;}\
310003131ecSBarry Smith while (nnz > 0) {\
311003131ecSBarry Smith sum -=  xv[0] * r[xi[0]] - xv[1] * r[xi[1]] -\
312003131ecSBarry Smith 	xv[2] * r[xi[2]] - xv[3] * r[xi[3]];\
313003131ecSBarry Smith xv  += 4; xi += 4; nnz -= 4; }}}
314003131ecSBarry Smith 
315003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
316003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
317003131ecSBarry Smith PetscInt __i,__i1,__i2;\
318003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
319003131ecSBarry Smith sum -= (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
320003131ecSBarry Smith if (nnz & 0x1) sum -= xv[__i] * r[xi[__i]];}
321003131ecSBarry Smith 
322003131ecSBarry Smith #else
323003131ecSBarry Smith #define PetscSparseDenseMinusDot(sum,r,xv,xi,nnz) {\
324003131ecSBarry Smith PetscInt __i;\
325003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum -= xv[__i] * r[xi[__i]];}
326003131ecSBarry Smith #endif
327003131ecSBarry Smith 
328003131ecSBarry Smith 
329003131ecSBarry Smith 
330003131ecSBarry Smith /*
331003131ecSBarry Smith     PetscSparseDensePlusDot - The inner kernel of matrix-vector product \sum_i xv[i] * r[xi[i]] for CSR storage
332003131ecSBarry Smith 
333003131ecSBarry Smith   Input Parameters:
334003131ecSBarry Smith +  nnz - the number of entries
335003131ecSBarry Smith .  r - the array of vector values
336003131ecSBarry Smith .  xv - the matrix values for the row
337003131ecSBarry Smith -  xi - the column indices of the nonzeros in the row
338003131ecSBarry Smith 
339003131ecSBarry Smith   Output Parameter:
340003131ecSBarry Smith .  sum - the sum of results
341003131ecSBarry Smith 
342003131ecSBarry Smith   PETSc compile flags:
343b2843cf1SBarry Smith +   PETSC_KERNEL_USE_UNROLL_4 -  don't use this; it changes nnz and hence is WRONG
344e0e43300SBarry Smith -   PETSC_KERNEL_USE_UNROLL_2 -
345003131ecSBarry Smith 
346003131ecSBarry Smith .seealso: PetscSparseDenseMinusDot()
347003131ecSBarry Smith 
348003131ecSBarry Smith */
349003131ecSBarry Smith #ifdef PETSC_KERNEL_USE_UNROLL_4
350003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
351003131ecSBarry Smith if (nnz > 0) {\
352003131ecSBarry Smith switch (nnz & 0x3) {\
353003131ecSBarry Smith case 3: sum += *xv++ * r[*xi++];\
354003131ecSBarry Smith case 2: sum += *xv++ * r[*xi++];\
355003131ecSBarry Smith case 1: sum += *xv++ * r[*xi++];\
356003131ecSBarry Smith nnz -= 4;}\
357003131ecSBarry Smith while (nnz > 0) {\
358003131ecSBarry Smith sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] +\
359003131ecSBarry Smith 	xv[2] * r[xi[2]] + xv[3] * r[xi[3]];\
360003131ecSBarry Smith xv  += 4; xi += 4; nnz -= 4; }}}
361003131ecSBarry Smith 
362003131ecSBarry Smith #elif defined(PETSC_KERNEL_USE_UNROLL_2)
363003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
364003131ecSBarry Smith PetscInt __i,__i1,__i2;\
365003131ecSBarry Smith for(__i=0;__i<nnz-1;__i+=2) {__i1 = xi[__i]; __i2=xi[__i+1];\
366003131ecSBarry Smith sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);}\
367003131ecSBarry Smith if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
368003131ecSBarry Smith 
369003131ecSBarry Smith #else
370003131ecSBarry Smith #define PetscSparseDensePlusDot(sum,r,xv,xi,nnz) {\
371003131ecSBarry Smith  PetscInt __i;\
372003131ecSBarry Smith for(__i=0;__i<nnz;__i++) sum += xv[__i] * r[xi[__i]];}
373003131ecSBarry Smith #endif
374003131ecSBarry Smith 
3752d40f771SBarry Smith #endif
376