1*47d993e7Ssuyashtn /* Portions of this code are under: 2*47d993e7Ssuyashtn Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved. 3*47d993e7Ssuyashtn */ 46524c165SJacob Faibussowitsch #ifndef __DENSE_H 5f6f390a9SLois Curfman McInnes #define __DENSE_H 6af0996ceSBarry Smith #include <petsc/private/matimpl.h> 76718818eSStefano Zampini /* TODO REMOVE */ 8afea5027SHong Zhang #include <../src/mat/impls/aij/seq/aij.h> /* Mat_MatTransMatMult is defined here */ 9f6f390a9SLois Curfman McInnes 10f6f390a9SLois Curfman McInnes /* 11f6f390a9SLois Curfman McInnes MATSEQDENSE format - conventional dense Fortran storage (by columns) 12f6f390a9SLois Curfman McInnes */ 13f6f390a9SLois Curfman McInnes 14f6f390a9SLois Curfman McInnes typedef struct { 15ea709b57SSatish Balay PetscScalar *v; /* matrix elements */ 16d3042a70SBarry Smith PetscScalar *unplacedarray; /* if one called MatDensePlaceArray(), this is where it stashed the original */ 17ace3abfcSBarry Smith PetscBool roworiented; /* if true, row oriented input (default) */ 1813f74950SBarry Smith PetscInt pad; /* padding */ 194ce68768SBarry Smith PetscBLASInt *pivots; /* pivots in LU factorization */ 20a49dc2a2SStefano Zampini PetscBLASInt lfwork; /* length of work array in factorization */ 21a49dc2a2SStefano Zampini PetscScalar *fwork; /* work array in factorization */ 224905a7bcSToby Isaac PetscScalar *tau; /* scalar factors of QR factorization */ 234396437dSToby Isaac Vec qrrhs; /* RHS for solving with QR (solution vector can't hold copy of RHS) */ 24284134d9SBarry Smith PetscBLASInt lda; /* Lapack leading dimension of data */ 254905a7bcSToby Isaac PetscBLASInt rank; /* numerical rank (of a QR factorized matrix) */ 26ace3abfcSBarry Smith PetscBool user_alloc; /* true if the user provided the dense data */ 27d3042a70SBarry Smith PetscBool unplaced_user_alloc; 28abc3b08eSStefano Zampini Mat ptapwork; /* workspace (SeqDense matrix) for PtAP */ 296718818eSStefano Zampini 305ea7661aSPierre Jolivet /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */ 315ea7661aSPierre Jolivet Mat cmat; /* matrix representation of a given subset of columns */ 326947451fSStefano Zampini Vec cvec; /* vector representation of a given column */ 336947451fSStefano Zampini const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */ 346947451fSStefano Zampini PetscInt vecinuse; /* if cvec is in use (col = vecinuse-1) */ 355ea7661aSPierre Jolivet PetscInt matinuse; /* if cmat is in use (cbegin = matinuse-1) */ 36f6f390a9SLois Curfman McInnes } Mat_SeqDense; 37f6f390a9SLois Curfman McInnes 384222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat, Mat, PetscReal, Mat); 39ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat, Mat, Mat); 404222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat, Mat, PetscReal, Mat); 41ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat, Mat, Mat); 424222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat, Mat, PetscReal, Mat); 43ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat, Mat, Mat); 444222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat, Mat, PetscReal, Mat); 45ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat, Mat, Mat); 464222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqBAIJ_SeqDense(Mat, Mat, PetscReal, Mat); 47ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqBAIJ_SeqDense(Mat, Mat, Mat); 484222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_SeqSBAIJ_SeqDense(Mat, Mat, PetscReal, Mat); 49c2916339SPierre Jolivet PETSC_INTERN PetscErrorCode MatMatMultNumeric_SeqSBAIJ_SeqDense(Mat, Mat, Mat); 5052c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat, Mat, PetscReal, Mat *); 5152c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat, Mat, Mat); 52ca15aa20SStefano Zampini 534222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat); 544222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat); 554222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat); 564222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat); 574222ddf1SHong Zhang 58ca15aa20SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat); 59ca15aa20SStefano Zampini 60*47d993e7Ssuyashtn /* Used by SeqDenseCUDA and seqDenseHIP */ 61ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat, Mat, MatDuplicateOption); 62ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatNorm_SeqDense(Mat, NormType, PetscReal *); 63ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatView_SeqDense(Mat, PetscViewer); 64ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDestroy_SeqDense(Mat); 65ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetArray_SeqDense(Mat, PetscScalar *[]); 66ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreArray_SeqDense(Mat, PetscScalar *[]); 67ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatAXPY_SeqDense(Mat, PetscScalar, Mat, MatStructure); 68ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTransposeAdd_SeqDense(Mat, Vec, Vec, Vec); 69ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat, Vec, Vec); 70ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat, Vec, Vec, Vec); 71ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat, Vec, Vec); 72ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatDuplicate_SeqDense(Mat, MatDuplicateOption, Mat *); 73ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat, PetscScalar *); 74ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactor_SeqDense(Mat, IS, const MatFactorInfo *); 75ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactor_SeqDense(Mat, IS, IS, const MatFactorInfo *); 76bf5a80bcSToby Isaac PETSC_INTERN PetscErrorCode MatQRFactor_SeqDense(Mat, IS, const MatFactorInfo *); 77ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat, Mat, IS, const MatFactorInfo *); 78ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat, Mat, IS, IS, const MatFactorInfo *); 79bf5a80bcSToby Isaac PETSC_INTERN PetscErrorCode MatQRFactorSymbolic_SeqDense(Mat, Mat, IS, const MatFactorInfo *); 80ca15aa20SStefano Zampini PETSC_INTERN PetscErrorCode MatSeqDenseSymmetrize_Private(Mat, PetscBool); 81637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatGetColumnVector_SeqDense(Mat, Vec, PetscInt); 82637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat, PetscScalar); 83126e4e93SJose E. Roman PETSC_INTERN PetscErrorCode MatShift_SeqDense(Mat, PetscScalar); 846947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat, PetscInt, Vec *); 856947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat, PetscInt, Vec *); 866947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat, PetscInt, Vec *); 876947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat, PetscInt, Vec *); 886947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat, PetscInt, Vec *); 896947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat, PetscInt, Vec *); 90a2748737SPierre Jolivet PETSC_INTERN PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat, PetscInt, PetscInt, PetscInt, PetscInt, Mat *); 915ea7661aSPierre Jolivet PETSC_INTERN PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat, Mat *); 9217359960SJose E. Roman PETSC_INTERN PetscErrorCode MatDenseSetLDA_SeqDense(Mat, PetscInt); 93a76f77c3SStefano Zampini PETSC_INTERN PetscErrorCode MatCopy_SeqDense(Mat, Mat, MatStructure); 943d8925e7SStefano Zampini PETSC_INTERN PetscErrorCode MatZeroEntries_SeqDense(Mat); 9575f6d85dSStefano Zampini PETSC_INTERN PetscErrorCode MatSetUp_SeqDense(Mat); 963faff063SStefano Zampini PETSC_INTERN PetscErrorCode MatSetRandom_SeqDense(Mat, PetscRandom); 9718476b09SLois Curfman McInnes 982bf066beSStefano Zampini #if defined(PETSC_HAVE_CUDA) 99126e4e93SJose E. Roman PETSC_INTERN PetscErrorCode MatShift_DenseCUDA_Private(PetscScalar *, PetscScalar, PetscInt, PetscInt, PetscInt, PetscInt); 1002bf066beSStefano Zampini PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat); 1012bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDenseCUDA_SeqDense(Mat, MatType, MatReuse, Mat *); 1022bf066beSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat, MatType, MatReuse, Mat *); 1032bf066beSStefano Zampini #endif 1042bf066beSStefano Zampini 105*47d993e7Ssuyashtn #if defined(PETSC_HAVE_HIP) 106*47d993e7Ssuyashtn PETSC_INTERN PetscErrorCode MatShift_DenseHIP_Private(PetscScalar *, PetscScalar, PetscInt, PetscInt, PetscInt, PetscInt); 107*47d993e7Ssuyashtn PETSC_EXTERN PetscErrorCode MatSeqDenseHIPInvertFactors_Private(Mat); 108*47d993e7Ssuyashtn PETSC_INTERN PetscErrorCode MatConvert_SeqDenseHIP_SeqDense(Mat, MatType, MatReuse, Mat *); 109*47d993e7Ssuyashtn PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqDenseHIP(Mat, MatType, MatReuse, Mat *); 110*47d993e7Ssuyashtn #endif 111*47d993e7Ssuyashtn 11205709791SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat); 113ec01deb9SMatthew Knepley 114d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 115d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 116d528f656SJakub Kruzik 1178491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat, PetscViewer); 1188491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat, PetscViewer); 119eb91f321SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat, PetscViewer); 120f6f390a9SLois Curfman McInnes #endif 121