1 /* Portions of this code are under: 2 Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved. 3 */ 4 #ifndef __DENSE_H 5 #define __DENSE_H 6 #include <petsc/private/matimpl.h> 7 /* TODO REMOVE */ 8 #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */ 9 10 /* 11 MATSEQDENSE format - conventional dense Fortran storage (by columns) 12 */ 13 14 typedef struct { 15 PetscScalar *v; /* matrix elements */ 16 PetscScalar *unplacedarray; /* if one called MatDensePlaceArray(), this is where it stashed the original */ 17 PetscBool roworiented; /* if true, row oriented input (default) */ 18 PetscInt pad; /* padding */ 19 PetscBLASInt *pivots; /* pivots in LU factorization */ 20 PetscBLASInt lfwork; /* length of work array in factorization */ 21 PetscScalar *fwork; /* work array in factorization */ 22 PetscScalar *tau; /* scalar factors of QR factorization */ 23 Vec qrrhs; /* RHS for solving with QR (solution vector can't hold copy of RHS) */ 24 PetscBLASInt lda; /* Lapack leading dimension of data */ 25 PetscBLASInt rank; /* numerical rank (of a QR factorized matrix) */ 26 PetscBool user_alloc; /* true if the user provided the dense data */ 27 PetscBool unplaced_user_alloc; 28 29 /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */ 30 Mat cmat; /* matrix representation of a given subset of columns */ 31 Vec cvec; /* vector representation of a given column */ 32 const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */ 33 PetscInt vecinuse; /* if cvec is in use (col = vecinuse-1) */ 34 PetscInt matinuse; /* if cmat is in use (cbegin = matinuse-1) */ 35 } Mat_SeqDense; 36 37 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat, Mat, PetscReal, Mat); 38 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat, Mat, Mat); 39 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat, Mat, PetscReal, Mat); 40 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat, Mat, Mat); 41 PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat, Mat, PetscReal, Mat); 42 PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat, Mat, Mat); 43 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat, Mat, PetscReal, Mat); 44 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat, Mat, Mat); 45 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat, Mat, PetscReal, Mat); 46 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat, Mat, Mat); 47 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat, Mat, PetscReal, Mat); 48 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat, Mat, Mat); 49 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat, Mat, PetscReal, Mat *); 50 PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat, Mat, Mat); 51 52 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat); 53 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat); 54 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat); 55 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat); 56 57 PETSC_INTERN PetscErrorCode MatCreate_SeqDense(Mat); 58 59 /* Used by SeqDenseCUDA and seqDenseHIP */ 60 PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat, Mat, MatDuplicateOption); 61 PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat, NormType, PetscReal *); 62 PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat, PetscViewer); 63 PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat); 64 PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat, PetscScalar *[]); 65 PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat, PetscScalar *[]); 66 PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat, PetscScalar, Mat, MatStructure); 67 PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat, Vec, Vec, Vec); 68 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat, Vec, Vec); 69 PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat, Vec, Vec, Vec); 70 PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat, Vec, Vec); 71 PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat, MatDuplicateOption, Mat *); 72 PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat, PetscScalar *); 73 PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat, IS, const MatFactorInfo *); 74 PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat, IS, IS, const MatFactorInfo *); 75 PETSC_INTERN PetscErrorCode MatQRFactor_SeqDense(Mat, IS, const MatFactorInfo *); 76 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat, Mat, IS, const MatFactorInfo *); 77 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat, Mat, IS, IS, const MatFactorInfo *); 78 PETSC_INTERN PetscErrorCode MatQRFactorSymbolic_SeqDense(Mat, Mat, IS, const MatFactorInfo *); 79 PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat, PetscBool); 80 PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat, Vec, PetscInt); 81 PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat, PetscScalar); 82 PETSC_INTERN PetscErrorCode MatShift_SeqDense(Mat, PetscScalar); 83 PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat, PetscInt, Vec *); 84 PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat, PetscInt, Vec *); 85 PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat, PetscInt, Vec *); 86 PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat, PetscInt, Vec *); 87 PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat, PetscInt, Vec *); 88 PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat, PetscInt, Vec *); 89 PETSC_INTERN PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat, PetscInt, PetscInt, PetscInt, PetscInt, Mat *); 90 PETSC_INTERN PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat, Mat *); 91 PETSC_INTERN PetscErrorCode MatDenseSetLDA_SeqDense(Mat, PetscInt); 92 PETSC_INTERN PetscErrorCode MatCopy_SeqDense(Mat, Mat, MatStructure); 93 PETSC_INTERN PetscErrorCode MatZeroEntries_SeqDense(Mat); 94 PETSC_INTERN PetscErrorCode MatSetUp_SeqDense(Mat); 95 PETSC_INTERN PetscErrorCode MatSetRandom_SeqDense(Mat, PetscRandom); 96 97 #if defined(PETSC_HAVE_CUDA) 98 PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Internal(Mat); 99 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDenseCUDA_SeqDenseCUDA_Internal(Mat, Mat, Mat, PetscBool, PetscBool); 100 PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat, MatType, MatReuse, Mat *); 101 #endif 102 103 #if defined(PETSC_HAVE_HIP) 104 PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqDenseHIPInvertFactors_Internal(Mat); 105 PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDenseHIP_SeqDenseHIP_Internal(Mat, Mat, Mat, PetscBool, PetscBool); 106 PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseHIP(Mat, MatType, MatReuse, Mat *); 107 #endif 108 109 PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat); 110 111 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 112 PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 113 114 PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat, PetscViewer); 115 PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat, PetscViewer); 116 PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat, PetscViewer); 117 #endif 118