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