xref: /petsc/src/mat/impls/dense/seq/dense.h (revision 1ebf93c6b7d760d592de6ebe6cdc0debaa3caf75) !
1 
2 #if !defined(__DENSE_H)
3 #define __DENSE_H
4 #include <petsc/private/matimpl.h>
5 #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */
6 
7 /*
8   MATSEQDENSE format - conventional dense Fortran storage (by columns)
9 */
10 
11 typedef struct {
12   PetscScalar  *v;                /* matrix elements */
13   PetscBool    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   PetscBool    changelda;         /* change lda on resize? Default unless user set lda */
18   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
19   PetscBool    user_alloc;        /* true if the user provided the dense data */
20   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */
21 
22   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMult_SeqAIJ_SeqDense */
23 } Mat_SeqDense;
24 
25 extern PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
26 extern PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
27 extern PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
28 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
29 extern PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
30 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
31 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
32 
33 PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
34 
35 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
36 PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
37 
38 #endif
39