Actual source code: dense.h
4: #include include/private/matimpl.h
7: /*
8: MATSEQDENSE format - conventional dense Fortran storage (by columns)
9: */
11: typedef struct {
12: PetscScalar *v; /* matrix elements */
13: PetscTruth roworiented; /* if true, row oriented input (default) */
14: PetscInt pad; /* padding */
15: PetscBLASInt *pivots; /* pivots in LU factorization */
16: PetscBLASInt lda; /* Lapack leading dimension of data */
17: PetscTruth changelda; /* change lda on resize? Default unless user set lda */
18: PetscBLASInt Mmax,Nmax; /* indicates the largest dimensions of data possible */
19: PetscTruth user_alloc; /* true if the user provided the dense data */
20: } Mat_SeqDense;
22: EXTERN PetscErrorCode MatMult_SeqDense(Mat A,Vec,Vec);
23: EXTERN PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec,Vec,Vec);
24: EXTERN PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec,Vec);
25: EXTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec,Vec,Vec);
26: EXTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
27: EXTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
28: EXTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
29: EXTERN PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
30: EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
31: EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
32: EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
33: EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
34: EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
36: #endif