xref: /petsc/src/mat/impls/baij/seq/baij.h (revision 6c6c5352cdaa6aa3a8bd2967d4dd73a4e9ba1251)
173f4d377SMatthew Knepley /* $Id: baij.h,v 1.35 2001/08/07 03:02:55 balay Exp $ */
22593348eSBarry Smith 
370f55243SBarry Smith #include "src/mat/matimpl.h"
42593348eSBarry Smith 
535aab85fSBarry Smith #if !defined(__BAIJ_H)
635aab85fSBarry Smith #define __BAIJ_H
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 
13*6c6c5352SBarry Smith /* This header is shared by the SeqSBAIJ matrix */
14*6c6c5352SBarry Smith #define SEQBAIJHEADER \
15*6c6c5352SBarry Smith   PetscTruth       sorted;       /* if true, rows are sorted by increasing columns */                \
16*6c6c5352SBarry Smith   PetscTruth       roworiented;  /* if true, row-oriented input, default */                          \
17*6c6c5352SBarry Smith   int              nonew;        /* 1 don't add new nonzeros, -1 generate error on new */            \
18*6c6c5352SBarry Smith   PetscTruth       singlemalloc; /* if true a, i, and j have been obtained with                      \
19*6c6c5352SBarry Smith                                         one big malloc */                                            \
20*6c6c5352SBarry Smith   int              bs,bs2;       /* block size, square of block size */                              \
21*6c6c5352SBarry Smith   int              mbs,nbs;      /* rows/bs, columns/bs */                                           \
22*6c6c5352SBarry Smith   int              nz,maxnz;     /* nonzeros, allocated nonzeros */                                  \
23*6c6c5352SBarry Smith   int              *diag;        /* pointers to diagonal elements */                                 \
24*6c6c5352SBarry Smith   int              *i;           /* pointer to beginning of each row */                              \
25*6c6c5352SBarry Smith   int              *imax;        /* maximum space allocated for each row */                          \
26*6c6c5352SBarry Smith   int              *ilen;        /* actual length of each row */                                     \
27*6c6c5352SBarry Smith   int              *j;           /* column values: j + i[k] - 1 is start of row k */                 \
28*6c6c5352SBarry Smith   MatScalar        *a;           /* nonzero elements */                                              \
29*6c6c5352SBarry Smith   IS               row,col,icol; /* index sets, used for reorderings */                              \
30*6c6c5352SBarry Smith   PetscScalar      *solve_work;  /* work space used in MatSolve */                                   \
31*6c6c5352SBarry Smith   int              reallocs;     /* number of mallocs done during MatSetValues()                     \
32*6c6c5352SBarry Smith                                     as more values are set then were preallocated */                 \
33*6c6c5352SBarry Smith   PetscScalar      *mult_work;   /* work array for matrix vector product*/                           \
34*6c6c5352SBarry Smith   PetscScalar      *saved_values;                                                                    \
35*6c6c5352SBarry Smith                                                                                                      \
36*6c6c5352SBarry Smith   PetscTruth       keepzeroedrows; /* if true, MatZeroRows() will not change nonzero structure */    \
37*6c6c5352SBarry Smith                                                                                                      \
38*6c6c5352SBarry Smith   int              setvalueslen;   /* only used for single precision */                              \
39*6c6c5352SBarry Smith   MatScalar        *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied   \
40*6c6c5352SBarry Smith                                       before calling MatSetValuesXXX_SeqBAIJ_MatScalar() */          \
41*6c6c5352SBarry Smith                                                                                                      \
42*6c6c5352SBarry Smith   PetscTruth       pivotinblocks;  /* pivot inside factorization of each diagonal block */           \
43*6c6c5352SBarry Smith                                                                                                      \
44*6c6c5352SBarry Smith   int              *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */      \
45c4319e64SHong Zhang   Mat              XtoY;             /* used by MatAXPY() */
46*6c6c5352SBarry Smith 
47*6c6c5352SBarry Smith typedef struct {
48*6c6c5352SBarry Smith   SEQBAIJHEADER
4935aab85fSBarry Smith } Mat_SeqBAIJ;
502593348eSBarry Smith 
51b380c88cSHong Zhang EXTERN int MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
52ca44d042SBarry Smith EXTERN int MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*);
53ca44d042SBarry Smith EXTERN int MatMarkDiagonal_SeqBAIJ(Mat);
54be3590efSBarry Smith 
55b380c88cSHong Zhang EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
56b380c88cSHong Zhang EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*);
5754138f6bSKris Buschelman EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
5854138f6bSKris Buschelman EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
59268466fbSBarry Smith EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
6054138f6bSKris Buschelman EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
6154138f6bSKris Buschelman EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
62268466fbSBarry Smith EXTERN int MatScale_SeqBAIJ(const PetscScalar*,Mat);
6354138f6bSKris Buschelman EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
6454138f6bSKris Buschelman EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
6554138f6bSKris Buschelman EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec);
6654138f6bSKris Buschelman EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
6754138f6bSKris Buschelman EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
6854138f6bSKris Buschelman EXTERN int MatZeroEntries_SeqBAIJ(Mat);
6954138f6bSKris Buschelman 
7054138f6bSKris Buschelman EXTERN int MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat);
7154138f6bSKris Buschelman EXTERN int MatSeqBAIJ_UpdateSolvers(Mat);
7254138f6bSKris Buschelman 
7354138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_Update(Mat,Vec,Vec);
7454138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
75e37c484bSBarry Smith EXTERN int MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
7654138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
77ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
7854138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
79ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
8054138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
81ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
8254138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
8354138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec);
844110d75bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec);
850959a075SKris Buschelman EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec);
864110d75bSKris Buschelman #endif
8754138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
88ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
8954138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
90ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
9154138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
92ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
9354138f6bSKris Buschelman EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
9454138f6bSKris Buschelman 
9554138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec);
9654138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
9754138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
9854138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
9954138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
10054138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
10154138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
10254138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
10354138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
10454138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
10554138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
10654138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
10754138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
10854138f6bSKris Buschelman EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
109ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
11054138f6bSKris Buschelman 
11154138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
11254138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
11354138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,Mat*);
11454138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
11554138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,Mat*);
11654138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
11754138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,Mat*);
11854138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
11954138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,Mat*);
1200959a075SKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,Mat*);
12154138f6bSKris Buschelman #else
12254138f6bSKris Buschelman #endif
12354138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
12454138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,Mat*);
12554138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
12654138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,Mat*);
12754138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_7(Mat,Mat*);
12854138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,Mat*);
12954138f6bSKris Buschelman EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
13054138f6bSKris Buschelman 
13154138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
13254138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
13354138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
13454138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
13554138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
13654138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
13754138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
13854138f6bSKris Buschelman EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
13954138f6bSKris Buschelman 
14054138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
14154138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
14254138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
14354138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
14454138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
14554138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
14654138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
14754138f6bSKris Buschelman EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
148f248c16bSBarry Smith EXTERN int MatLoad_SeqBAIJ(PetscViewer,const MatType,Mat*);
1492593348eSBarry Smith 
1502593348eSBarry Smith #endif
151