xref: /petsc/src/mat/impls/dense/seq/dense.h (revision 4ae313f42201736c22765e1d4ea08232bbac4af6)
1f6f390a9SLois Curfman McInnes 
2f6f390a9SLois Curfman McInnes #if !defined(__DENSE_H)
3f6f390a9SLois Curfman McInnes #define __DENSE_H
4da33ede1SBarry Smith #include "src/mat/matimpl.h"
5da33ede1SBarry Smith 
6f6f390a9SLois Curfman McInnes 
7f6f390a9SLois Curfman McInnes /*
8f6f390a9SLois Curfman McInnes   MATSEQDENSE format - conventional dense Fortran storage (by columns)
9f6f390a9SLois Curfman McInnes */
10f6f390a9SLois Curfman McInnes 
11f6f390a9SLois Curfman McInnes typedef struct {
12ea709b57SSatish Balay   PetscScalar  *v;                /* matrix elements */
13273d9f13SBarry Smith   PetscTruth   roworiented;       /* if true, row oriented input (default) */
1413f74950SBarry Smith   PetscInt     pad;               /* padding */
154ce68768SBarry Smith   PetscBLASInt *pivots;           /* pivots in LU factorization */
16284134d9SBarry Smith   PetscBLASInt lda;               /* Lapack leading dimension of data */
1721a2c019SBarry Smith   PetscTruth   changelda;         /* change lda on resize? Default unless user set lda */
1821a2c019SBarry Smith   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
19273d9f13SBarry Smith   PetscTruth   user_alloc;        /* true if the user provided the dense data */
20f6f390a9SLois Curfman McInnes } Mat_SeqDense;
21f6f390a9SLois Curfman McInnes 
22dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqDense(Mat A,Vec,Vec);
23dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec,Vec,Vec);
24dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec,Vec);
25dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec,Vec,Vec);
26a9fe9ddaSSatish Balay EXTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
27a9fe9ddaSSatish Balay EXTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
28a9fe9ddaSSatish Balay EXTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
29a9fe9ddaSSatish Balay EXTERN PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
30a9fe9ddaSSatish Balay EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
31a9fe9ddaSSatish Balay EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
32*4ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
33*4ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
34*4ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
3518476b09SLois Curfman McInnes 
36f6f390a9SLois Curfman McInnes #endif
37