Actual source code: petscmat_kokkos.hpp

  1: #pragma once

  3: #include <petscconf.h>

  5: /* SUBMANSEC = Mat */

  7: #if defined(PETSC_HAVE_KOKKOS)

  9:   #include <Kokkos_Core.hpp>
 10: #include <petscmat.h>

 12: /*@C
 13:    MatCreateSeqAIJKokkosWithKokkosViews - Creates a MATSEQAIJKOKKOS matrix with Kokkos views of the aij data

 15:    Synopsis:
 16: #include <petscmat_kokkos.hpp>
 17:    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);

 19:    Logically Collective, No Fortran Support

 21:    Input Parameter:
 22: +  comm  - the MPI communicator
 23: -  m     - row size
 24: -  n     - the column size
 25: -  i     - the Kokkos view of row data (in Kokkos::DefaultExecutionSpace)
 26: -  j     - the Kokkos view of the column data (in Kokkos::DefaultExecutionSpace)
 27: -  a     - the Kokkos view of the values (in Kokkos::DefaultExecutionSpace)

 29:    Output Parameter:
 30: .  A  - the `MATSEQAIJKOKKOS` matrix

 32:    Level: intermediate

 34:    Notes:
 35:    Creates a Mat given the csr data input as Kokkos views. This routine allows a Mat
 36:    to be built without involving the host. Don't modify entries in the views after this routine.
 37:    There should be no outstanding asynchronous operations on the views (ie this routine does not call fence()
 38:    before using the views)

 40: .seealso:
 41: @*/
 42: PetscErrorCode MatCreateSeqAIJKokkosWithKokkosViews(MPI_Comm, PetscInt, PetscInt, Kokkos::View<PetscInt *> &, Kokkos::View<PetscInt *> &, Kokkos::View<PetscScalar *> &, Mat *);

 44: #endif