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