xref: /petsc/src/mat/impls/dense/seq/dense.h (revision 71383f7570b59a18517995bd8cb04ed04c3be442)
147d993e7Ssuyashtn /* Portions of this code are under:
247d993e7Ssuyashtn    Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved.
347d993e7Ssuyashtn */
4a4963045SJacob Faibussowitsch #pragma once
5af0996ceSBarry Smith #include <petsc/private/matimpl.h>
66718818eSStefano Zampini /* TODO REMOVE */
7cc1eb50dSBarry Smith #include <../src/mat/impls/aij/seq/aij.h> /* MatProductCtx_MatTransMatMult is defined here */
8f6f390a9SLois Curfman McInnes 
9f6f390a9SLois Curfman McInnes /*
10f6f390a9SLois Curfman McInnes   MATSEQDENSE format - conventional dense Fortran storage (by columns)
11f6f390a9SLois Curfman McInnes */
12f6f390a9SLois Curfman McInnes 
13f6f390a9SLois Curfman McInnes typedef struct {
14ea709b57SSatish Balay   PetscScalar  *v;             /* matrix elements */
15d3042a70SBarry Smith   PetscScalar  *unplacedarray; /* if one called MatDensePlaceArray(), this is where it stashed the original */
165e116b59SBarry Smith   PetscBool     roworiented;   /* if true, row-oriented input (default) */
1713f74950SBarry Smith   PetscInt      pad;           /* padding */
184ce68768SBarry Smith   PetscBLASInt *pivots;        /* pivots in LU factorization */
19a49dc2a2SStefano Zampini   PetscBLASInt  lfwork;        /* length of work array in factorization */
20a49dc2a2SStefano Zampini   PetscScalar  *fwork;         /* work array in factorization */
214905a7bcSToby Isaac   PetscScalar  *tau;           /* scalar factors of QR factorization */
224396437dSToby Isaac   Vec           qrrhs;         /* RHS for solving with QR (solution vector can't hold copy of RHS) */
23284134d9SBarry Smith   PetscBLASInt  lda;           /* Lapack leading dimension of data */
244905a7bcSToby Isaac   PetscBLASInt  rank;          /* numerical rank (of a QR factorized matrix) */
25ace3abfcSBarry Smith   PetscBool     user_alloc;    /* true if the user provided the dense data */
26d3042a70SBarry Smith   PetscBool     unplaced_user_alloc;
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);
46*1766d9c3SPierre Jolivet PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_SeqBAIJ_SeqDense(Mat, Mat, PetscReal, Mat);
47*1766d9c3SPierre Jolivet PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_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);
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 
564742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatCreate_SeqDense(Mat);
57ca15aa20SStefano Zampini 
5847d993e7Ssuyashtn /* Used by SeqDenseCUDA and seqDenseHIP */
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);
66459e8d23SBlanca Mellado Pinto PETSC_INTERN PetscErrorCode MatMultHermitianTransposeAdd_SeqDense(Mat, Vec, Vec, Vec);
67459e8d23SBlanca Mellado Pinto PETSC_INTERN PetscErrorCode MatMultHermitianTranspose_SeqDense(Mat, Vec, Vec);
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);
823853def2SToby Isaac PETSC_INTERN PetscErrorCode MatConjugate_SeqDense(Mat);
83637a0070SStefano Zampini PETSC_INTERN PetscErrorCode MatScale_SeqDense(Mat, PetscScalar);
84126e4e93SJose E. Roman PETSC_INTERN PetscErrorCode MatShift_SeqDense(Mat, PetscScalar);
856947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat, PetscInt, Vec *);
866947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat, PetscInt, Vec *);
876947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat, PetscInt, Vec *);
886947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat, PetscInt, Vec *);
896947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat, PetscInt, Vec *);
906947451fSStefano Zampini PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat, PetscInt, Vec *);
91a2748737SPierre Jolivet PETSC_INTERN PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat, PetscInt, PetscInt, PetscInt, PetscInt, Mat *);
925ea7661aSPierre Jolivet PETSC_INTERN PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat, Mat *);
9317359960SJose E. Roman PETSC_INTERN PetscErrorCode MatDenseSetLDA_SeqDense(Mat, PetscInt);
94a76f77c3SStefano Zampini PETSC_INTERN PetscErrorCode MatCopy_SeqDense(Mat, Mat, MatStructure);
953d8925e7SStefano Zampini PETSC_INTERN PetscErrorCode MatZeroEntries_SeqDense(Mat);
9675f6d85dSStefano Zampini PETSC_INTERN PetscErrorCode MatSetUp_SeqDense(Mat);
973faff063SStefano Zampini PETSC_INTERN PetscErrorCode MatSetRandom_SeqDense(Mat, PetscRandom);
9814277c92SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetDiagonal_SeqDense(Mat, Vec);
9918476b09SLois Curfman McInnes 
100d016bddeSToby Isaac PETSC_INTERN PetscErrorCode MatMultAddColumnRangeKernel_SeqDense(Mat, Vec, Vec, Vec, PetscInt, PetscInt, PetscBool, PetscBool);
101d016bddeSToby Isaac PETSC_INTERN PetscErrorCode MatMultColumnRangeKernel_SeqDense(Mat, Vec, Vec, PetscInt, PetscInt, PetscBool, PetscBool);
102d016bddeSToby Isaac PETSC_INTERN PetscErrorCode MatMultColumnRange_SeqDense(Mat, Vec, Vec, PetscInt, PetscInt);
1030be0d8bdSHansol Suh PETSC_INTERN PetscErrorCode MatMultAddColumnRange_SeqDense(Mat, Vec, Vec, Vec, PetscInt, PetscInt);
1040be0d8bdSHansol Suh PETSC_INTERN PetscErrorCode MatMultHermitianTransposeColumnRange_SeqDense(Mat, Vec, Vec, PetscInt, PetscInt);
1050be0d8bdSHansol Suh PETSC_INTERN PetscErrorCode MatMultHermitianTransposeAddColumnRange_SeqDense(Mat, Vec, Vec, Vec, PetscInt, PetscInt);
1060be0d8bdSHansol Suh 
1072bf066beSStefano Zampini #if defined(PETSC_HAVE_CUDA)
1084742e46bSJacob Faibussowitsch PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Internal(Mat);
1094742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode                MatMatMultNumeric_SeqDenseCUDA_SeqDenseCUDA_Internal(Mat, Mat, Mat, PetscBool, PetscBool);
1102bf066beSStefano Zampini PETSC_INTERN PetscErrorCode                MatConvert_SeqDense_SeqDenseCUDA(Mat, MatType, MatReuse, Mat *);
1112bf066beSStefano Zampini #endif
1122bf066beSStefano Zampini 
11347d993e7Ssuyashtn #if defined(PETSC_HAVE_HIP)
1144742e46bSJacob Faibussowitsch PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqDenseHIPInvertFactors_Internal(Mat);
1154742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode                MatMatMultNumeric_SeqDenseHIP_SeqDenseHIP_Internal(Mat, Mat, Mat, PetscBool, PetscBool);
11647d993e7Ssuyashtn PETSC_INTERN PetscErrorCode                MatConvert_SeqDense_SeqDenseHIP(Mat, MatType, MatReuse, Mat *);
11747d993e7Ssuyashtn #endif
11847d993e7Ssuyashtn 
119d6acfc2dSPierre Jolivet PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
120ec01deb9SMatthew Knepley 
121d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
122d528f656SJakub Kruzik PETSC_INTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIDense(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
123d528f656SJakub Kruzik 
1248491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatView_Dense_Binary(Mat, PetscViewer);
1258491ab44SLisandro Dalcin PETSC_INTERN PetscErrorCode MatLoad_Dense_Binary(Mat, PetscViewer);
126eb91f321SVaclav Hapla PETSC_INTERN PetscErrorCode MatLoad_Dense_HDF5(Mat, PetscViewer);
127d16ceb75SStefano Zampini 
128d16ceb75SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseCreateColumnVec_Private(Mat, Vec *);
129