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