xref: /petsc/src/ts/utils/dmplexlandau/kokkos/landau.kokkos.cxx (revision f13dfd9ea68e0ddeee984e65c377a1819eab8a8a)
1 /*
2   Implements the Kokkos kernel
3 */
4 #include <petscconf.h>
5 // On Sunspot, as of 05/04/2023, the default SYCL compiler oneapi/eng-compiler/2022.12.30.003
6 // failed to link this file, so we add "defined(PETSC_HAVE_SYCL)" to disable landau with SYCL
7 // to facilitate petsc users' experiments on Sunspot.
8 // TODO: remove it when Intel fixed the bug. See MR !6412
9 #if defined(PETSC_HAVE_CUDA_CLANG) || defined(PETSC_HAVE_SYCL)
10   #include <petsclandau.h>
11   #define LANDAU_NOT_IMPLEMENTED SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not supported with CLANG or SYCL")
12 PetscErrorCode LandauKokkosJacobian(DM[], const PetscInt, const PetscInt, const PetscInt, const PetscInt, const PetscInt[], PetscReal[], PetscScalar[], const PetscScalar[], const LandauStaticData *, const PetscReal, const PetscLogEvent[], const PetscInt[], const PetscInt[], Mat[], Mat)
13 {
14   LANDAU_NOT_IMPLEMENTED;
15 }
16 PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps *, pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE], PetscInt[], PetscInt)
17 {
18   LANDAU_NOT_IMPLEMENTED;
19 }
20 PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps *, PetscInt)
21 {
22   LANDAU_NOT_IMPLEMENTED;
23 }
24 PetscErrorCode LandauKokkosStaticDataSet(DM, const PetscInt, const PetscInt, const PetscInt, const PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], PetscReal[], LandauStaticData *)
25 {
26   LANDAU_NOT_IMPLEMENTED;
27 }
28 PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData *)
29 {
30   LANDAU_NOT_IMPLEMENTED;
31 }
32 #else
33   #include <petscvec_kokkos.hpp>
34   #include <petsclandau.h>
35   #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
36   #include <petsc/private/deviceimpl.h>
37   #include <petscts.h>
38 
39   #include <Kokkos_Core.hpp>
40   #include <cstdio>
41 typedef Kokkos::TeamPolicy<>::member_type team_member;
42   #include "../land_tensors.h"
43 
44 namespace landau_inner_red
45 { // namespace helps with name resolution in reduction identity
46 template <class ScalarType>
47 struct array_type {
48   ScalarType gg2[LANDAU_DIM];
49   ScalarType gg3[LANDAU_DIM][LANDAU_DIM];
50 
51   KOKKOS_INLINE_FUNCTION // Default constructor - Initialize to 0's
52   array_type()
53   {
54     for (int j = 0; j < LANDAU_DIM; j++) {
55       gg2[j] = 0;
56       for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] = 0;
57     }
58   }
59   KOKKOS_INLINE_FUNCTION // Copy Constructor
60   array_type(const array_type &rhs)
61   {
62     for (int j = 0; j < LANDAU_DIM; j++) {
63       gg2[j] = rhs.gg2[j];
64       for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] = rhs.gg3[j][k];
65     }
66   }
67   KOKKOS_INLINE_FUNCTION // add operator
68     array_type &
69     operator+=(const array_type &src)
70   {
71     for (int j = 0; j < LANDAU_DIM; j++) {
72       gg2[j] += src.gg2[j];
73       for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] += src.gg3[j][k];
74     }
75     return *this;
76   }
77   KOKKOS_INLINE_FUNCTION // volatile add operator
78     void
79     operator+=(const volatile array_type &src) volatile
80   {
81     for (int j = 0; j < LANDAU_DIM; j++) {
82       gg2[j] += src.gg2[j];
83       for (int k = 0; k < LANDAU_DIM; k++) gg3[j][k] += src.gg3[j][k];
84     }
85   }
86 };
87 typedef array_type<PetscReal> TensorValueType; // used to simplify code below
88 } // namespace landau_inner_red
89 
90 namespace Kokkos
91 { //reduction identity must be defined in Kokkos namespace
92 template <>
93 struct reduction_identity<landau_inner_red::TensorValueType> {
94   KOKKOS_FORCEINLINE_FUNCTION static landau_inner_red::TensorValueType sum() { return landau_inner_red::TensorValueType(); }
95 };
96 } // namespace Kokkos
97 
98 extern "C" {
99 PetscErrorCode LandauKokkosCreateMatMaps(P4estVertexMaps maps[], pointInterpolationP4est (*pointMaps)[LANDAU_MAX_Q_FACE], PetscInt Nf[], PetscInt grid)
100 {
101   P4estVertexMaps                                                                                                                                       h_maps; /* host container */
102   const Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>>     h_points((pointInterpolationP4est *)pointMaps, maps[grid].num_reduced);
103   const Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQND], Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_gidx((LandauIdx *)maps[grid].gIdx, maps[grid].num_elements);
104   Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight>     *d_points = new Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight>("points", maps[grid].num_reduced);
105   Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQND], Kokkos::LayoutRight> *d_gidx   = new Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQND], Kokkos::LayoutRight>("gIdx", maps[grid].num_elements);
106 
107   PetscFunctionBegin;
108   Kokkos::deep_copy(*d_gidx, h_gidx);
109   Kokkos::deep_copy(*d_points, h_points);
110   h_maps.num_elements = maps[grid].num_elements;
111   h_maps.num_face     = maps[grid].num_face;
112   h_maps.num_reduced  = maps[grid].num_reduced;
113   h_maps.deviceType   = maps[grid].deviceType;
114   h_maps.numgrids     = maps[grid].numgrids;
115   h_maps.Nf           = Nf[grid];
116   h_maps.c_maps       = (pointInterpolationP4est(*)[LANDAU_MAX_Q_FACE])d_points->data();
117   maps[grid].vp1      = (void *)d_points;
118   h_maps.gIdx         = (LandauIdx(*)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQND])d_gidx->data();
119   maps[grid].vp2      = (void *)d_gidx;
120   {
121     Kokkos::View<P4estVertexMaps, Kokkos::HostSpace> h_maps_k(&h_maps);
122     Kokkos::View<P4estVertexMaps>                   *d_maps_k = new Kokkos::View<P4estVertexMaps>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(), h_maps_k));
123     Kokkos::deep_copy(*d_maps_k, h_maps_k);
124     maps[grid].d_self = d_maps_k->data();
125     maps[grid].vp3    = (void *)d_maps_k;
126   }
127   PetscFunctionReturn(PETSC_SUCCESS);
128 }
129 PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps maps[], PetscInt num_grids)
130 {
131   PetscFunctionBegin;
132   for (PetscInt grid = 0; grid < num_grids; grid++) {
133     auto a = static_cast<Kokkos::View<pointInterpolationP4est *[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight> *>(maps[grid].vp1);
134     auto b = static_cast<Kokkos::View<LandauIdx *[LANDAU_MAX_SPECIES][LANDAU_MAX_NQND], Kokkos::LayoutRight> *>(maps[grid].vp2);
135     auto c = static_cast<Kokkos::View<P4estVertexMaps *> *>(maps[grid].vp3);
136     delete a;
137     delete b;
138     delete c;
139   }
140   PetscFunctionReturn(PETSC_SUCCESS);
141 }
142 
143 PetscErrorCode LandauKokkosStaticDataSet(DM plex, const PetscInt Nq, const PetscInt Nb, const PetscInt batch_sz, const PetscInt num_grids, PetscInt a_numCells[], PetscInt a_species_offset[], PetscInt a_mat_offset[], PetscReal a_nu_alpha[], PetscReal a_nu_beta[], PetscReal a_invMass[], PetscReal a_lambdas[], PetscReal a_invJ[], PetscReal a_x[], PetscReal a_y[], PetscReal a_z[], PetscReal a_w[], LandauStaticData *SData_d)
144 {
145   PetscReal       *BB, *DD;
146   PetscTabulation *Tf;
147   PetscInt         dim;
148   PetscInt         ip_offset[LANDAU_MAX_GRIDS + 1], ipf_offset[LANDAU_MAX_GRIDS + 1], elem_offset[LANDAU_MAX_GRIDS + 1], nip, IPf_sz, Nftot;
149   PetscDS          prob;
150 
151   PetscFunctionBegin;
152   PetscCall(DMGetDimension(plex, &dim));
153   PetscCall(DMGetDS(plex, &prob));
154   PetscCheck(LANDAU_DIM == dim, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM);
155   PetscCall(PetscDSGetTabulation(prob, &Tf));
156   BB           = Tf[0]->T[0];
157   DD           = Tf[0]->T[1];
158   ip_offset[0] = ipf_offset[0] = elem_offset[0] = 0;
159   nip                                           = 0;
160   IPf_sz                                        = 0;
161   for (PetscInt grid = 0; grid < num_grids; grid++) {
162     PetscInt nfloc        = a_species_offset[grid + 1] - a_species_offset[grid];
163     elem_offset[grid + 1] = elem_offset[grid] + a_numCells[grid];
164     nip += a_numCells[grid] * Nq;
165     ip_offset[grid + 1] = nip;
166     IPf_sz += Nq * nfloc * a_numCells[grid];
167     ipf_offset[grid + 1] = IPf_sz;
168   }
169   Nftot = a_species_offset[num_grids];
170   PetscCall(PetscKokkosInitializeCheck());
171   {
172     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_alpha(a_nu_alpha, Nftot);
173     auto                                                                                                            alpha = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("alpha", Nftot);
174     SData_d->alpha                                                                                                        = static_cast<void *>(alpha);
175     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_beta(a_nu_beta, Nftot);
176     auto                                                                                                            beta = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("beta", Nftot);
177     SData_d->beta                                                                                                        = static_cast<void *>(beta);
178     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_invMass(a_invMass, Nftot);
179     auto                                                                                                            invMass = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("invMass", Nftot);
180     SData_d->invMass                                                                                                        = static_cast<void *>(invMass);
181     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_lambdas(a_lambdas, LANDAU_MAX_GRIDS * LANDAU_MAX_GRIDS);
182     auto                                                                                                            lambdas = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("lambdas", LANDAU_MAX_GRIDS * LANDAU_MAX_GRIDS);
183     SData_d->lambdas                                                                                                        = static_cast<void *>(lambdas);
184     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_BB(BB, Nq * Nb);
185     auto                                                                                                            B = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("B", Nq * Nb);
186     SData_d->B                                                                                                        = static_cast<void *>(B);
187     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_DD(DD, Nq * Nb * dim);
188     auto                                                                                                            D = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("D", Nq * Nb * dim);
189     SData_d->D                                                                                                        = static_cast<void *>(D);
190     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_invJ(a_invJ, nip * dim * dim);
191     auto                                                                                                            invJ = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("invJ", nip * dim * dim);
192     SData_d->invJ                                                                                                        = static_cast<void *>(invJ);
193     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_x(a_x, nip);
194     auto                                                                                                            x = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("x", nip);
195     SData_d->x                                                                                                        = static_cast<void *>(x);
196     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_y(a_y, nip);
197     auto                                                                                                            y = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("y", nip);
198     SData_d->y                                                                                                        = static_cast<void *>(y);
199     const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_w(a_w, nip);
200     auto                                                                                                            w = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("w", nip);
201     SData_d->w                                                                                                        = static_cast<void *>(w);
202 
203     Kokkos::deep_copy(*alpha, h_alpha);
204     Kokkos::deep_copy(*beta, h_beta);
205     Kokkos::deep_copy(*invMass, h_invMass);
206     Kokkos::deep_copy(*lambdas, h_lambdas);
207     Kokkos::deep_copy(*B, h_BB);
208     Kokkos::deep_copy(*D, h_DD);
209     Kokkos::deep_copy(*invJ, h_invJ);
210     Kokkos::deep_copy(*x, h_x);
211     Kokkos::deep_copy(*y, h_y);
212     Kokkos::deep_copy(*w, h_w);
213 
214     if (dim == 3) {
215       const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_z(a_z, nip);
216       auto                                                                                                            z = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("z", nip);
217       SData_d->z                                                                                                        = static_cast<void *>(z);
218       Kokkos::deep_copy(*z, h_z);
219     } else SData_d->z = NULL;
220 
221     //
222     const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_NCells(a_numCells, num_grids);
223     auto                                                                                                           nc = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("NCells", num_grids);
224     SData_d->NCells                                                                                                   = static_cast<void *>(nc);
225     Kokkos::deep_copy(*nc, h_NCells);
226 
227     const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_species_offset(a_species_offset, num_grids + 1);
228     auto                                                                                                           soff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("species_offset", num_grids + 1);
229     SData_d->species_offset                                                                                             = static_cast<void *>(soff);
230     Kokkos::deep_copy(*soff, h_species_offset);
231 
232     const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_mat_offset(a_mat_offset, num_grids + 1);
233     auto                                                                                                           moff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("mat_offset", num_grids + 1);
234     SData_d->mat_offset                                                                                                 = static_cast<void *>(moff);
235     Kokkos::deep_copy(*moff, h_mat_offset);
236 
237     const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_ip_offset(ip_offset, num_grids + 1);
238     auto                                                                                                           ipoff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("ip_offset", num_grids + 1);
239     SData_d->ip_offset                                                                                                   = static_cast<void *>(ipoff);
240     Kokkos::deep_copy(*ipoff, h_ip_offset);
241 
242     const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_elem_offset(elem_offset, num_grids + 1);
243     auto                                                                                                           eoff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("elem_offset", num_grids + 1);
244     SData_d->elem_offset                                                                                                = static_cast<void *>(eoff);
245     Kokkos::deep_copy(*eoff, h_elem_offset);
246 
247     const Kokkos::View<PetscInt *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_ipf_offset(ipf_offset, num_grids + 1);
248     auto                                                                                                           ipfoff = new Kokkos::View<PetscInt *, Kokkos::LayoutLeft>("ipf_offset", num_grids + 1);
249     SData_d->ipf_offset                                                                                                   = static_cast<void *>(ipfoff);
250     Kokkos::deep_copy(*ipfoff, h_ipf_offset);
251   #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
252     auto ipfdf_data = new Kokkos::View<PetscReal ***, Kokkos::LayoutLeft>("fdf", batch_sz, dim + 1, IPf_sz);
253   #else
254     auto ipfdf_data = new Kokkos::View<PetscReal ***, Kokkos::LayoutRight>("fdf", batch_sz, dim + 1, IPf_sz);
255   #endif
256     SData_d->ipfdf_data = static_cast<void *>(ipfdf_data);
257     auto Eq_m           = new Kokkos::View<PetscReal *, Kokkos::LayoutLeft>("Eq_m", Nftot); // allocate but do not set, same for whole batch
258     SData_d->Eq_m       = static_cast<void *>(Eq_m);
259     const Kokkos::View<LandauIdx *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_coo_elem_offsets((LandauIdx *)SData_d->coo_elem_offsets, SData_d->coo_n_cellsTot + 1);
260     auto                                                                                                            coo_elem_offsets = new Kokkos::View<LandauIdx *, Kokkos::LayoutLeft>("coo_elem_offsets", SData_d->coo_n_cellsTot + 1);
261     const Kokkos::View<LandauIdx *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_coo_elem_fullNb((LandauIdx *)SData_d->coo_elem_fullNb, SData_d->coo_n_cellsTot);
262     auto                                                                                                            coo_elem_fullNb = new Kokkos::View<LandauIdx *, Kokkos::LayoutLeft>("coo_elem_offsets", SData_d->coo_n_cellsTot);
263     const Kokkos::View<LandauIdx *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_coo_elem_point_offsets((LandauIdx *)SData_d->coo_elem_point_offsets, SData_d->coo_n_cellsTot * (LANDAU_MAX_NQND + 1));
264     auto coo_elem_point_offsets = new Kokkos::View<LandauIdx *, Kokkos::LayoutLeft>("coo_elem_point_offsets", SData_d->coo_n_cellsTot * (LANDAU_MAX_NQND + 1));
265     // assign & copy
266     Kokkos::deep_copy(*coo_elem_offsets, h_coo_elem_offsets);
267     Kokkos::deep_copy(*coo_elem_fullNb, h_coo_elem_fullNb);
268     Kokkos::deep_copy(*coo_elem_point_offsets, h_coo_elem_point_offsets);
269     // need to free this now and use pointer space
270     PetscCall(PetscFree3(SData_d->coo_elem_offsets, SData_d->coo_elem_fullNb, SData_d->coo_elem_point_offsets));
271     SData_d->coo_elem_offsets       = static_cast<void *>(coo_elem_offsets);
272     SData_d->coo_elem_fullNb        = static_cast<void *>(coo_elem_fullNb);
273     SData_d->coo_elem_point_offsets = static_cast<void *>(coo_elem_point_offsets);
274     auto coo_vals                   = new Kokkos::View<PetscScalar *, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace>("coo_vals", SData_d->coo_size);
275     SData_d->coo_vals               = static_cast<void *>(coo_vals);
276   }
277   SData_d->maps = NULL; // not used
278   PetscFunctionReturn(PETSC_SUCCESS);
279 }
280 
281 PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData *SData_d)
282 {
283   PetscFunctionBegin;
284   if (SData_d->alpha) {
285     auto alpha = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->alpha);
286     delete alpha;
287     SData_d->alpha = NULL;
288     auto beta      = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->beta);
289     delete beta;
290     auto invMass = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invMass);
291     delete invMass;
292     auto lambdas = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->lambdas);
293     delete lambdas;
294     auto B = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->B);
295     delete B;
296     auto D = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->D);
297     delete D;
298     auto invJ = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invJ);
299     delete invJ;
300     auto x = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->x);
301     delete x;
302     auto y = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->y);
303     delete y;
304     if (SData_d->z) {
305       auto z = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->z);
306       delete z;
307     }
308   #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
309     auto z = static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutLeft> *>(SData_d->ipfdf_data);
310   #else
311     auto z = static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutRight> *>(SData_d->ipfdf_data);
312   #endif
313     delete z;
314     auto w = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->w);
315     delete w;
316     auto Eq_m = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->Eq_m);
317     delete Eq_m;
318     // offset
319     auto nc = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->NCells);
320     delete nc;
321     auto soff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->species_offset);
322     delete soff;
323     auto moff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->mat_offset);
324     delete moff;
325     auto ipoff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ip_offset);
326     delete ipoff;
327     auto eoff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->elem_offset);
328     delete eoff;
329     auto ipfoff = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ipf_offset);
330     delete ipfoff;
331     auto coo_elem_offsets = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>((void *)SData_d->coo_elem_offsets);
332     delete coo_elem_offsets;
333     auto coo_elem_fullNb = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>((void *)SData_d->coo_elem_fullNb);
334     delete coo_elem_fullNb;
335     auto coo_elem_point_offsets = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>((void *)SData_d->coo_elem_point_offsets);
336     delete coo_elem_point_offsets;
337     SData_d->coo_elem_offsets       = NULL;
338     SData_d->coo_elem_point_offsets = NULL;
339     SData_d->coo_elem_fullNb        = NULL;
340     auto coo_vals                   = static_cast<Kokkos::View<PetscScalar *, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace> *>((void *)SData_d->coo_vals);
341     delete coo_vals;
342   }
343   PetscFunctionReturn(PETSC_SUCCESS);
344 }
345 
346   #define KOKKOS_SHARED_LEVEL 0   // 0 is shared, 1 is global
347 
348 KOKKOS_INLINE_FUNCTION
349 PetscErrorCode landau_mat_assemble(PetscScalar *coo_vals, const PetscScalar Aij, const PetscInt f, const PetscInt g, const PetscInt Nb, PetscInt moffset, const PetscInt elem, const PetscInt fieldA, const P4estVertexMaps *d_maps, const LandauIdx coo_elem_offsets[], const LandauIdx coo_elem_fullNb[], const LandauIdx (*coo_elem_point_offsets)[LANDAU_MAX_NQND + 1], const PetscInt glb_elem_idx, const PetscInt bid_coo_sz_batch)
350 {
351   PetscInt               idx, q, nr, nc;
352   const LandauIdx *const Idxs                         = &d_maps->gIdx[elem][fieldA][0];
353   PetscScalar            row_scale[LANDAU_MAX_Q_FACE] = {0}, col_scale[LANDAU_MAX_Q_FACE] = {0};
354   { // mirror (i,j) in CreateStaticGPUData
355     const int fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb;
356     idx = Idxs[f];
357     if (idx >= 0) {
358       nr           = 1;
359       row_scale[0] = 1.;
360     } else {
361       idx = -idx - 1;
362       for (q = 0, nr = 0; q < d_maps->num_face; q++, nr++) {
363         if (d_maps->c_maps[idx][q].gid < 0) break;
364         row_scale[q] = d_maps->c_maps[idx][q].scale;
365       }
366     }
367     idx = Idxs[g];
368     if (idx >= 0) {
369       nc           = 1;
370       col_scale[0] = 1.;
371     } else {
372       idx = -idx - 1;
373       nc  = d_maps->num_face;
374       for (q = 0, nc = 0; q < d_maps->num_face; q++, nc++) {
375         if (d_maps->c_maps[idx][q].gid < 0) break;
376         col_scale[q] = d_maps->c_maps[idx][q].scale;
377       }
378     }
379     const int idx0 = bid_coo_sz_batch + coo_elem_offsets[glb_elem_idx] + fieldA * fullNb2 + fullNb * coo_elem_point_offsets[glb_elem_idx][f] + nr * coo_elem_point_offsets[glb_elem_idx][g];
380     for (int q = 0, idx2 = idx0; q < nr; q++) {
381       for (int d = 0; d < nc; d++, idx2++) coo_vals[idx2] = row_scale[q] * col_scale[d] * Aij;
382     }
383   }
384   return PETSC_SUCCESS;
385 }
386 
387 PetscErrorCode LandauKokkosJacobian(DM plex[], const PetscInt Nq, const PetscInt Nb, const PetscInt batch_sz, const PetscInt num_grids, const PetscInt a_numCells[], PetscReal a_Eq_m[], PetscScalar a_elem_closure[], const PetscScalar a_xarray[], const LandauStaticData *SData_d, const PetscReal shift, const PetscLogEvent events[], const PetscInt a_mat_offset[], const PetscInt a_species_offset[], Mat subJ[], Mat JacP)
388 {
389   using scr_mem_t   = Kokkos::DefaultExecutionSpace::scratch_memory_space;
390   using real2_scr_t = Kokkos::View<PetscScalar **, Kokkos::LayoutRight, scr_mem_t>;
391   using g2_scr_t    = Kokkos::View<PetscReal ***, Kokkos::LayoutRight, scr_mem_t>;
392   using g3_scr_t    = Kokkos::View<PetscReal ****, Kokkos::LayoutRight, scr_mem_t>;
393   PetscInt         dim, num_cells_max, Nf_max, num_cells_batch;
394   int              nfaces = 0, vector_size = 512 / Nq;
395   LandauCtx       *ctx;
396   PetscReal       *d_Eq_m     = NULL;
397   PetscScalar     *d_vertex_f = NULL;
398   P4estVertexMaps *maps[LANDAU_MAX_GRIDS]; // this gets captured
399   PetscContainer   container;
400   const int        conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp == 0) ? Nq : 1;
401   const PetscInt   coo_sz_batch = SData_d->coo_size / batch_sz;                                                 // capture
402   auto             d_alpha_k    = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->alpha); //static data
403   const PetscReal *d_alpha      = d_alpha_k->data();
404   const PetscInt   Nftot        = d_alpha_k->size(); // total number of species
405   auto             d_beta_k     = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->beta);
406   const PetscReal *d_beta       = d_beta_k->data();
407   auto             d_invMass_k  = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invMass);
408   const PetscReal *d_invMass    = d_invMass_k->data();
409   auto             d_lambdas_k  = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->lambdas);
410   const PetscReal *d_lambdas    = d_lambdas_k->data();
411   auto             d_B_k        = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->B);
412   const PetscReal *d_BB         = d_B_k->data();
413   auto             d_D_k        = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->D);
414   const PetscReal *d_DD         = d_D_k->data();
415   auto             d_invJ_k     = *static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->invJ); // use Kokkos vector in kernels
416   auto             d_x_k        = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->x);     //static data
417   const PetscReal *d_x          = d_x_k->data();
418   auto             d_y_k        = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->y); //static data
419   const PetscReal *d_y          = d_y_k->data();
420   auto             d_z_k        = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->z); //static data
421   const PetscReal *d_z          = (LANDAU_DIM == 3) ? d_z_k->data() : NULL;
422   auto             d_w_k        = *static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->w); //static data
423   const PetscReal *d_w          = d_w_k.data();
424   // grid offsets - single vertex grid data
425   auto            d_numCells_k       = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->NCells);
426   const PetscInt *d_numCells         = d_numCells_k->data();
427   auto            d_species_offset_k = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->species_offset);
428   const PetscInt *d_species_offset   = d_species_offset_k->data();
429   auto            d_mat_offset_k     = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->mat_offset);
430   const PetscInt *d_mat_offset       = d_mat_offset_k->data();
431   auto            d_ip_offset_k      = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ip_offset);
432   const PetscInt *d_ip_offset        = d_ip_offset_k->data();
433   auto            d_ipf_offset_k     = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->ipf_offset);
434   const PetscInt *d_ipf_offset       = d_ipf_offset_k->data();
435   auto            d_elem_offset_k    = static_cast<Kokkos::View<PetscInt *, Kokkos::LayoutLeft> *>(SData_d->elem_offset);
436   const PetscInt *d_elem_offset      = d_elem_offset_k->data();
437   #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, including batched vertices
438   Kokkos::View<PetscReal ***, Kokkos::LayoutLeft> d_fdf_k = *static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutLeft> *>(SData_d->ipfdf_data);
439   #else
440   Kokkos::View<PetscReal ***, Kokkos::LayoutRight> d_fdf_k = *static_cast<Kokkos::View<PetscReal ***, Kokkos::LayoutRight> *>(SData_d->ipfdf_data);
441   #endif
442   auto d_Eq_m_k = static_cast<Kokkos::View<PetscReal *, Kokkos::LayoutLeft> *>(SData_d->Eq_m); // static storage, dynamic data - E(t), copy later, single vertex
443   // COO
444   auto       d_coo_elem_offsets_k                           = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>(SData_d->coo_elem_offsets);
445   LandauIdx *d_coo_elem_offsets                             = (SData_d->coo_size == 0) ? NULL : d_coo_elem_offsets_k->data();
446   auto       d_coo_elem_fullNb_k                            = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>(SData_d->coo_elem_fullNb);
447   LandauIdx *d_coo_elem_fullNb                              = (SData_d->coo_size == 0) ? NULL : d_coo_elem_fullNb_k->data();
448   auto       d_coo_elem_point_offsets_k                     = static_cast<Kokkos::View<LandauIdx *, Kokkos::LayoutLeft> *>(SData_d->coo_elem_point_offsets);
449   LandauIdx(*d_coo_elem_point_offsets)[LANDAU_MAX_NQND + 1] = (SData_d->coo_size == 0) ? NULL : (LandauIdx(*)[LANDAU_MAX_NQND + 1]) d_coo_elem_point_offsets_k->data();
450   auto         d_coo_vals_k                                 = static_cast<Kokkos::View<PetscScalar *, Kokkos::LayoutRight> *>(SData_d->coo_vals);
451   PetscScalar *d_coo_vals                                   = (SData_d->coo_size == 0) ? NULL : d_coo_vals_k->data();
452 
453   PetscFunctionBegin;
454   while (vector_size & (vector_size - 1)) vector_size = vector_size & (vector_size - 1);
455   if (vector_size > 16) vector_size = 16; // printf("DEBUG\n");
456   PetscCall(PetscLogEventBegin(events[3], 0, 0, 0, 0));
457   PetscCall(DMGetApplicationContext(plex[0], &ctx));
458   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
459   PetscCall(DMGetDimension(plex[0], &dim));
460   PetscCheck(LANDAU_DIM == dim, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM);
461   if (ctx->gpu_assembly) {
462     PetscCall(PetscObjectQuery((PetscObject)JacP, "assembly_maps", (PetscObject *)&container));
463     if (container) {
464       P4estVertexMaps *h_maps = NULL;
465       PetscCall(PetscContainerGetPointer(container, (void **)&h_maps));
466       for (PetscInt grid = 0; grid < num_grids; grid++) {
467         if (h_maps[grid].d_self) {
468           maps[grid] = h_maps[grid].d_self;
469           nfaces     = h_maps[grid].num_face; // nface=0 for CPU assembly
470         } else {
471           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container");
472         }
473       }
474       PetscCheck(d_coo_vals, PETSC_COMM_SELF, PETSC_ERR_PLIB, "d_coo_vals==0");
475     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container");
476   } else {
477     for (PetscInt grid = 0; grid < num_grids; grid++) maps[grid] = NULL;
478     nfaces    = 0;
479     container = NULL;
480   }
481   num_cells_batch = Nf_max = num_cells_max = 0;
482   for (PetscInt grid = 0; grid < num_grids; grid++) {
483     int Nfloc = a_species_offset[grid + 1] - a_species_offset[grid];
484     if (Nfloc > Nf_max) Nf_max = Nfloc;
485     if (a_numCells[grid] > num_cells_max) num_cells_max = a_numCells[grid];
486     num_cells_batch += a_numCells[grid]; // we don't have a host element offset here (but in ctx)
487   }
488   const int elem_mat_num_cells_max_grid = container ? 0 : num_cells_max;
489   #ifdef LAND_SUPPORT_CPU_ASS
490   const int                                           totDim_max = Nf_max * Nb, elem_mat_size_max = totDim_max * totDim_max;
491   Kokkos::View<PetscScalar ****, Kokkos::LayoutRight> d_elem_mats("element matrices", batch_sz, num_grids, elem_mat_num_cells_max_grid, elem_mat_size_max); // not used (cpu assembly)
492   #endif
493   const Kokkos::View<PetscReal *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_Eq_m_k(a_Eq_m, Nftot);
494   if (a_elem_closure || a_xarray) {
495     Kokkos::deep_copy(*d_Eq_m_k, h_Eq_m_k);
496     d_Eq_m = d_Eq_m_k->data();
497   } else d_Eq_m = NULL;
498   PetscCall(PetscKokkosInitializeCheck());
499   PetscCall(PetscLogEventEnd(events[3], 0, 0, 0, 0));
500   if (a_elem_closure || a_xarray) { // Jacobian, create f & df
501     Kokkos::View<PetscScalar *, Kokkos::LayoutLeft> *d_vertex_f_k = NULL;
502     PetscCall(PetscLogEventBegin(events[1], 0, 0, 0, 0));
503     if (a_elem_closure) {
504       PetscInt closure_sz = 0; // argh, don't have this on the host!!!
505       for (PetscInt grid = 0; grid < num_grids; grid++) {
506         PetscInt nfloc = a_species_offset[grid + 1] - a_species_offset[grid];
507         closure_sz += Nb * nfloc * a_numCells[grid];
508       }
509       closure_sz *= batch_sz;
510       d_vertex_f_k = new Kokkos::View<PetscScalar *, Kokkos::LayoutLeft>("closure", closure_sz);
511       const Kokkos::View<PetscScalar *, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_closure_k(a_elem_closure, closure_sz); // Vertex data for each element
512       Kokkos::deep_copy(*d_vertex_f_k, h_closure_k);
513       d_vertex_f = d_vertex_f_k->data();
514     } else {
515       d_vertex_f = (PetscScalar *)a_xarray;
516     }
517     PetscCall(PetscLogEventEnd(events[1], 0, 0, 0, 0));
518     PetscCall(PetscLogEventBegin(events[8], 0, 0, 0, 0));
519     PetscCall(PetscLogGpuTimeBegin());
520 
521     const int scr_bytes_fdf = real2_scr_t::shmem_size(Nf_max, Nb);
522     PetscCall(PetscInfo(plex[0], "df & f shared memory size: %d bytes in level, %d num cells total=%d, team size=%d, vector size=%d, #face=%d, Nf_max=%d\n", scr_bytes_fdf, KOKKOS_SHARED_LEVEL, (int)(num_cells_batch * batch_sz), team_size, vector_size, nfaces, (int)Nf_max));
523     Kokkos::parallel_for(
524       "f, df", Kokkos::TeamPolicy<>(num_cells_batch * batch_sz, team_size, /* Kokkos::AUTO */ vector_size).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes_fdf)), KOKKOS_LAMBDA(const team_member team) {
525         const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank() % b_Nelem, b_id = team.league_rank() / b_Nelem, IPf_sz_glb = d_ipf_offset[num_grids];
526         // find my grid
527         PetscInt grid = 0;
528         while (b_elem_idx >= d_elem_offset[grid + 1]) grid++;
529         {
530           const PetscInt loc_nip = d_numCells[grid] * Nq, loc_Nf = d_species_offset[grid + 1] - d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid];
531           const PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, d_mat_offset);
532           {
533             real2_scr_t      s_coef_k(team.team_scratch(KOKKOS_SHARED_LEVEL), Nf_max, Nb);
534             PetscScalar     *coef;
535             const PetscReal *invJe = &d_invJ_k((d_ip_offset[grid] + loc_elem * Nq) * dim * dim);
536             // un pack IPData
537             if (!maps[grid]) {
538               coef = &d_vertex_f[b_id * IPf_sz_glb + d_ipf_offset[grid] + loc_elem * Nb * loc_Nf]; // closure and IP indexing are the same
539             } else {
540               coef = s_coef_k.data();
541               Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, (int)loc_Nf), [=](int f) {
542                 //for (int f = 0; f < loc_Nf; ++f) {
543                 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, (int)Nb), [=](int b) {
544                   //for (int b = 0; b < Nb; ++b) {
545                   LandauIdx idx = maps[grid]->gIdx[loc_elem][f][b];
546                   if (idx >= 0) {
547                     coef[f * Nb + b] = d_vertex_f[idx + moffset]; // xarray data, not IP, need raw array to deal with CPU assemble case (not used)
548                   } else {
549                     idx              = -idx - 1;
550                     coef[f * Nb + b] = 0;
551                     for (int q = 0; q < maps[grid]->num_face; q++) {
552                       PetscInt    id    = maps[grid]->c_maps[idx][q].gid;
553                       PetscScalar scale = maps[grid]->c_maps[idx][q].scale;
554                       if (id >= 0) coef[f * Nb + b] += scale * d_vertex_f[id + moffset];
555                     }
556                   }
557                 });
558               });
559             }
560             team.team_barrier();
561             Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nq), [=](int myQi) {
562               const PetscReal *const invJ = &invJe[myQi * dim * dim]; // b_elem_idx: batch element index
563               const PetscReal       *Bq = &d_BB[myQi * Nb], *Dq = &d_DD[myQi * Nb * dim];
564               Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, (int)loc_Nf), [=](int f) {
565                 PetscInt       b, e, d;
566                 PetscReal      refSpaceDer[LANDAU_DIM];
567                 const PetscInt idx    = d_ipf_offset[grid] + f * loc_nip + loc_elem * Nq + myQi;
568                 d_fdf_k(b_id, 0, idx) = 0.0;
569                 for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0;
570                 for (b = 0; b < Nb; ++b) {
571                   d_fdf_k(b_id, 0, idx) += Bq[b] * PetscRealPart(coef[f * Nb + b]);
572                   for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[b * dim + d] * PetscRealPart(coef[f * Nb + b]);
573                 }
574                 for (d = 0; d < dim; ++d) {
575                   for (e = 0, d_fdf_k(b_id, d + 1, idx) = 0.0; e < dim; ++e) d_fdf_k(b_id, d + 1, idx) += invJ[e * dim + d] * refSpaceDer[e];
576                 }
577               }); // f
578             });   // myQi
579           }
580           team.team_barrier();
581         } // 'grid' scope
582       }); // global elems - fdf
583     Kokkos::fence();
584     PetscCall(PetscLogGpuTimeEnd());
585     PetscCall(PetscLogEventEnd(events[8], 0, 0, 0, 0));
586     // Jacobian
587     size_t      maximum_shared_mem_size = 64000;
588     PetscDevice device;
589 
590     PetscCall(PetscDeviceGetDefault_Internal(&device));
591     PetscCall(PetscDeviceGetAttribute(device, PETSC_DEVICE_ATTR_SIZE_T_SHARED_MEM_PER_BLOCK, &maximum_shared_mem_size));
592     size_t    jac_scr_bytes    = (size_t)2 * (g2_scr_t::shmem_size(dim, Nf_max, Nq) + g3_scr_t::shmem_size(dim, dim, Nf_max, Nq));
593     const int jac_shared_level = (jac_scr_bytes > maximum_shared_mem_size) ? 1 : KOKKOS_SHARED_LEVEL;
594     // device function/lambda
595     auto jac_lambda = KOKKOS_LAMBDA(const team_member team)
596     {
597       const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank() % b_Nelem, b_id = team.league_rank() / b_Nelem;
598       // find my grid
599       PetscInt grid = 0;
600       while (b_elem_idx >= d_elem_offset[grid + 1]) grid++;
601       {
602         const PetscInt loc_Nf = d_species_offset[grid + 1] - d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid];
603         const PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, d_mat_offset);
604         const PetscInt f_off   = d_species_offset[grid];
605   #ifdef LAND_SUPPORT_CPU_ASS
606         const PetscInt totDim = loc_Nf * Nb;
607   #endif
608         g2_scr_t g2(team.team_scratch(jac_shared_level), dim, loc_Nf, Nq);
609         g3_scr_t g3(team.team_scratch(jac_shared_level), dim, dim, loc_Nf, Nq);
610         g2_scr_t gg2(team.team_scratch(jac_shared_level), dim, loc_Nf, Nq);
611         g3_scr_t gg3(team.team_scratch(jac_shared_level), dim, dim, loc_Nf, Nq);
612         // get g2[] & g3[]
613         Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nq), [=](int myQi) {
614           using Kokkos::parallel_reduce;
615           const PetscInt                    jpidx_glb = d_ip_offset[grid] + loc_elem * Nq + myQi;
616           const PetscReal                  *invJ      = &d_invJ_k(jpidx_glb * dim * dim);
617           const PetscReal                   vj[3] = {d_x[jpidx_glb], d_y[jpidx_glb], d_z ? d_z[jpidx_glb] : 0}, wj = d_w[jpidx_glb];
618           landau_inner_red::TensorValueType gg_temp; // reduce on part of gg2 and g33 for IP jpidx_g
619           Kokkos::parallel_reduce(
620             Kokkos::ThreadVectorRange(team, (int)d_ip_offset[num_grids]),
621             [=](const int &ipidx, landau_inner_red::TensorValueType &ggg) {
622               const PetscReal wi = d_w[ipidx], x = d_x[ipidx], y = d_y[ipidx];
623               PetscReal       temp1[3] = {0, 0, 0}, temp2 = 0;
624               PetscInt        fieldB, d2, d3, f_off_r, grid_r, ipidx_g, nip_loc_r, loc_Nf_r;
625   #if LANDAU_DIM == 2
626               PetscReal Ud[2][2], Uk[2][2], mask = (PetscAbs(vj[0] - x) < 100 * PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[1] - y) < 100 * PETSC_SQRT_MACHINE_EPSILON) ? 0. : 1.;
627               LandauTensor2D(vj, x, y, Ud, Uk, mask);
628   #else
629                 PetscReal U[3][3], z = d_z[ipidx], mask = (PetscAbs(vj[0]-x) < 100*PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[1]-y) < 100*PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[2]-z) < 100*PETSC_SQRT_MACHINE_EPSILON) ? 0. : 1.;
630                 LandauTensor3D(vj, x, y, z, U, mask);
631   #endif
632               grid_r = 0;
633               while (ipidx >= d_ip_offset[grid_r + 1]) grid_r++; // yuck search for grid
634               f_off_r   = d_species_offset[grid_r];
635               ipidx_g   = ipidx - d_ip_offset[grid_r];
636               nip_loc_r = d_numCells[grid_r] * Nq;
637               loc_Nf_r  = d_species_offset[grid_r + 1] - d_species_offset[grid_r];
638               for (fieldB = 0; fieldB < loc_Nf_r; ++fieldB) { // fieldB is \beta  d_lambdas[grid][grid_r]
639                 const PetscInt  idx = d_ipf_offset[grid_r] + fieldB * nip_loc_r + ipidx_g;
640                 const PetscReal ff1 = d_beta[fieldB + f_off_r] * d_lambdas[LANDAU_MAX_GRIDS * grid + grid_r], ff2 = d_invMass[fieldB + f_off_r] * ff1;
641                 temp1[0] += d_fdf_k(b_id, 1, idx) * ff2;
642                 temp1[1] += d_fdf_k(b_id, 2, idx) * ff2;
643   #if LANDAU_DIM == 3
644                 temp1[2] += d_fdf_k(b_id, 3, idx) * ff2;
645   #endif
646                 temp2 += d_fdf_k(b_id, 0, idx) * ff1;
647               }
648               temp1[0] *= wi;
649               temp1[1] *= wi;
650   #if LANDAU_DIM == 3
651               temp1[2] *= wi;
652   #endif
653               temp2 *= wi;
654   #if LANDAU_DIM == 2
655               for (d2 = 0; d2 < 2; d2++) {
656                 for (d3 = 0; d3 < 2; ++d3) {
657                   /* K = U * grad(f): g2=e: i,A */
658                   ggg.gg2[d2] += Uk[d2][d3] * temp1[d3];
659                   /* D = -U * (I \kron (fx)): g3=f: i,j,A */
660                   ggg.gg3[d2][d3] += Ud[d2][d3] * temp2;
661                 }
662               }
663   #else
664                 for (d2 = 0; d2 < 3; ++d2) {
665                   for (d3 = 0; d3 < 3; ++d3) {
666                     /* K = U * grad(f): g2 = e: i,A */
667                     ggg.gg2[d2] += U[d2][d3]*temp1[d3];
668                     /* D = -U * (I \kron (fx)): g3 = f: i,j,A */
669                     ggg.gg3[d2][d3] += U[d2][d3]*temp2;
670                   }
671                 }
672   #endif
673             },
674             Kokkos::Sum<landau_inner_red::TensorValueType>(gg_temp));
675           // add alpha and put in gg2/3
676           Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, (int)loc_Nf), [&](const int &fieldA) { // \alpha
677             PetscInt d2, d3;
678             for (d2 = 0; d2 < dim; d2++) {
679               gg2(d2, fieldA, myQi) = gg_temp.gg2[d2] * d_alpha[fieldA + f_off];
680               for (d3 = 0; d3 < dim; d3++) gg3(d2, d3, fieldA, myQi) = -gg_temp.gg3[d2][d3] * d_alpha[fieldA + f_off] * d_invMass[fieldA + f_off];
681             }
682           });
683           /* add electric field term once per IP */
684           Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, (int)loc_Nf), [&](const int &fieldA) { gg2(dim - 1, fieldA, myQi) += d_Eq_m[fieldA + f_off]; });
685           Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, (int)loc_Nf), [=](const int &fieldA) {
686             int d, d2, d3, dp;
687             /* Jacobian transform - g2, g3 - per thread (2D) */
688             for (d = 0; d < dim; ++d) {
689               g2(d, fieldA, myQi) = 0;
690               for (d2 = 0; d2 < dim; ++d2) {
691                 g2(d, fieldA, myQi) += invJ[d * dim + d2] * gg2(d2, fieldA, myQi);
692                 g3(d, d2, fieldA, myQi) = 0;
693                 for (d3 = 0; d3 < dim; ++d3) {
694                   for (dp = 0; dp < dim; ++dp) g3(d, d2, fieldA, myQi) += invJ[d * dim + d3] * gg3(d3, dp, fieldA, myQi) * invJ[d2 * dim + dp];
695                 }
696                 g3(d, d2, fieldA, myQi) *= wj;
697               }
698               g2(d, fieldA, myQi) *= wj;
699             }
700           });
701         }); // Nq
702         team.team_barrier();
703         { /* assemble */
704           for (PetscInt fieldA = 0; fieldA < loc_Nf; fieldA++) {
705             /* assemble */
706             Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nb), [=](int f) {
707               Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, Nb), [=](int g) {
708                 PetscScalar t = 0;
709                 for (int qj = 0; qj < Nq; qj++) { // look at others integration points
710                   const PetscReal *BJq = &d_BB[qj * Nb], *DIq = &d_DD[qj * Nb * dim];
711                   for (int d = 0; d < dim; ++d) {
712                     t += DIq[f * dim + d] * g2(d, fieldA, qj) * BJq[g];
713                     for (int d2 = 0; d2 < dim; ++d2) t += DIq[f * dim + d] * g3(d, d2, fieldA, qj) * DIq[g * dim + d2];
714                   }
715                 }
716                 if (elem_mat_num_cells_max_grid) { // CPU assembly
717   #ifdef LAND_SUPPORT_CPU_ASS
718                   const PetscInt fOff                     = (fieldA * Nb + f) * totDim + fieldA * Nb + g;
719                   d_elem_mats(b_id, grid, loc_elem, fOff) = t;
720   #endif
721                 } else {
722                   static_cast<void>(landau_mat_assemble(d_coo_vals, t, f, g, Nb, moffset, loc_elem, fieldA, maps[grid], d_coo_elem_offsets, d_coo_elem_fullNb, d_coo_elem_point_offsets, b_elem_idx, b_id * coo_sz_batch));
723                 }
724               });
725             });
726           }
727         }
728       } // scope with 'grid'
729     };
730   #if defined(PETSC_HAVE_HIP)
731     const int lbound2 = 1;
732   #else
733     const int lbound2 = 2;
734   #endif
735     PetscCall(PetscLogEventBegin(events[4], 0, 0, 0, 0));
736     PetscCall(PetscLogGpuTimeBegin());
737     PetscCall(PetscInfo(plex[0], "Jacobian shared memory size: %zu bytes in level %s (max shared=%zu), num cells total=%d, team size=%d, vector size=%d, #face=%d, Nf_max=%d\n", jac_scr_bytes, jac_shared_level == 0 ? "local" : "global", maximum_shared_mem_size, (int)(num_cells_batch * batch_sz), team_size, vector_size, nfaces, (int)Nf_max));
738     Kokkos::parallel_for("Jacobian", Kokkos::TeamPolicy<Kokkos::LaunchBounds<256, lbound2>>(num_cells_batch * batch_sz, team_size, vector_size).set_scratch_size(jac_shared_level, Kokkos::PerTeam(jac_scr_bytes)), jac_lambda); // Kokkos::LaunchBounds<512,2>
739     Kokkos::fence();
740     PetscCall(PetscLogGpuTimeEnd());
741     PetscCall(PetscLogEventEnd(events[4], 0, 0, 0, 0));
742     if (d_vertex_f_k) delete d_vertex_f_k;
743   } else { // mass
744     PetscCall(PetscLogEventBegin(events[16], 0, 0, 0, 0));
745     PetscCall(PetscLogGpuTimeBegin());
746     PetscCall(PetscInfo(plex[0], "Mass team size=%d vector size=%d #face=%d Nb=%" PetscInt_FMT ", %s assembly\n", team_size, vector_size, nfaces, Nb, d_coo_vals ? "COO" : "CSR"));
747     Kokkos::parallel_for(
748       "Mass", Kokkos::TeamPolicy<>(num_cells_batch * batch_sz, team_size, vector_size), KOKKOS_LAMBDA(const team_member team) { // Kokkos::LaunchBounds<512,4>
749         const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank() % b_Nelem, b_id = team.league_rank() / b_Nelem;
750         // find my grid
751         PetscInt grid = 0;
752         while (b_elem_idx >= d_elem_offset[grid + 1]) grid++;
753         {
754           const PetscInt loc_Nf = d_species_offset[grid + 1] - d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid], jpidx_0 = d_ip_offset[grid] + loc_elem * Nq;
755   #ifdef LAND_SUPPORT_CPU_ASS
756           const PetscInt totDim = loc_Nf * Nb;
757   #endif
758           const PetscInt moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, d_mat_offset);
759           for (int fieldA = 0; fieldA < loc_Nf; fieldA++) {
760             /* assemble */
761             Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, Nb), [=](int f) {
762               Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, Nb), [=](int g) {
763                 PetscScalar t = 0;
764                 for (int qj = 0; qj < Nq; qj++) { // look at others integration points
765                   const PetscReal *BJq       = &d_BB[qj * Nb];
766                   const PetscInt   jpidx_glb = jpidx_0 + qj;
767                   if (dim == 2) {
768                     t += BJq[f] * d_w_k(jpidx_glb) * shift * BJq[g] * 2. * PETSC_PI;
769                   } else {
770                     t += BJq[f] * d_w_k(jpidx_glb) * shift * BJq[g];
771                   }
772                 }
773                 if (elem_mat_num_cells_max_grid) {
774   #ifdef LAND_SUPPORT_CPU_ASS
775                   const PetscInt fOff                     = (fieldA * Nb + f) * totDim + fieldA * Nb + g;
776                   d_elem_mats(b_id, grid, loc_elem, fOff) = t;
777   #endif
778                 } else {
779                   static_cast<void>(landau_mat_assemble(d_coo_vals, t, f, g, Nb, moffset, loc_elem, fieldA, maps[grid], d_coo_elem_offsets, d_coo_elem_fullNb, d_coo_elem_point_offsets, b_elem_idx, b_id * coo_sz_batch));
780                 }
781               });
782             });
783           } // field
784         } // grid
785       });
786     Kokkos::fence();
787     PetscCall(PetscLogGpuTimeEnd());
788     PetscCall(PetscLogEventEnd(events[16], 0, 0, 0, 0));
789   }
790   if (d_coo_vals) {
791     PetscCall(MatSetValuesCOO(JacP, d_coo_vals, ADD_VALUES));
792   #ifndef LAND_SUPPORT_CPU_ASS
793     Kokkos::fence(); // for timers
794   #endif
795   } else if (elem_mat_num_cells_max_grid) { // CPU assembly
796   #ifdef LAND_SUPPORT_CPU_ASS
797     Kokkos::View<PetscScalar ****, Kokkos::LayoutRight>::HostMirror h_elem_mats = Kokkos::create_mirror_view(d_elem_mats);
798     Kokkos::deep_copy(h_elem_mats, d_elem_mats);
799     PetscCheck(!container, PETSC_COMM_SELF, PETSC_ERR_PLIB, "?????");
800     for (PetscInt b_id = 0; b_id < batch_sz; b_id++) { // OpenMP (once)
801       for (PetscInt grid = 0; grid < num_grids; grid++) {
802         PetscSection       section, globalSection;
803         PetscInt           cStart, cEnd;
804         Mat                B       = subJ[LAND_PACK_IDX(b_id, grid)];
805         PetscInt           moffset = LAND_MOFFSET(b_id, grid, batch_sz, num_grids, a_mat_offset), nloc, nzl, colbuf[1024], row;
806         const PetscInt    *cols;
807         const PetscScalar *vals;
808         PetscCall(PetscLogEventBegin(events[5], 0, 0, 0, 0));
809         PetscCall(DMPlexGetHeightStratum(plex[grid], 0, &cStart, &cEnd));
810         PetscCall(DMGetLocalSection(plex[grid], &section));
811         PetscCall(DMGetGlobalSection(plex[grid], &globalSection));
812         PetscCall(PetscLogEventEnd(events[5], 0, 0, 0, 0));
813         PetscCall(PetscLogEventBegin(events[6], 0, 0, 0, 0));
814         for (PetscInt ej = cStart; ej < cEnd; ++ej) {
815           const PetscScalar *elMat = &h_elem_mats(b_id, grid, ej - cStart, 0);
816           PetscCall(DMPlexMatSetClosure(plex[grid], section, globalSection, B, ej, elMat, ADD_VALUES));
817           if (grid == 0 && ej == -1) {
818             const PetscInt loc_Nf = a_species_offset[grid + 1] - a_species_offset[grid], totDim = loc_Nf * Nb;
819             int            d, f;
820             PetscPrintf(PETSC_COMM_SELF, "Kokkos Element matrix %d/%d\n", 1, (int)a_numCells[grid]);
821             for (d = 0; d < totDim; ++d) {
822               for (f = 0; f < totDim; ++f) PetscPrintf(PETSC_COMM_SELF, " %12.5e", PetscRealPart(elMat[d * totDim + f]));
823               PetscPrintf(PETSC_COMM_SELF, "\n");
824             }
825             exit(14);
826           }
827         }
828         // move nest matrix to global JacP
829         PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
830         PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
831         PetscCall(MatGetSize(B, &nloc, NULL));
832         for (int i = 0; i < nloc; i++) {
833           PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
834           PetscCheck(nzl <= 1024, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "Row too big: %" PetscInt_FMT, nzl);
835           for (int j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset;
836           row = i + moffset;
837           PetscCall(MatSetValues(JacP, 1, &row, nzl, colbuf, vals, ADD_VALUES));
838           PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
839         }
840         PetscCall(MatDestroy(&subJ[LAND_PACK_IDX(b_id, grid)]));
841         PetscCall(PetscLogEventEnd(events[6], 0, 0, 0, 0));
842       } // grids
843     }
844   #else
845     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "CPU assembly not supported");
846   #endif
847   }
848   PetscFunctionReturn(PETSC_SUCCESS);
849 }
850 } // extern "C"
851 #endif
852