xref: /petsc/src/mat/impls/dense/seq/dense.h (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
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 lfwork;            /* length of work array in factorization */
17   PetscScalar  *fwork;            /* work array in factorization */
18   PetscBLASInt lda;               /* Lapack leading dimension of data */
19   PetscBool    changelda;         /* change lda on resize? Default unless user set lda */
20   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
21   PetscBool    user_alloc;        /* true if the user provided the dense data */
22   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */
23 
24   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMult_SeqAIJ_SeqDense */
25 } Mat_SeqDense;
26 
27 extern PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
28 extern PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
29 extern PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
30 extern PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
31 extern PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
32 extern PetscErrorCode MatMatTransposeMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
33 extern PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
34 extern PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
35 extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
36 extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
37 
38 PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
39 
40 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
41 PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
42 PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
43 
44 #endif
45