xref: /petsc/src/mat/impls/baij/seq/baij.h (revision 8bc8193efbc389280f83b3d41dffa9e2d23e2ace)
1 
2 #if !defined(__BAIJ_H)
3 #define __BAIJ_H
4 #include "src/mat/matimpl.h"
5 
6 
7 /*
8   MATSEQBAIJ format - Block compressed row storage. The i[] and j[]
9   arrays start at 0.
10 */
11 
12 /* This header is shared by the SeqSBAIJ matrix */
13 #define SEQBAIJHEADER \
14   PetscTruth       sorted;       /* if true, rows are sorted by increasing columns */                \
15   PetscTruth       roworiented;  /* if true, row-oriented input, default */                          \
16   PetscInt         nonew;        /* 1 don't add new nonzeros, -1 generate error on new */            \
17   PetscTruth       singlemalloc; /* if true a, i, and j have been obtained with                      \
18                                         one big malloc */                                            \
19   PetscInt         bs2;          /*  square of block size */                                         \
20   PetscInt         mbs,nbs;      /* rows/bs, columns/bs */                                           \
21   PetscInt         nz,maxnz;     /* nonzeros, allocated nonzeros */                                  \
22   PetscInt         *diag;        /* pointers to diagonal elements */                                 \
23   PetscInt         *i;           /* pointer to beginning of each row */                              \
24   PetscInt         *imax;        /* maximum space allocated for each row */                          \
25   PetscInt         *ilen;        /* actual length of each row */                                     \
26   PetscInt         *j;           /* column values: j + i[k] - 1 is start of row k */                 \
27   MatScalar        *a;           /* nonzero elements */                                              \
28   IS               row,col,icol; /* index sets, used for reorderings */                              \
29   PetscScalar      *solve_work;  /* work space used in MatSolve */                                   \
30   PetscInt         reallocs;     /* number of mallocs done during MatSetValues()                     \
31                                     as more values are set then were preallocated */                 \
32   PetscScalar      *mult_work;   /* work array for matrix vector product*/                           \
33   PetscScalar      *saved_values;                                                                    \
34                                                                                                      \
35   PetscTruth       keepzeroedrows; /* if true, MatZeroRows() will not change nonzero structure */    \
36                                                                                                      \
37   PetscInt         setvalueslen;   /* only used for single precision */                              \
38   MatScalar        *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied   \
39                                       before calling MatSetValuesXXX_SeqBAIJ_MatScalar() */          \
40                                                                                                      \
41   PetscTruth       pivotinblocks;  /* pivot inside factorization of each diagonal block */           \
42                                                                                                      \
43   PetscInt         *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */      \
44   Mat              XtoY;             /* used by MatAXPY() */                                         \
45   PetscScalar      *idiag;           /* inverse of block diagonal  */                                \
46   PetscTruth       idiagvalid;       /* if above has correct/current values */
47 
48 typedef struct {
49   SEQBAIJHEADER
50 } Mat_SeqBAIJ;
51 
52 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
53 EXTERN PetscErrorCode MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*);
54 EXTERN PetscErrorCode MatMissingDiagonal_SeqBAIJ(Mat);
55 EXTERN PetscErrorCode MatMarkDiagonal_SeqBAIJ(Mat);
56 
57 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
58 EXTERN PetscErrorCode MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*);
59 EXTERN PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat,PetscInt,IS*,PetscInt);
60 EXTERN PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
61 EXTERN PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
62 EXTERN PetscErrorCode MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
63 EXTERN PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
64 EXTERN PetscErrorCode MatScale_SeqBAIJ(const PetscScalar*,Mat);
65 EXTERN PetscErrorCode MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
66 EXTERN PetscErrorCode MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
67 EXTERN PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat,Vec);
68 EXTERN PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
69 EXTERN PetscErrorCode MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
70 EXTERN PetscErrorCode MatZeroEntries_SeqBAIJ(Mat);
71 
72 EXTERN PetscErrorCode MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat);
73 EXTERN PetscErrorCode MatSeqBAIJ_UpdateSolvers(Mat);
74 
75 EXTERN PetscErrorCode MatSolve_SeqBAIJ_Update(Mat,Vec,Vec);
76 EXTERN PetscErrorCode MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
77 EXTERN PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
78 EXTERN PetscErrorCode MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
79 EXTERN PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
80 EXTERN PetscErrorCode MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
81 EXTERN PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
82 EXTERN PetscErrorCode MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
83 EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
84 #if defined(PETSC_HAVE_SSE)
85 EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec);
86 EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec);
87 EXTERN PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec);
88 #endif
89 EXTERN PetscErrorCode MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
90 EXTERN PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
91 EXTERN PetscErrorCode MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
92 EXTERN PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
93 EXTERN PetscErrorCode MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
94 EXTERN PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
95 EXTERN PetscErrorCode MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
96 
97 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec);
98 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
99 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
100 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
101 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
102 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
103 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
104 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
105 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
106 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
107 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
108 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
109 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
110 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
111 EXTERN PetscErrorCode MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
112 
113 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
114 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
115 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,Mat*);
116 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
117 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,Mat*);
118 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
119 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,Mat*);
120 #if defined(PETSC_HAVE_SSE)
121 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,Mat*);
122 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,Mat*);
123 #else
124 #endif
125 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
126 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,Mat*);
127 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
128 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,Mat*);
129 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7(Mat,Mat*);
130 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,Mat*);
131 EXTERN PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
132 
133 EXTERN PetscErrorCode MatMult_SeqBAIJ_1(Mat,Vec,Vec);
134 EXTERN PetscErrorCode MatMult_SeqBAIJ_2(Mat,Vec,Vec);
135 EXTERN PetscErrorCode MatMult_SeqBAIJ_3(Mat,Vec,Vec);
136 EXTERN PetscErrorCode MatMult_SeqBAIJ_4(Mat,Vec,Vec);
137 EXTERN PetscErrorCode MatMult_SeqBAIJ_5(Mat,Vec,Vec);
138 EXTERN PetscErrorCode MatMult_SeqBAIJ_6(Mat,Vec,Vec);
139 EXTERN PetscErrorCode MatMult_SeqBAIJ_7(Mat,Vec,Vec);
140 EXTERN PetscErrorCode MatMult_SeqBAIJ_N(Mat,Vec,Vec);
141 
142 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
143 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
144 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
145 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
146 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
147 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
148 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
149 EXTERN PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
150 EXTERN PetscErrorCode MatLoad_SeqBAIJ(PetscViewer,const MatType,Mat*);
151 
152 #endif
153