xref: /petsc/src/mat/impls/baij/seq/baij.h (revision a1d92eed66694224bc6faa1ab6e88dc1bb90150d)
173f4d377SMatthew Knepley /* $Id: baij.h,v 1.35 2001/08/07 03:02:55 balay Exp $ */
22593348eSBarry Smith 
335aab85fSBarry Smith #if !defined(__BAIJ_H)
435aab85fSBarry Smith #define __BAIJ_H
5da33ede1SBarry Smith #include "src/mat/matimpl.h"
6da33ede1SBarry Smith 
72593348eSBarry Smith 
82593348eSBarry Smith /*
935aab85fSBarry Smith   MATSEQBAIJ format - Block compressed row storage. The i[] and j[]
10e24b481bSBarry Smith   arrays start at 0.
112593348eSBarry Smith */
122593348eSBarry Smith 
136c6c5352SBarry Smith /* This header is shared by the SeqSBAIJ matrix */
146c6c5352SBarry Smith #define SEQBAIJHEADER \
156c6c5352SBarry Smith   PetscTruth       sorted;       /* if true, rows are sorted by increasing columns */                \
166c6c5352SBarry Smith   PetscTruth       roworiented;  /* if true, row-oriented input, default */                          \
176c6c5352SBarry Smith   int              nonew;        /* 1 don't add new nonzeros, -1 generate error on new */            \
186c6c5352SBarry Smith   PetscTruth       singlemalloc; /* if true a, i, and j have been obtained with                      \
196c6c5352SBarry Smith                                         one big malloc */                                            \
206c6c5352SBarry Smith   int              bs,bs2;       /* block size, square of block size */                              \
216c6c5352SBarry Smith   int              mbs,nbs;      /* rows/bs, columns/bs */                                           \
226c6c5352SBarry Smith   int              nz,maxnz;     /* nonzeros, allocated nonzeros */                                  \
236c6c5352SBarry Smith   int              *diag;        /* pointers to diagonal elements */                                 \
246c6c5352SBarry Smith   int              *i;           /* pointer to beginning of each row */                              \
256c6c5352SBarry Smith   int              *imax;        /* maximum space allocated for each row */                          \
266c6c5352SBarry Smith   int              *ilen;        /* actual length of each row */                                     \
276c6c5352SBarry Smith   int              *j;           /* column values: j + i[k] - 1 is start of row k */                 \
286c6c5352SBarry Smith   MatScalar        *a;           /* nonzero elements */                                              \
296c6c5352SBarry Smith   IS               row,col,icol; /* index sets, used for reorderings */                              \
306c6c5352SBarry Smith   PetscScalar      *solve_work;  /* work space used in MatSolve */                                   \
316c6c5352SBarry Smith   int              reallocs;     /* number of mallocs done during MatSetValues()                     \
326c6c5352SBarry Smith                                     as more values are set then were preallocated */                 \
336c6c5352SBarry Smith   PetscScalar      *mult_work;   /* work array for matrix vector product*/                           \
346c6c5352SBarry Smith   PetscScalar      *saved_values;                                                                    \
356c6c5352SBarry Smith                                                                                                      \
366c6c5352SBarry Smith   PetscTruth       keepzeroedrows; /* if true, MatZeroRows() will not change nonzero structure */    \
376c6c5352SBarry Smith                                                                                                      \
386c6c5352SBarry Smith   int              setvalueslen;   /* only used for single precision */                              \
396c6c5352SBarry Smith   MatScalar        *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied   \
406c6c5352SBarry Smith                                       before calling MatSetValuesXXX_SeqBAIJ_MatScalar() */          \
416c6c5352SBarry Smith                                                                                                      \
426c6c5352SBarry Smith   PetscTruth       pivotinblocks;  /* pivot inside factorization of each diagonal block */           \
436c6c5352SBarry Smith                                                                                                      \
446c6c5352SBarry Smith   int              *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */      \
45*a1d92eedSBarry Smith   Mat              XtoY;             /* used by MatAXPY() */                                         \
46*a1d92eedSBarry Smith   PetscScalar      *idiag;           /* inverse of block diagonal  */
476c6c5352SBarry Smith 
486c6c5352SBarry Smith typedef struct {
496c6c5352SBarry Smith   SEQBAIJHEADER
5035aab85fSBarry Smith } Mat_SeqBAIJ;
512593348eSBarry Smith 
52b380c88cSHong Zhang EXTERN int MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
53ca44d042SBarry Smith EXTERN int MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*);
54ca44d042SBarry Smith EXTERN int MatMarkDiagonal_SeqBAIJ(Mat);
55be3590efSBarry Smith 
56b380c88cSHong Zhang EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
57b380c88cSHong Zhang EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*);
5854138f6bSKris Buschelman EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
5954138f6bSKris Buschelman EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
60268466fbSBarry Smith EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
6154138f6bSKris Buschelman EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
6254138f6bSKris Buschelman EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
63268466fbSBarry Smith EXTERN int MatScale_SeqBAIJ(const PetscScalar*,Mat);
6454138f6bSKris Buschelman EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
6554138f6bSKris Buschelman EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
6654138f6bSKris Buschelman EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec);
6754138f6bSKris Buschelman EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
6854138f6bSKris Buschelman EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
6954138f6bSKris Buschelman EXTERN int MatZeroEntries_SeqBAIJ(Mat);
7054138f6bSKris Buschelman 
7154138f6bSKris Buschelman EXTERN int MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat);
7254138f6bSKris Buschelman EXTERN int MatSeqBAIJ_UpdateSolvers(Mat);
7354138f6bSKris Buschelman 
7454138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_Update(Mat,Vec,Vec);
7554138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
76e37c484bSBarry Smith EXTERN int MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
7754138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
78ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
7954138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
80ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
8154138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
82ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
8354138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
8454138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec);
854110d75bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec);
860959a075SKris Buschelman EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec);
874110d75bSKris Buschelman #endif
8854138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
89ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
9054138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
91ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
9254138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
93ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
9454138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
9554138f6bSKris Buschelman 
9654138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec);
9754138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
9854138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
9954138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
10054138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
10154138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
10254138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
10354138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
10454138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
10554138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
10654138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
10754138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
10854138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
10954138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
110ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
11154138f6bSKris Buschelman 
11254138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
11354138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
11454138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,Mat*);
11554138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
11654138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,Mat*);
11754138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
11854138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,Mat*);
11954138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
12054138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,Mat*);
1210959a075SKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,Mat*);
12254138f6bSKris Buschelman #else
12354138f6bSKris Buschelman #endif
12454138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
12554138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,Mat*);
12654138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
12754138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,Mat*);
12854138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_7(Mat,Mat*);
12954138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,Mat*);
13054138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
13154138f6bSKris Buschelman 
13254138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
13354138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
13454138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
13554138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
13654138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
13754138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
13854138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
13954138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
14054138f6bSKris Buschelman 
14154138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
14254138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
14354138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
14454138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
14554138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
14654138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
14754138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
14854138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
149f248c16bSBarry Smith EXTERN int MatLoad_SeqBAIJ(PetscViewer,const MatType,Mat*);
1502593348eSBarry Smith 
1512593348eSBarry Smith #endif
152