xref: /petsc/src/mat/impls/baij/seq/baij.h (revision e6b907ac9e74227fac7d27fb51c2636d4aaee38b)
12593348eSBarry Smith 
235aab85fSBarry Smith #if !defined(__BAIJ_H)
335aab85fSBarry Smith #define __BAIJ_H
4da33ede1SBarry Smith #include "src/mat/matimpl.h"
5*e6b907acSBarry Smith #include "src/mat/impls/aij/seq/aij.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 \
15521d7252SBarry Smith   PetscInt         bs2;              /*  square of block size */                                     \
16690b6cddSBarry Smith   PetscInt         mbs,nbs;          /* rows/bs, columns/bs */                                       \
176c6c5352SBarry Smith   PetscScalar      *mult_work;       /* work array for matrix vector product*/                       \
186c6c5352SBarry Smith   PetscScalar      *saved_values;                                                                    \
196c6c5352SBarry Smith                                                                                                      \
206ad2eaddSHong Zhang   Mat              sbaijMat;         /* mat in sbaij format */                                       \
216c6c5352SBarry Smith                                                                                                      \
22690b6cddSBarry Smith   PetscInt         setvalueslen;     /* only used for single precision */                            \
236c6c5352SBarry Smith   MatScalar        *setvaluescopy;   /* area double precision values in MatSetValuesXXX() are copied \
246c6c5352SBarry Smith                                       before calling MatSetValuesXXX_SeqBAIJ_MatScalar() */          \
256c6c5352SBarry Smith                                                                                                      \
266c6c5352SBarry Smith   PetscTruth       pivotinblocks;    /* pivot inside factorization of each diagonal block */         \
276c6c5352SBarry Smith                                                                                                      \
28b01c7715SBarry Smith   PetscScalar      *idiag;           /* inverse of block diagonal  */                                \
29*e6b907acSBarry Smith   PetscTruth       idiagvalid       /* if above has correct/current values */
30*e6b907acSBarry Smith 
316c6c5352SBarry Smith 
326c6c5352SBarry Smith typedef struct {
33*e6b907acSBarry Smith   SEQAIJHEADER;
34*e6b907acSBarry Smith   SEQBAIJHEADER;
3535aab85fSBarry Smith } Mat_SeqBAIJ;
362593348eSBarry Smith 
37f349cc81SSatish Balay EXTERN_C_BEGIN
38ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqBAIJSetPreallocation_SeqBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
39f349cc81SSatish Balay EXTERN_C_END
40dfbe8321SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
41c05c3958SHong Zhang EXTERN PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat,IS,MatFactorInfo*,Mat *);
42c05c3958SHong Zhang EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat,IS,MatFactorInfo*,Mat*);
43af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat,MatFactorInfo*,Mat*);
44af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
45dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*);
46dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqBAIJ(Mat);
47dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqBAIJ(Mat);
48be3590efSBarry Smith 
49dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
50dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*);
51690b6cddSBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat,PetscInt,IS*,PetscInt);
52690b6cddSBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
53690b6cddSBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
54dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
55dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
56f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatScale_SeqBAIJ(Mat,PetscScalar);
57dfbe8321SBarry Smith EXTERN PetscErrorCode MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
58dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
59dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat,Vec);
60dfbe8321SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
61dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
62dfbe8321SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqBAIJ(Mat);
6354138f6bSKris Buschelman 
64dfbe8321SBarry Smith EXTERN PetscErrorCode MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat);
65dfbe8321SBarry Smith EXTERN PetscErrorCode MatSeqBAIJ_UpdateSolvers(Mat);
6654138f6bSKris Buschelman 
67dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_Update(Mat,Vec,Vec);
68dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
69dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
70dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
71dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
72dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
73dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
74dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
75dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
7654138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
77dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec);
78dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec);
79dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec);
804110d75bSKris Buschelman #endif
81dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
82dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
83dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
84dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
85dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
86dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
87dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
8854138f6bSKris Buschelman 
89dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec);
90dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
91dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
92dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
93dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
94dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
95dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
96dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
97dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
98dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
99dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
100dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
101dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
102dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
103dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
10454138f6bSKris Buschelman 
105af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat,MatFactorInfo*,Mat*);
106af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat,MatFactorInfo*,Mat*);
107af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
108af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat,MatFactorInfo*,Mat*);
109af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
110af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat,MatFactorInfo*,Mat*);
111af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
11254138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
113af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,MatFactorInfo*,Mat*);
114af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,MatFactorInfo*,Mat*);
11554138f6bSKris Buschelman #else
11654138f6bSKris Buschelman #endif
117af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat,MatFactorInfo*,Mat*);
118af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
119af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat,MatFactorInfo*,Mat*);
120af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
121af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7(Mat,MatFactorInfo*,Mat*);
122af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
123af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat,MatFactorInfo*,Mat*);
12454138f6bSKris Buschelman 
125dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_1(Mat,Vec,Vec);
126dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_2(Mat,Vec,Vec);
127dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_3(Mat,Vec,Vec);
128dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_4(Mat,Vec,Vec);
129dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_5(Mat,Vec,Vec);
130dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_6(Mat,Vec,Vec);
131dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_7(Mat,Vec,Vec);
132dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_N(Mat,Vec,Vec);
13354138f6bSKris Buschelman 
134dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
135dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
136dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
137dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
138dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
139dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
140dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
141dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
142f69a0ea3SMatthew Knepley EXTERN PetscErrorCode MatLoad_SeqBAIJ(PetscViewer, MatType,Mat*);
1432593348eSBarry Smith 
1442593348eSBarry Smith #endif
145