xref: /petsc/src/mat/impls/dense/seq/dense.h (revision abc3b08ea96a1ec1fce68b6dd23e12ba18b3f617)
1f6f390a9SLois Curfman McInnes 
2f6f390a9SLois Curfman McInnes #if !defined(__DENSE_H)
3f6f390a9SLois Curfman McInnes #define __DENSE_H
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
5afea5027SHong Zhang #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */
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 */
13ace3abfcSBarry Smith   PetscBool    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 */
17ace3abfcSBarry Smith   PetscBool    changelda;         /* change lda on resize? Default unless user set lda */
1821a2c019SBarry Smith   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
19ace3abfcSBarry Smith   PetscBool    user_alloc;        /* true if the user provided the dense data */
20*abc3b08eSStefano Zampini   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */
21afea5027SHong Zhang 
22afea5027SHong Zhang   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMult_SeqAIJ_SeqDense */
23f6f390a9SLois Curfman McInnes } Mat_SeqDense;
24f6f390a9SLois Curfman McInnes 
2509573ac7SBarry Smith extern PetscErrorCode MatMult_SeqDense(Mat A,Vec,Vec);
2609573ac7SBarry Smith extern PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec,Vec,Vec);
2709573ac7SBarry Smith extern PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec,Vec);
2809573ac7SBarry Smith extern PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec,Vec,Vec);
2909573ac7SBarry Smith extern PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
3009573ac7SBarry Smith extern PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
3175648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
3275648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
3375648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
3409573ac7SBarry Smith extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
3509573ac7SBarry Smith extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
3618476b09SLois Curfman McInnes 
37afea5027SHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
38afea5027SHong Zhang 
398cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
408cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
41ec01deb9SMatthew Knepley 
42f6f390a9SLois Curfman McInnes #endif
43