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