Actual source code: matseqdensecuda.cu
1: #include "../matseqdensecupm.hpp"
3: using namespace Petsc::mat::cupm;
4: using Petsc::device::cupm::DeviceType;
6: static constexpr impl::MatDense_Seq_CUPM<DeviceType::CUDA> cupm_mat{};
8: /*MC
9: MATSEQDENSECUDA - "seqdensecuda" - A matrix type to be used for sequential dense matrices on
10: GPUs.
12: Options Database Keys:
13: . -mat_type seqdensecuda - sets the matrix type to `MATSEQDENSECUDA` during a call to `MatSetFromOptions()`
15: Level: beginner
17: .seealso: `Mat`, `MatType`, `MATSEQDENSE`
18: M*/
19: PETSC_INTERN PetscErrorCode MatCreate_SeqDenseCUDA(Mat A)
20: {
21: PetscFunctionBegin;
22: PetscCall(cupm_mat.Create(A));
23: PetscFunctionReturn(PETSC_SUCCESS);
24: }
26: PETSC_INTERN PetscErrorCode MatSolverTypeRegister_DENSECUDA(void)
27: {
28: PetscFunctionBegin;
29: PetscCall(impl::MatSolverTypeRegister_DENSECUPM<DeviceType::CUDA>());
30: PetscFunctionReturn(PETSC_SUCCESS);
31: }
33: PetscErrorCode MatConvert_SeqDense_SeqDenseCUDA(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
34: {
35: PetscFunctionBegin;
36: PetscCall(cupm_mat.Convert_SeqDense_SeqDenseCUPM(A, newtype, reuse, newmat));
37: PetscFunctionReturn(PETSC_SUCCESS);
38: }
40: PetscErrorCode MatMatMultNumeric_SeqDenseCUDA_SeqDenseCUDA_Internal(Mat A, Mat B, Mat C, PetscBool TA, PetscBool TB)
41: {
42: PetscFunctionBegin;
43: PetscCall(impl::MatMatMultNumeric_SeqDenseCUPM_SeqDenseCUPM<DeviceType::CUDA>(A, B, C, TA, TB));
44: PetscFunctionReturn(PETSC_SUCCESS);
45: }
47: PetscErrorCode MatSeqDenseCUDAInvertFactors_Internal(Mat A)
48: {
49: PetscFunctionBegin;
50: PetscCall(cupm_mat.InvertFactors(A));
51: PetscFunctionReturn(PETSC_SUCCESS);
52: }
54: /*@C
55: MatCreateSeqDenseCUDA - Creates a sequential matrix in dense format using CUDA.
57: Collective
59: Input Parameters:
60: + comm - MPI communicator
61: . m - number of rows
62: . n - number of columns
63: - data - optional location of GPU matrix data. Pass `NULL` to let PETSc to control matrix
64: memory allocation
66: Output Parameter:
67: . A - the matrix
69: Level: intermediate
71: .seealso: `Mat`, `MatType`, `MATSEQDENSECUDA`, `MatCreate()`, `MatCreateSeqDense()`
72: @*/
73: PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm comm, PetscInt m, PetscInt n, PetscScalar *data, Mat *A)
74: {
75: PetscFunctionBegin;
76: PetscCall(MatCreateSeqDenseCUPM<DeviceType::CUDA>(comm, m, n, data, A));
77: PetscFunctionReturn(PETSC_SUCCESS);
78: }