xref: /petsc/src/mat/impls/dense/seq/dense.h (revision c2916339d3c5437df97526bab4e3cf11b3acd91c)
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 */
13d3042a70SBarry 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 */
17a49dc2a2SStefano Zampini   PetscBLASInt lfwork;            /* length of work array in factorization */
18a49dc2a2SStefano Zampini   PetscScalar  *fwork;            /* work array in factorization */
19284134d9SBarry Smith   PetscBLASInt lda;               /* Lapack leading dimension of data */
20ace3abfcSBarry Smith   PetscBool    changelda;         /* change lda on resize? Default unless user set lda */
2121a2c019SBarry Smith   PetscBLASInt Mmax,Nmax;         /* indicates the largest dimensions of data possible */
22ace3abfcSBarry Smith   PetscBool    user_alloc;        /* true if the user provided the dense data */
23d3042a70SBarry Smith   PetscBool    unplaced_user_alloc;
24abc3b08eSStefano Zampini   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */
25afea5027SHong Zhang 
26afea5027SHong Zhang   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMult_SeqAIJ_SeqDense */
27f6f390a9SLois Curfman McInnes } Mat_SeqDense;
28f6f390a9SLois Curfman McInnes 
29ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
30ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
31ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
32ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
33ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
34ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatTransposeMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
35ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat*);
36ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
37ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
38ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
39ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
40ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat,Mat,Mat);
41*c2916339SPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat,Mat,PetscReal,Mat*);
42*c2916339SPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat,Mat,Mat);
4352c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat,Mat,PetscReal,Mat*);
4452c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat,Mat,Mat);
45ca15aa20SStefano Zampini 
46ca15aa20SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);
47ca15aa20SStefano Zampini 
48ca15aa20SStefano Zampini /* Used by SeqDenseCUDA */
49ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption);
50ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*);
51ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
52ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat);
53ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]);
54ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]);
55ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure);
56ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
57ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
58ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
59ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
60ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*);
61ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*);
62ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*);
63ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
64ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
65ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*);
66ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool);
6718476b09SLois Curfman McInnes 
682bf066beSStefano Zampini #if defined(PETSC_HAVE_CUDA)
692bf066beSStefano Zampini PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat);
702bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*);
712bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*);
722bf066beSStefano Zampini #endif
732bf066beSStefano Zampini 
74afea5027SHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
75afea5027SHong Zhang 
76150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
77a001520aSPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMult_SeqBAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
78*c2916339SPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMult_SeqSBAIJ_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
7952c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMult_Nest_Dense(Mat,Mat,MatReuse,PetscReal,Mat*);
80150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat,Mat,MatReuse,PetscReal,Mat*);
8105709791SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
82ec01deb9SMatthew Knepley 
83d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
84d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
85d528f656SJakub Kruzik 
86eb91f321SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer);
87f6f390a9SLois Curfman McInnes #endif
88