Actual source code: dense.h
petsc-3.8.4 2018-03-24
4: #include <petsc/private/matimpl.h>
5: #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */
7: /*
8: MATSEQDENSE format - conventional dense Fortran storage (by columns)
9: */
11: typedef struct {
12: PetscScalar *v; /* matrix elements */
13: PetscScalar *unplacedarray; /* if one called MatDensePlaceArray(), this is where it stashed the original */
14: PetscBool roworiented; /* if true, row oriented input (default) */
15: PetscInt pad; /* padding */
16: PetscBLASInt *pivots; /* pivots in LU factorization */
17: PetscBLASInt lfwork; /* length of work array in factorization */
18: PetscScalar *fwork; /* work array in factorization */
19: PetscBLASInt lda; /* Lapack leading dimension of data */
20: PetscBool changelda; /* change lda on resize? Default unless user set lda */
21: PetscBLASInt Mmax,Nmax; /* indicates the largest dimensions of data possible */
22: PetscBool user_alloc; /* true if the user provided the dense data */
23: PetscBool unplaced_user_alloc;
24: Mat ptapwork; /* workspace (SeqDense matrix) for PtAP */
26: Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult_SeqAIJ_SeqDense */
27: } Mat_SeqDense;
29: extern PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
30: extern PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
31: extern PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
32: extern PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
33: extern PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
34: extern PetscErrorCode MatMatTransposeMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
35: extern PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
36: extern PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
37: extern PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
38: extern PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
40: PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
42: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
43: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
44: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
46: #endif