xref: /petsc/src/mat/impls/baij/seq/baij.h (revision f349cc8147fd8be7fc2bb8bc6274b020a4c1628f)
12593348eSBarry Smith 
235aab85fSBarry Smith #if !defined(__BAIJ_H)
335aab85fSBarry Smith #define __BAIJ_H
4da33ede1SBarry Smith #include "src/mat/matimpl.h"
5da33ede1SBarry Smith 
62593348eSBarry Smith 
72593348eSBarry Smith /*
835aab85fSBarry Smith   MATSEQBAIJ format - Block compressed row storage. The i[] and j[]
9e24b481bSBarry Smith   arrays start at 0.
102593348eSBarry Smith */
112593348eSBarry Smith 
126c6c5352SBarry Smith /* This header is shared by the SeqSBAIJ matrix */
136c6c5352SBarry Smith #define SEQBAIJHEADER \
146c6c5352SBarry Smith   PetscTruth       sorted;       /* if true, rows are sorted by increasing columns */                \
156c6c5352SBarry Smith   PetscTruth       roworiented;  /* if true, row-oriented input, default */                          \
16690b6cddSBarry Smith   PetscInt         nonew;        /* 1 don't add new nonzeros, -1 generate error on new */            \
176c6c5352SBarry Smith   PetscTruth       singlemalloc; /* if true a, i, and j have been obtained with                      \
186c6c5352SBarry Smith                                         one big malloc */                                            \
19521d7252SBarry Smith   PetscInt         bs2;          /*  square of block size */                                         \
20690b6cddSBarry Smith   PetscInt         mbs,nbs;      /* rows/bs, columns/bs */                                           \
21690b6cddSBarry Smith   PetscInt         nz,maxnz;     /* nonzeros, allocated nonzeros */                                  \
22690b6cddSBarry Smith   PetscInt         *diag;        /* pointers to diagonal elements */                                 \
23690b6cddSBarry Smith   PetscInt         *i;           /* pointer to beginning of each row */                              \
24690b6cddSBarry Smith   PetscInt         *imax;        /* maximum space allocated for each row */                          \
25690b6cddSBarry Smith   PetscInt         *ilen;        /* actual length of each row */                                     \
26690b6cddSBarry Smith   PetscInt         *j;           /* column values: j + i[k] - 1 is start of row k */                 \
276c6c5352SBarry Smith   MatScalar        *a;           /* nonzero elements */                                              \
286c6c5352SBarry Smith   IS               row,col,icol; /* index sets, used for reorderings */                              \
296c6c5352SBarry Smith   PetscScalar      *solve_work;  /* work space used in MatSolve */                                   \
30690b6cddSBarry Smith   PetscInt         reallocs;     /* number of mallocs done during MatSetValues()                     \
316c6c5352SBarry Smith                                     as more values are set then were preallocated */                 \
326c6c5352SBarry Smith   PetscScalar      *mult_work;   /* work array for matrix vector product*/                           \
336c6c5352SBarry Smith   PetscScalar      *saved_values;                                                                    \
346c6c5352SBarry Smith                                                                                                      \
356c6c5352SBarry Smith   PetscTruth       keepzeroedrows; /* if true, MatZeroRows() will not change nonzero structure */    \
366ad2eaddSHong Zhang   Mat              sbaijMat;         /* mat in sbaij format */                                       \
376c6c5352SBarry Smith                                                                                                      \
38690b6cddSBarry Smith   PetscInt         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                                                                                                      \
44690b6cddSBarry Smith   PetscInt         *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */      \
45a1d92eedSBarry Smith   Mat              XtoY;             /* used by MatAXPY() */                                         \
46b01c7715SBarry Smith   PetscScalar      *idiag;           /* inverse of block diagonal  */                                \
4773e7a558SHong Zhang   PetscTruth       idiagvalid;       /* if above has correct/current values */                       \
4873e7a558SHong Zhang   Mat_CompressedRow compressedrow;   /* use compressed row format */
496c6c5352SBarry Smith 
506c6c5352SBarry Smith typedef struct {
516c6c5352SBarry Smith   SEQBAIJHEADER
5235aab85fSBarry Smith } Mat_SeqBAIJ;
532593348eSBarry Smith 
54*f349cc81SSatish Balay EXTERN_C_BEGIN
55ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqBAIJSetPreallocation_SeqBAIJ(Mat,PetscInt,PetscInt,PetscInt*);
56*f349cc81SSatish Balay EXTERN_C_END
57dfbe8321SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
58c05c3958SHong Zhang EXTERN PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat,IS,MatFactorInfo*,Mat *);
59c05c3958SHong Zhang EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat,IS,MatFactorInfo*,Mat*);
60af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat,MatFactorInfo*,Mat*);
61af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
62dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*);
63dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqBAIJ(Mat);
64dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqBAIJ(Mat);
65be3590efSBarry Smith 
66dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
67dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*);
68690b6cddSBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat,PetscInt,IS*,PetscInt);
69690b6cddSBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
70690b6cddSBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
71dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
72dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
73dfbe8321SBarry Smith EXTERN PetscErrorCode MatScale_SeqBAIJ(const PetscScalar*,Mat);
74dfbe8321SBarry Smith EXTERN PetscErrorCode MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
75dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
76dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat,Vec);
77dfbe8321SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
78dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
79dfbe8321SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqBAIJ(Mat);
8054138f6bSKris Buschelman 
81dfbe8321SBarry Smith EXTERN PetscErrorCode MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat);
82dfbe8321SBarry Smith EXTERN PetscErrorCode MatSeqBAIJ_UpdateSolvers(Mat);
8354138f6bSKris Buschelman 
84dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_Update(Mat,Vec,Vec);
85dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
86dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
87dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
88dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
89dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
90dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
91dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
92dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
9354138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
94dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec);
95dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec);
96dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec);
974110d75bSKris Buschelman #endif
98dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
99dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
100dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
101dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
102dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
103dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
104dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10554138f6bSKris Buschelman 
106dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec);
107dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
108dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
109dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
110dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
111dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
112dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
113dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
114dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
115dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
116dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
117dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
118dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
119dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
120dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
12154138f6bSKris Buschelman 
122af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat,MatFactorInfo*,Mat*);
123af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat,MatFactorInfo*,Mat*);
124af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
125af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat,MatFactorInfo*,Mat*);
126af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
127af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat,MatFactorInfo*,Mat*);
128af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
12954138f6bSKris Buschelman #if defined(PETSC_HAVE_SSE)
130af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,MatFactorInfo*,Mat*);
131af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,MatFactorInfo*,Mat*);
13254138f6bSKris Buschelman #else
13354138f6bSKris Buschelman #endif
134af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat,MatFactorInfo*,Mat*);
135af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
136af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat,MatFactorInfo*,Mat*);
137af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
138af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7(Mat,MatFactorInfo*,Mat*);
139af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,MatFactorInfo*,Mat*);
140af281ebdSHong Zhang EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat,MatFactorInfo*,Mat*);
14154138f6bSKris Buschelman 
142dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_1(Mat,Vec,Vec);
143dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_2(Mat,Vec,Vec);
144dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_3(Mat,Vec,Vec);
145dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_4(Mat,Vec,Vec);
146dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_5(Mat,Vec,Vec);
147dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_6(Mat,Vec,Vec);
148dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_7(Mat,Vec,Vec);
149dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqBAIJ_N(Mat,Vec,Vec);
15054138f6bSKris Buschelman 
151dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
152dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
153dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
154dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
155dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
156dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
157dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
158dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
159dfbe8321SBarry Smith EXTERN PetscErrorCode MatLoad_SeqBAIJ(PetscViewer,const MatType,Mat*);
1602593348eSBarry Smith 
1612593348eSBarry Smith #endif
162