xref: /petsc/src/mat/impls/dense/seq/dense.h (revision 6947451f174971b5cbba4c8ab49dd392302c23da)
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 */
25*6947451fSStefano Zampini   /* Support for MatDenseGetColumnVec */
26*6947451fSStefano Zampini   Vec               cvec;      /* vector representation of a given column */
27*6947451fSStefano Zampini   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
28*6947451fSStefano Zampini   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
29afea5027SHong Zhang 
304222ddf1SHong Zhang   Mat_MatTransMatMult *atb;       /* used by MatTransposeMatMultxxx_SeqAIJ_SeqDense */
31f6f390a9SLois Curfman McInnes } Mat_SeqDense;
32f6f390a9SLois Curfman McInnes 
334222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
34ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
354222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
36ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
374222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
38ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
394222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
40ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
414222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
42ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat,Mat,Mat);
434222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
44c2916339SPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat,Mat,Mat);
4552c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat,Mat,PetscReal,Mat*);
4652c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat,Mat,Mat);
47ca15aa20SStefano Zampini 
484222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat);
494222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat);
504222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat);
514222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat);
524222ddf1SHong Zhang 
53ca15aa20SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);
54ca15aa20SStefano Zampini 
55ca15aa20SStefano Zampini /* Used by SeqDenseCUDA */
56ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption);
57ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*);
58ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
59ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat);
60ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]);
61ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]);
62ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure);
63ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
64ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
65ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
66ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
67ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*);
68ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*);
69ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*);
70ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
71ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
72ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*);
73ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool);
74637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat,Vec,PetscInt);
75637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat,PetscScalar);
76*6947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat,PetscInt,Vec*);
77*6947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat,PetscInt,Vec*);
78*6947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
79*6947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
80*6947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
81*6947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
8218476b09SLois Curfman McInnes 
832bf066beSStefano Zampini #if defined(PETSC_HAVE_CUDA)
842bf066beSStefano Zampini PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat);
852bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*);
862bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*);
872bf066beSStefano Zampini #endif
882bf066beSStefano Zampini 
89afea5027SHong Zhang PETSC_INTERN PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat);
9005709791SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
91ec01deb9SMatthew Knepley 
92d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
93d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
94d528f656SJakub Kruzik 
958491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat,PetscViewer);
968491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat,PetscViewer);
97eb91f321SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer);
98f6f390a9SLois Curfman McInnes #endif
99