Actual source code: bdiag.h
4: #include include/private/matimpl.h
7: /*
8: Mat_SeqBDiag (MATSEQBDIAG) - block-diagonal format, where each diagonal
9: element consists of a square block of size bs x bs. Dense storage
10: within each block is in column-major order. The diagonals are the
11: full length of the matrix. As a special case, blocks of size bs=1
12: (scalars) are supported as well.
13: */
15: typedef struct {
16: PetscInt mblock,nblock; /* block rows and columns */
17: PetscInt nonew; /* if true, no new nonzeros allowed in matrix */
18: PetscInt nonew_diag; /* if true, no new diagonals allowed in matrix */
19: PetscInt nz,maxnz; /* nonzeros, allocated nonzeros */
20: PetscInt nd; /* number of block diagonals */
21: PetscInt mainbd; /* the number of the main block diagonal */
22: PetscInt *diag; /* value of (row-col)/bs for each diagonal */
23: PetscInt *bdlen; /* block-length of each diagonal */
24: PetscInt ndim; /* diagonals come from an ndim pde (if 0, ignore) */
25: PetscInt ndims[3]; /* sizes of the mesh if ndim > 0 */
26: PetscTruth user_alloc; /* true if the user provided the diagonals */
27: PetscInt *colloc; /* holds column locations if using MatGetRow */
28: PetscScalar **diagv; /* The actual diagonals */
29: PetscScalar *dvalue; /* Used to hold a row if MatGetRow is used */
30: PetscInt *pivot; /* pivots for LU factorization (temporary loc) */
31: PetscTruth roworiented; /* inputs to MatSetValue() are row oriented (default = 1) */
32: PetscInt reallocs; /* number of allocations during MatSetValues */
33: PetscScalar *solvework; /* work space for triangular solves for large block sizes */
34: } Mat_SeqBDiag;
36: EXTERN PetscErrorCode MatNorm_SeqBDiag_Columns(Mat,PetscReal*,PetscInt);
37: EXTERN PetscErrorCode MatMult_SeqBDiag_1(Mat,Vec,Vec);
38: EXTERN PetscErrorCode MatMult_SeqBDiag_2(Mat,Vec,Vec);
39: EXTERN PetscErrorCode MatMult_SeqBDiag_3(Mat,Vec,Vec);
40: EXTERN PetscErrorCode MatMult_SeqBDiag_4(Mat,Vec,Vec);
41: EXTERN PetscErrorCode MatMult_SeqBDiag_5(Mat,Vec,Vec);
42: EXTERN PetscErrorCode MatMult_SeqBDiag_N(Mat A,Vec,Vec);
43: EXTERN PetscErrorCode MatMultAdd_SeqBDiag_1(Mat,Vec,Vec,Vec);
44: EXTERN PetscErrorCode MatMultAdd_SeqBDiag_2(Mat,Vec,Vec,Vec);
45: EXTERN PetscErrorCode MatMultAdd_SeqBDiag_3(Mat,Vec,Vec,Vec);
46: EXTERN PetscErrorCode MatMultAdd_SeqBDiag_4(Mat,Vec,Vec,Vec);
47: EXTERN PetscErrorCode MatMultAdd_SeqBDiag_5(Mat,Vec,Vec,Vec);
48: EXTERN PetscErrorCode MatMultAdd_SeqBDiag_N(Mat A,Vec,Vec,Vec);
49: EXTERN PetscErrorCode MatMultTranspose_SeqBDiag_1(Mat,Vec,Vec);
50: EXTERN PetscErrorCode MatMultTranspose_SeqBDiag_N(Mat A,Vec,Vec);
51: EXTERN PetscErrorCode MatMultTransposeAdd_SeqBDiag_1(Mat,Vec,Vec,Vec);
52: EXTERN PetscErrorCode MatMultTransposeAdd_SeqBDiag_N(Mat A,Vec,Vec,Vec);
53: EXTERN PetscErrorCode MatSetValues_SeqBDiag_1(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
54: EXTERN PetscErrorCode MatSetValues_SeqBDiag_N(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
55: EXTERN PetscErrorCode MatGetValues_SeqBDiag_1(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],PetscScalar []);
56: EXTERN PetscErrorCode MatGetValues_SeqBDiag_N(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],PetscScalar []);
57: EXTERN PetscErrorCode MatRelax_SeqBDiag_1(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
58: EXTERN PetscErrorCode MatRelax_SeqBDiag_N(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
59: EXTERN PetscErrorCode MatView_SeqBDiag(Mat,PetscViewer);
60: EXTERN PetscErrorCode MatGetInfo_SeqBDiag(Mat,MatInfoType,MatInfo*);
61: EXTERN PetscErrorCode MatGetRow_SeqBDiag(Mat,PetscInt,PetscInt *,PetscInt **,PetscScalar **);
62: EXTERN PetscErrorCode MatRestoreRow_SeqBDiag(Mat,PetscInt,PetscInt *,PetscInt **,PetscScalar **);
63: EXTERN PetscErrorCode MatTranspose_SeqBDiag(Mat,Mat *);
64: EXTERN PetscErrorCode MatNorm_SeqBDiag(Mat,NormType,PetscReal *);
65: EXTERN PetscErrorCode MatLUFactorSymbolic_SeqBDiag(Mat,IS,IS,MatFactorInfo*,Mat*);
66: EXTERN PetscErrorCode MatILUFactorSymbolic_SeqBDiag(Mat,IS,IS,MatFactorInfo*,Mat*);
67: EXTERN PetscErrorCode MatILUFactor_SeqBDiag(Mat,IS,IS,MatFactorInfo*);
68: EXTERN PetscErrorCode MatLUFactorNumeric_SeqBDiag_N(Mat,MatFactorInfo*,Mat*);
69: EXTERN PetscErrorCode MatLUFactorNumeric_SeqBDiag_1(Mat,MatFactorInfo*,Mat*);
70: EXTERN PetscErrorCode MatSolve_SeqBDiag_1(Mat,Vec,Vec);
71: EXTERN PetscErrorCode MatSolve_SeqBDiag_2(Mat,Vec,Vec);
72: EXTERN PetscErrorCode MatSolve_SeqBDiag_3(Mat,Vec,Vec);
73: EXTERN PetscErrorCode MatSolve_SeqBDiag_4(Mat,Vec,Vec);
74: EXTERN PetscErrorCode MatSolve_SeqBDiag_5(Mat,Vec,Vec);
75: EXTERN PetscErrorCode MatSolve_SeqBDiag_N(Mat,Vec,Vec);
76: EXTERN PetscErrorCode MatLoad_SeqBDiag(PetscViewer, MatType,Mat*);
77: EXTERN PetscErrorCode MatGetRow_MPIBDiag(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
78: EXTERN PetscErrorCode MatRestoreRow_MPIBDiag(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
79: EXTERN PetscErrorCode MatScale_SeqBDiag(Mat,PetscScalar);
81: #endif