1 2 #if !defined(__DENSE_H) 3 #define __DENSE_H 4 #include <petsc/private/matimpl.h> 5 /* TODO REMOVE */ 6 #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */ 7 8 /* 9 MATSEQDENSE format - conventional dense Fortran storage (by columns) 10 */ 11 12 typedef struct { 13 PetscScalar *v; /* matrix elements */ 14 PetscScalar *unplacedarray; /* if one called MatDensePlaceArray(), this is where it stashed the original */ 15 PetscBool roworiented; /* if true, row oriented input (default) */ 16 PetscInt pad; /* padding */ 17 PetscBLASInt *pivots; /* pivots in LU factorization */ 18 PetscBLASInt lfwork; /* length of work array in factorization */ 19 PetscScalar *fwork; /* work array in factorization */ 20 PetscBLASInt lda; /* Lapack leading dimension of data */ 21 PetscBool user_alloc; /* true if the user provided the dense data */ 22 PetscBool unplaced_user_alloc; 23 Mat ptapwork; /* workspace (SeqDense matrix) for PtAP */ 24 25 /* Support for MatDenseGetColumnVec */ 26 Vec cvec; /* vector representation of a given column */ 27 const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */ 28 PetscInt vecinuse; /* if cvec is in use (col = vecinuse-1) */ 29 } Mat_SeqDense; 30 31 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat); 32 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat); 33 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat); 34 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat); 35 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat,Mat,PetscReal,Mat); 36 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat,Mat,Mat); 37 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat); 38 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat); 39 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat,Mat,PetscReal,Mat); 40 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat,Mat,Mat); 41 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat,Mat,PetscReal,Mat); 42 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat,Mat,Mat); 43 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat,Mat,PetscReal,Mat*); 44 PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat,Mat,Mat); 45 46 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat); 47 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat); 48 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat); 49 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat); 50 51 PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat); 52 53 /* Used by SeqDenseCUDA */ 54 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat,Mat,MatDuplicateOption); 55 PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat,NormType,PetscReal*); 56 PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat,PetscViewer); 57 PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat); 58 PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat,PetscScalar*[]); 59 PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat,PetscScalar*[]); 60 PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat,PetscScalar,Mat,MatStructure); 61 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat,Vec,Vec,Vec); 62 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec); 63 PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec); 64 PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec); 65 PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat,MatDuplicateOption,Mat*); 66 PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat,PetscScalar*); 67 PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat,IS,const MatFactorInfo*); 68 PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*); 69 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat,Mat,IS,const MatFactorInfo*); 70 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat,Mat,IS,IS,const MatFactorInfo*); 71 PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat,PetscBool); 72 PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat,Vec,PetscInt); 73 PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat,PetscScalar); 74 PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat,PetscInt,Vec*); 75 PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat,PetscInt,Vec*); 76 PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat,PetscInt,Vec*); 77 PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat,PetscInt,Vec*); 78 PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat,PetscInt,Vec*); 79 PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat,PetscInt,Vec*); 80 81 #if defined(PETSC_HAVE_CUDA) 82 PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat); 83 PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat,MatType,MatReuse,Mat*); 84 PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat,MatType,MatReuse,Mat*); 85 #endif 86 87 PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat); 88 89 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 90 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 91 92 PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat,PetscViewer); 93 PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat,PetscViewer); 94 PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat,PetscViewer); 95 #endif 96