xref: /petsc/include/petscmat_kokkos.hpp (revision 23386071a3c0abcd54fc2a6f38cff7dacf103a83)
1*c0c276a7Ssdargavi #pragma once
2*c0c276a7Ssdargavi 
3*c0c276a7Ssdargavi #include <petscconf.h>
4*c0c276a7Ssdargavi 
5*c0c276a7Ssdargavi /* SUBMANSEC = Mat */
6*c0c276a7Ssdargavi 
7*c0c276a7Ssdargavi #if defined(PETSC_HAVE_KOKKOS)
8*c0c276a7Ssdargavi 
9*c0c276a7Ssdargavi   #include <Kokkos_Core.hpp>
10*c0c276a7Ssdargavi   #include <petscmat.h>
11*c0c276a7Ssdargavi 
12*c0c276a7Ssdargavi /*@C
13*c0c276a7Ssdargavi    MatCreateSeqAIJKokkosWithKokkosViews - Creates a MATSEQAIJKOKKOS matrix with Kokkos views of the aij data
14*c0c276a7Ssdargavi 
15*c0c276a7Ssdargavi    Synopsis:
16*c0c276a7Ssdargavi    #include <petscmat_kokkos.hpp>
17*c0c276a7Ssdargavi    PetscErrorCode MatCreateSeqAIJKokkosWithKokkosViews  (MPI_Comm comm, PetscInt m, PetscInt n, Kokkos::View<PetscInt *>& i_d, Kokkos::View<PetscInt *>& j_d, Kokkos::View<PetscScalar *>& a_d, Mat *A);
18*c0c276a7Ssdargavi 
19*c0c276a7Ssdargavi    Logically Collective, No Fortran Support
20*c0c276a7Ssdargavi 
21*c0c276a7Ssdargavi    Input Parameter:
22*c0c276a7Ssdargavi +  comm  - the MPI communicator
23*c0c276a7Ssdargavi -  m     - row size
24*c0c276a7Ssdargavi -  n     - the column size
25*c0c276a7Ssdargavi -  i     - the Kokkos view of row data (in Kokkos::DefaultExecutionSpace)
26*c0c276a7Ssdargavi -  j     - the Kokkos view of the column data (in Kokkos::DefaultExecutionSpace)
27*c0c276a7Ssdargavi -  a     - the Kokkos view of the values (in Kokkos::DefaultExecutionSpace)
28*c0c276a7Ssdargavi 
29*c0c276a7Ssdargavi    Output Parameter:
30*c0c276a7Ssdargavi .  A  - the `MATSEQAIJKOKKOS` matrix
31*c0c276a7Ssdargavi 
32*c0c276a7Ssdargavi    Level: intermediate
33*c0c276a7Ssdargavi 
34*c0c276a7Ssdargavi    Notes:
35*c0c276a7Ssdargavi    Creates a Mat given the csr data input as Kokkos views. This routine allows a Mat
36*c0c276a7Ssdargavi    to be built without involving the host. Don't modify entries in the views after this routine.
37*c0c276a7Ssdargavi    There should be no outstanding asynchronous operations on the views (ie this routine does not call fence()
38*c0c276a7Ssdargavi    before using the views)
39*c0c276a7Ssdargavi 
40*c0c276a7Ssdargavi .seealso:
41*c0c276a7Ssdargavi @*/
42*c0c276a7Ssdargavi PetscErrorCode MatCreateSeqAIJKokkosWithKokkosViews(MPI_Comm, PetscInt, PetscInt, Kokkos::View<PetscInt *> &, Kokkos::View<PetscInt *> &, Kokkos::View<PetscScalar *> &, Mat *);
43*c0c276a7Ssdargavi 
44*c0c276a7Ssdargavi #endif
45