Actual source code: baij.h
1: /* $Id: baij.h,v 1.35 2001/08/07 03:02:55 balay Exp $ */
3: #include src/mat/matimpl.h
8: /*
9: MATSEQBAIJ format - Block compressed row storage. The i[] and j[]
10: arrays start at 0.
11: */
13: typedef struct {
14: PetscTruth sorted; /* if true, rows are sorted by increasing columns */
15: PetscTruth roworiented; /* if true, row-oriented input, default */
16: int 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: int bs,bs2; /* block size, square of block size */
20: int mbs,nbs; /* rows/bs, columns/bs */
21: int nz,maxnz; /* nonzeros, allocated nonzeros */
22: int *diag; /* pointers to diagonal elements */
23: int *i; /* pointer to beginning of each row */
24: int *imax; /* maximum space allocated for each row */
25: int *ilen; /* actual length of each row */
26: int *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: int 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;
35: PetscTruth keepzeroedrows; /* if true, MatZeroRows() will not change nonzero structure */
37: #if defined(PETSC_USE_MAT_SINGLE)
38: int setvalueslen;
39: MatScalar *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied
40: before calling MatSetValuesXXX_SeqBAIJ_MatScalar() */
41: #endif
42: PetscTruth pivotinblocks; /* pivot inside factorization of each diagonal block */
44: int *xtoy,*xtoyB; /* map nonzero pattern of X into Y's, used by MatAXPY() */
45: Mat XtoY; /* used by MatAXPY() */
46: } Mat_SeqBAIJ;
48: EXTERN int MatILUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
49: EXTERN int MatDuplicate_SeqBAIJ(Mat,MatDuplicateOption,Mat*);
50: EXTERN int MatMarkDiagonal_SeqBAIJ(Mat);
52: EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
53: EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatFactorInfo*);
54: EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
55: EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
56: EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
57: EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec);
58: EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
59: EXTERN int MatScale_SeqBAIJ(const PetscScalar*,Mat);
60: EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *);
61: EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*);
62: EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec);
63: EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
64: EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
65: EXTERN int MatZeroEntries_SeqBAIJ(Mat);
67: EXTERN int MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(Mat);
68: EXTERN int MatSeqBAIJ_UpdateSolvers(Mat);
70: EXTERN int MatSolve_SeqBAIJ_Update(Mat,Vec,Vec);
71: EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
72: EXTERN int MatSolve_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
73: EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
74: EXTERN int MatSolve_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
75: EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
76: EXTERN int MatSolve_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
77: EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
78: EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
79: #if defined(PETSC_HAVE_SSE)
80: EXTERN int MatSolve_SeqBAIJ_4_SSE_Demotion(Mat,Vec,Vec);
81: EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion(Mat,Vec,Vec);
82: EXTERN int MatSolve_SeqBAIJ_4_NaturalOrdering_SSE_Demotion_usj(Mat,Vec,Vec);
83: #endif
84: EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
85: EXTERN int MatSolve_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
86: EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
87: EXTERN int MatSolve_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
88: EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
89: EXTERN int MatSolve_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
90: EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
92: EXTERN int MatSolveTranspose_SeqBAIJ_Update(Mat,Vec,Vec);
93: EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec);
94: EXTERN int MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
95: EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec);
96: EXTERN int MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
97: EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec);
98: EXTERN int MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
99: EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec);
100: EXTERN int MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
101: EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec);
102: EXTERN int MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
103: EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec);
104: EXTERN int MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
105: EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec);
106: EXTERN int MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
108: EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
109: EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
110: EXTERN int MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat,Mat*);
111: EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
112: EXTERN int MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat,Mat*);
113: EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
114: EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat,Mat*);
115: #if defined(PETSC_HAVE_SSE)
116: EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE(Mat,Mat*);
117: EXTERN int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat,Mat*);
118: #else
119: #endif
120: EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
121: EXTERN int MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat,Mat*);
122: EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
123: EXTERN int MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat,Mat*);
124: EXTERN int MatLUFactorNumeric_SeqBAIJ_7(Mat,Mat*);
125: EXTERN int MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat,Mat*);
126: EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
128: EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
129: EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
130: EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
131: EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
132: EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
133: EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
134: EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
135: EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
137: EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
138: EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
139: EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
140: EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
141: EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
142: EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
143: EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
144: EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
145: EXTERN int MatLoad_SeqBAIJ(PetscViewer,MatType,Mat*);
147: #endif