xref: /petsc/src/mat/impls/dense/seq/dense.h (revision 4396437d7f67f9fe41f15190ecbfe4af75f5cbcd)
1f6f390a9SLois Curfman McInnes 
2f6f390a9SLois Curfman McInnes #if !defined(__DENSE_H)
3f6f390a9SLois Curfman McInnes #define __DENSE_H
4af0996ceSBarry Smith #include <petsc/private/matimpl.h>
56718818eSStefano Zampini /* TODO REMOVE */
6afea5027SHong Zhang #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */
7f6f390a9SLois Curfman McInnes 
8f6f390a9SLois Curfman McInnes /*
9f6f390a9SLois Curfman McInnes   MATSEQDENSE format - conventional dense Fortran storage (by columns)
10f6f390a9SLois Curfman McInnes */
11f6f390a9SLois Curfman McInnes 
12f6f390a9SLois Curfman McInnes typedef struct {
13ea709b57SSatish Balay   PetscScalar  *v;                /* matrix elements */
14d3042a70SBarry Smith   PetscScalar  *unplacedarray;    /* if one called MatDensePlaceArray(), this is where it stashed the original */
15ace3abfcSBarry Smith   PetscBool    roworiented;       /* if true, row oriented input (default) */
1613f74950SBarry Smith   PetscInt     pad;               /* padding */
174ce68768SBarry Smith   PetscBLASInt *pivots;           /* pivots in LU factorization */
18a49dc2a2SStefano Zampini   PetscBLASInt lfwork;            /* length of work array in factorization */
19a49dc2a2SStefano Zampini   PetscScalar  *fwork;            /* work array in factorization */
204905a7bcSToby Isaac   PetscScalar  *tau;              /* scalar factors of QR factorization */
21*4396437dSToby Isaac   Vec          qrrhs;            /* RHS for solving with QR (solution vector can't hold copy of RHS) */
22284134d9SBarry Smith   PetscBLASInt lda;               /* Lapack leading dimension of data */
234905a7bcSToby Isaac   PetscBLASInt rank;              /* numerical rank (of a QR factorized matrix) */
24ace3abfcSBarry Smith   PetscBool    user_alloc;        /* true if the user provided the dense data */
25d3042a70SBarry Smith   PetscBool    unplaced_user_alloc;
26abc3b08eSStefano Zampini   Mat          ptapwork;          /* workspace (SeqDense matrix) for PtAP */
276718818eSStefano Zampini 
285ea7661aSPierre Jolivet   /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
295ea7661aSPierre Jolivet   Mat               cmat;      /* matrix representation of a given subset of columns */
306947451fSStefano Zampini   Vec               cvec;      /* vector representation of a given column */
316947451fSStefano Zampini   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
326947451fSStefano Zampini   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
335ea7661aSPierre Jolivet   PetscInt          matinuse;  /* if cmat is in use (cbegin = matinuse-1) */
34f6f390a9SLois Curfman McInnes } Mat_SeqDense;
35f6f390a9SLois Curfman McInnes 
364222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
37ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
384222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
39ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
404222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat);
41ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat);
424222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
43ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
444222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
45ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat,Mat,Mat);
464222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
47c2916339SPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat,Mat,Mat);
4852c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat,Mat,PetscReal,Mat*);
4952c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat,Mat,Mat);
50ca15aa20SStefano Zampini 
514222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat);
524222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat);
534222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat);
544222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat);
554222ddf1SHong Zhang 
56ca15aa20SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat);
57ca15aa20SStefano Zampini 
58ca15aa20SStefano Zampini /* Used by SeqDenseCUDA */
59ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption);
60ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*);
61ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer);
62ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat);
63ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]);
64ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]);
65ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure);
66ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec);
67ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec);
68ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec);
69ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec);
70ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*);
71ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*);
72ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*);
73ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*);
74ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*);
75ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*);
76ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool);
77637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat,Vec,PetscInt);
78637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat,PetscScalar);
796947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat,PetscInt,Vec*);
806947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat,PetscInt,Vec*);
816947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
826947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat,PetscInt,Vec*);
836947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
846947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat,PetscInt,Vec*);
855ea7661aSPierre Jolivet PETSC_INTERN PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat,PetscInt,PetscInt,Mat*);
865ea7661aSPierre Jolivet PETSC_INTERN PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat,Mat*);
8717359960SJose E. Roman PETSC_INTERN PetscErrorCode MatDenseSetLDA_SeqDense(Mat,PetscInt);
88a76f77c3SStefano Zampini PETSC_INTERN PetscErrorCode MatCopy_SeqDense(Mat,Mat,MatStructure);
8918476b09SLois Curfman McInnes 
902bf066beSStefano Zampini #if defined(PETSC_HAVE_CUDA)
912bf066beSStefano Zampini PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat);
922bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*);
932bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*);
942bf066beSStefano Zampini #endif
952bf066beSStefano Zampini 
9605709791SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
97ec01deb9SMatthew Knepley 
98d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
99d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
100d528f656SJakub Kruzik 
1018491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat,PetscViewer);
1028491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat,PetscViewer);
103eb91f321SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer);
104f6f390a9SLois Curfman McInnes #endif
105