xref: /petsc/src/mat/impls/dense/seq/dense.h (revision d3042a70f471d64055a67fbedc03cbed6f0fa57d)
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 */
13*d3042a70SBarry Smith   PetscScalar  *unplacedarray;    /* if one called MatDensePlaceArray(), this is where it stashed the original */
14ace3abfcSBarry Smith   PetscBool    roworiented;       /* if true, row oriented input (default) */
1513f74950SBarry Smith   PetscInt     pad;               /* padding */
164ce68768SBarry Smith   PetscBLASInt *pivots;           /* pivots in LU factorization */
17284134d9SBarry Smith   PetscBLASInt lda;               /* Lapack leading dimension of data */
18ace3abfcSBarry Smith   PetscBool    changelda;         /* change lda on resize? Default unless user set lda */
1921a2c019SBarry Smith   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
20ace3abfcSBarry Smith   PetscBool    user_alloc;        /* true if the user provided the dense data */
21*d3042a70SBarry Smith   PetscBool    unplaced_user_alloc;
22abc3b08eSStefano Zampini   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */
23afea5027SHong Zhang 
24afea5027SHong Zhang   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMult_SeqAIJ_SeqDense */
25f6f390a9SLois Curfman McInnes } Mat_SeqDense;
26f6f390a9SLois Curfman McInnes 
2709573ac7SBarry Smith extern PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
2809573ac7SBarry Smith extern PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
2975648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
3075648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
3175648e8dSHong Zhang extern PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
3209573ac7SBarry Smith extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
3309573ac7SBarry Smith extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
3418476b09SLois Curfman McInnes 
35afea5027SHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
36afea5027SHong Zhang 
37150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
38150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
39ec01deb9SMatthew Knepley 
40f6f390a9SLois Curfman McInnes #endif
41