xref: /petsc/src/ksp/pc/impls/pbjacobi/kokkos/pbjacobi_kok.kokkos.cxx (revision ffeef943c8ee50edff320d8a3135bb0c94853e4c)
1 #include <petscvec_kokkos.hpp>
2 #include <../src/vec/vec/impls/seq/kokkos/veckokkosimpl.hpp>
3 #include <petscdevice.h>
4 #include <../src/ksp/pc/impls/pbjacobi/pbjacobi.h>
5 
6 struct PC_PBJacobi_Kokkos {
7   PetscScalarKokkosDualView diag_dual;
8 
9   PC_PBJacobi_Kokkos(PetscInt len, PetscScalar *diag_ptr_h)
10   {
11     PetscScalarKokkosViewHost diag_h(diag_ptr_h, len);
12     auto                      diag_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), diag_h);
13     diag_dual                        = PetscScalarKokkosDualView(diag_d, diag_h);
14   }
15 
16   PetscErrorCode Update(const PetscScalar *diag_ptr_h)
17   {
18     PetscFunctionBegin;
19     PetscCheck(diag_dual.view_host().data() == diag_ptr_h, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Host pointer has changed since last call");
20     PetscCallCXX(diag_dual.modify_host()); /* mark the host has newer data */
21     PetscCallCXX(diag_dual.sync_device());
22     PetscFunctionReturn(PETSC_SUCCESS);
23   }
24 };
25 
26 /* Make 'transpose' a template parameter instead of a function input parameter, so that
27  it will be a const in template instantiation and gets optimized out.
28 */
29 template <PetscBool transpose>
30 static PetscErrorCode PCApplyOrTranspose_PBJacobi_Kokkos(PC pc, Vec x, Vec y)
31 {
32   PC_PBJacobi               *jac   = (PC_PBJacobi *)pc->data;
33   PC_PBJacobi_Kokkos        *pckok = static_cast<PC_PBJacobi_Kokkos *>(jac->spptr);
34   ConstPetscScalarKokkosView xv;
35   PetscScalarKokkosView      yv;
36   PetscScalarKokkosView      Av = pckok->diag_dual.view_device();
37   const PetscInt             bs = jac->bs, mbs = jac->mbs, bs2 = bs * bs;
38   const char                *label = transpose ? "PCApplyTranspose_PBJacobi_Kokkos" : "PCApply_PBJacobi_Kokkos";
39 
40   PetscFunctionBegin;
41   PetscCall(PetscLogGpuTimeBegin());
42   VecErrorIfNotKokkos(x);
43   VecErrorIfNotKokkos(y);
44   PetscCall(VecGetKokkosView(x, &xv));
45   PetscCall(VecGetKokkosViewWrite(y, &yv));
46   PetscCallCXX(Kokkos::parallel_for(
47     label, bs * mbs, KOKKOS_LAMBDA(PetscInt row) {
48       const PetscScalar *Ap, *xp;
49       PetscScalar       *yp;
50       PetscInt           i, j, k;
51 
52       k  = row / bs;                                /* k-th block */
53       i  = row % bs;                                /* this thread deals with i-th row of the block */
54       Ap = &Av(bs2 * k + i * (transpose ? bs : 1)); /* Ap points to the first entry of i-th row */
55       xp = &xv(bs * k);
56       yp = &yv(bs * k);
57       /* multiply i-th row (column) with x */
58       yp[i] = 0.0;
59       for (j = 0; j < bs; j++) {
60         yp[i] += Ap[0] * xp[j];
61         Ap += (transpose ? 1 : bs); /* block is in column major order */
62       }
63     }));
64   PetscCall(VecRestoreKokkosView(x, &xv));
65   PetscCall(VecRestoreKokkosViewWrite(y, &yv));
66   PetscCall(PetscLogGpuFlops(bs * bs * mbs * 2)); /* FMA on entries in all blocks */
67   PetscCall(PetscLogGpuTimeEnd());
68   PetscFunctionReturn(PETSC_SUCCESS);
69 }
70 
71 static PetscErrorCode PCDestroy_PBJacobi_Kokkos(PC pc)
72 {
73   PC_PBJacobi *jac = (PC_PBJacobi *)pc->data;
74 
75   PetscFunctionBegin;
76   PetscCallCXX(delete static_cast<PC_PBJacobi_Kokkos *>(jac->spptr));
77   PetscCall(PCDestroy_PBJacobi(pc));
78   PetscFunctionReturn(PETSC_SUCCESS);
79 }
80 
81 PETSC_INTERN PetscErrorCode PCSetUp_PBJacobi_Kokkos(PC pc)
82 {
83   PC_PBJacobi *jac = (PC_PBJacobi *)pc->data;
84   PetscInt     len;
85 
86   PetscFunctionBegin;
87   PetscCall(PCSetUp_PBJacobi_Host(pc)); /* Compute the inverse on host now. Might worth doing it on device directly */
88   len = jac->bs * jac->bs * jac->mbs;
89   if (!jac->spptr) {
90     PetscCallCXX(jac->spptr = new PC_PBJacobi_Kokkos(len, const_cast<PetscScalar *>(jac->diag)));
91   } else {
92     PC_PBJacobi_Kokkos *pckok = static_cast<PC_PBJacobi_Kokkos *>(jac->spptr);
93     PetscCall(pckok->Update(jac->diag));
94   }
95   PetscCall(PetscLogCpuToGpu(sizeof(PetscScalar) * len));
96 
97   pc->ops->apply          = PCApplyOrTranspose_PBJacobi_Kokkos<PETSC_FALSE>;
98   pc->ops->applytranspose = PCApplyOrTranspose_PBJacobi_Kokkos<PETSC_TRUE>;
99   pc->ops->destroy        = PCDestroy_PBJacobi_Kokkos;
100   PetscFunctionReturn(PETSC_SUCCESS);
101 }
102