xref: /petsc/src/ts/utils/dmplexlandau/kokkos/landau.kokkos.cxx (revision 7d5fd1e4d9337468ad3f05b65b7facdcd2dfd2a4)
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.Nq = Nq;
93   h_maps.c_maps = (pointInterpolationP4est (*)[LANDAU_MAX_Q_FACE]) d_points->data();
94   maps[grid].vp1 = (void*)d_points;
95   h_maps.gIdx = (LandauIdx (*)[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ]) d_gidx->data();
96   maps[grid].vp2 = (void*)d_gidx;
97   {
98     Kokkos::View<P4estVertexMaps, Kokkos::HostSpace> h_maps_k(&h_maps);
99     Kokkos::View<P4estVertexMaps>                    *d_maps_k = new Kokkos::View<P4estVertexMaps>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(),h_maps_k));
100     Kokkos::deep_copy (*d_maps_k, h_maps_k);
101     maps[grid].d_self = d_maps_k->data();
102     maps[grid].vp3 = (void*)d_maps_k;
103   }
104   PetscFunctionReturn(0);
105 }
106 PetscErrorCode LandauKokkosDestroyMatMaps(P4estVertexMaps maps[], PetscInt num_grids)
107 {
108   PetscFunctionBegin;
109   for (PetscInt grid=0;grid<num_grids;grid++) {
110     auto a = static_cast<Kokkos::View<pointInterpolationP4est*[LANDAU_MAX_Q_FACE], Kokkos::LayoutRight>*>(maps[grid].vp1);
111     auto b = static_cast<Kokkos::View<LandauIdx*[LANDAU_MAX_SPECIES][LANDAU_MAX_NQ], Kokkos::LayoutRight>*>(maps[grid].vp2);
112     auto c = static_cast<Kokkos::View<P4estVertexMaps*>*>(maps[grid].vp3);
113     delete a;  delete b;  delete c;
114   }
115   PetscFunctionReturn(0);
116 }
117 
118 PetscErrorCode LandauKokkosStaticDataSet(DM plex, const PetscInt Nq, const PetscInt num_grids, 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   PetscErrorCode  ierr;
124   PetscTabulation *Tf;
125   PetscInt        dim;
126   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;
127   PetscDS         prob;
128 
129   PetscFunctionBegin;
130   ierr = DMGetDimension(plex, &dim);CHKERRQ(ierr);
131   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
132   if (LANDAU_DIM != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM);
133   ierr = PetscDSGetTabulation(prob, &Tf);CHKERRQ(ierr);
134   BB   = Tf[0]->T[0]; DD = Tf[0]->T[1];
135   ip_offset[0] = ipf_offset[0] = elem_offset[0] = 0;
136   nip = 0;
137   IPf_sz = 0;
138   for (PetscInt grid=0 ; grid<num_grids ; grid++) {
139     PetscInt nfloc = a_species_offset[grid+1] - a_species_offset[grid];
140     elem_offset[grid+1] = elem_offset[grid] + a_numCells[grid];
141     nip += a_numCells[grid]*Nq;
142     ip_offset[grid+1] = nip;
143     IPf_sz += Nq*nfloc*a_numCells[grid];
144     ipf_offset[grid+1] = IPf_sz;
145   }
146   Nftot = a_species_offset[num_grids];
147   ierr = PetscKokkosInitializeCheck();CHKERRQ(ierr);
148   {
149     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_alpha (a_nu_alpha, Nftot);
150     auto alpha = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("alpha", Nftot);
151     SData_d->alpha = static_cast<void*>(alpha);
152     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_beta (a_nu_beta, Nftot);
153     auto beta = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("beta", Nftot);
154     SData_d->beta = static_cast<void*>(beta);
155     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_invMass (a_invMass,Nftot);
156     auto invMass = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("invMass", Nftot);
157     SData_d->invMass = static_cast<void*>(invMass);
158     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_BB (BB,Nq*Nb);
159     auto B = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("B", Nq*Nb);
160     SData_d->B = static_cast<void*>(B);
161     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_DD (DD,Nq*Nb*dim);
162     auto D = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("D", Nq*Nb*dim);
163     SData_d->D = static_cast<void*>(D);
164     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_invJ (a_invJ, nip*dim*dim);
165     auto invJ = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("invJ", nip*dim*dim);
166     SData_d->invJ = static_cast<void*>(invJ);
167     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_x (a_x, nip);
168     auto x = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("x", nip);
169     SData_d->x = static_cast<void*>(x);
170     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_y (a_y, nip);
171     auto y = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("y", nip);
172     SData_d->y = static_cast<void*>(y);
173     const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_w (a_w, nip);
174     auto w = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("w", nip);
175     SData_d->w = static_cast<void*>(w);
176 
177     Kokkos::deep_copy (*alpha, h_alpha);
178     Kokkos::deep_copy (*beta, h_beta);
179     Kokkos::deep_copy (*invMass, h_invMass);
180     Kokkos::deep_copy (*B, h_BB);
181     Kokkos::deep_copy (*D, h_DD);
182     Kokkos::deep_copy (*invJ, h_invJ);
183     Kokkos::deep_copy (*x, h_x);
184     Kokkos::deep_copy (*y, h_y);
185     Kokkos::deep_copy (*w, h_w);
186 
187     if (dim==3) {
188       const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_z (a_z , dim==3 ? nip : 0);
189       auto z = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("z", nip);
190       SData_d->z = static_cast<void*>(z);
191       Kokkos::deep_copy (*z, h_z);
192     } else SData_d->z = NULL;
193 
194     //
195     const Kokkos::View<PetscInt*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_NCells (a_numCells, num_grids);
196     auto nc = new Kokkos::View<PetscInt*, Kokkos::LayoutLeft> ("NCells",num_grids);
197     SData_d->NCells = static_cast<void*>(nc);
198     Kokkos::deep_copy (*nc, h_NCells);
199 
200     const Kokkos::View<PetscInt*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_species_offset (a_species_offset, num_grids+1);
201     auto soff = new Kokkos::View<PetscInt*, Kokkos::LayoutLeft> ("species_offset",num_grids+1);
202     SData_d->species_offset = static_cast<void*>(soff);
203     Kokkos::deep_copy (*soff, h_species_offset);
204 
205     const Kokkos::View<PetscInt*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_mat_offset (a_mat_offset, num_grids+1);
206     auto moff = new Kokkos::View<PetscInt*, Kokkos::LayoutLeft> ("mat_offset",num_grids+1);
207     SData_d->mat_offset = static_cast<void*>(moff);
208     Kokkos::deep_copy (*moff, h_mat_offset);
209 
210     const Kokkos::View<PetscInt*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_ip_offset (ip_offset, num_grids+1);
211     auto ipoff = new Kokkos::View<PetscInt*, Kokkos::LayoutLeft> ("ip_offset",num_grids+1);
212     SData_d->ip_offset = static_cast<void*>(ipoff);
213     Kokkos::deep_copy (*ipoff,  h_ip_offset);
214 
215     const Kokkos::View<PetscInt*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_elem_offset (elem_offset, num_grids+1);
216     auto eoff = new Kokkos::View<PetscInt*, Kokkos::LayoutLeft> ("elem_offset",num_grids+1);
217     SData_d->elem_offset = static_cast<void*>(eoff);
218     Kokkos::deep_copy (*eoff,  h_elem_offset);
219 
220     const Kokkos::View<PetscInt*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_ipf_offset (ipf_offset, num_grids+1);
221     auto ipfoff = new Kokkos::View<PetscInt*, Kokkos::LayoutLeft> ("ipf_offset",num_grids+1);
222     SData_d->ipf_offset = static_cast<void*>(ipfoff);
223     Kokkos::deep_copy (*ipfoff,  h_ipf_offset);
224 
225 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
226     auto ipfdf_data =  new Kokkos::View<PetscReal**, Kokkos::LayoutLeft > ("fdf", dim+1, IPf_sz);
227 #else
228     auto ipfdf_data =  new Kokkos::View<PetscReal**, Kokkos::LayoutRight > ("fdf", dim+1, IPf_sz);
229 #endif
230     SData_d->ipfdf_data = static_cast<void*>(ipfdf_data);
231 
232     auto Eq_m = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("Eq_m",Nftot); // allocate but do not set
233     SData_d->Eq_m = static_cast<void*>(Eq_m);
234   }
235   SData_d->maps = NULL; // not used
236   PetscFunctionReturn(0);
237 }
238 
239 PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData *SData_d)
240 {
241   PetscFunctionBegin;
242   if (SData_d->alpha) {
243     auto alpha = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->alpha);
244     delete alpha;
245     SData_d->alpha = NULL;
246     auto beta = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->beta);
247     delete beta;
248     auto invMass = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invMass);
249     delete invMass;
250     auto B = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->B);
251     delete B;
252     auto D = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->D);
253     delete D;
254     auto invJ = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invJ);
255     delete invJ;
256     auto x = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->x);
257     delete x;
258     auto y = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->y);
259     delete y;
260     if (SData_d->z) {
261       auto z = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->z);
262       delete z;
263     }
264 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
265     auto z = static_cast<Kokkos::View<PetscReal**, Kokkos::LayoutLeft>*>(SData_d->ipfdf_data);
266 #else
267     auto z = static_cast<Kokkos::View<PetscReal**, Kokkos::LayoutRight>*>(SData_d->ipfdf_data);
268 #endif
269     delete z;
270     auto w = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->w);
271     delete w;
272     auto Eq_m = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->Eq_m);
273     delete Eq_m;
274     // offset
275     auto nc = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->NCells);
276     delete nc;
277     auto soff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->species_offset);
278     delete soff;
279     auto moff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->mat_offset);
280     delete moff;
281     auto ipoff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ip_offset);
282     delete ipoff;
283     auto eoff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->elem_offset);
284     delete eoff;
285     auto ipfoff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ipf_offset);
286     delete ipfoff;
287   }
288   PetscFunctionReturn(0);
289 }
290 
291 #define KOKKOS_SHARED_LEVEL 1
292 KOKKOS_INLINE_FUNCTION
293 PetscErrorCode landau_mat_assemble(PetscSplitCSRDataStructure d_mat, const team_member team,
294                                    Kokkos::View<PetscScalar**, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace::scratch_memory_space> s_fieldMats,
295                                    Kokkos::View<PetscInt**, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace::scratch_memory_space>  s_idx,
296                                    Kokkos::View<PetscReal**, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace::scratch_memory_space> s_scale,
297                                    const PetscInt Nb, const PetscInt Nq, const PetscInt nfaces, const PetscInt moffset, const PetscInt elem, const PetscInt fieldA, const P4estVertexMaps *d_maps)
298 {
299   const LandauIdx *const Idxs = &d_maps->gIdx[elem][fieldA][0];
300   team.team_barrier();
301   Kokkos::parallel_for(Kokkos::TeamVectorRange(team,0,Nb), [=] (int f) {
302       PetscInt q, idx = Idxs[f];
303       if (idx >= 0) {
304         s_idx(f,0) = idx + moffset;
305         s_scale(f,0) = 1.;
306       } else {
307         idx = -idx - 1;
308         for (q = 0; q < nfaces; q++) {
309           s_idx(f,q) = d_maps->c_maps[idx][q].gid + moffset;
310           s_scale(f,q) = d_maps->c_maps[idx][q].scale;
311         }
312       }
313     });
314   team.team_barrier();
315   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) {
316       PetscInt nr,idx = Idxs[f];
317       if (idx >= 0) {
318         nr = 1;
319       } else {
320         nr = nfaces;
321       }
322       Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) {
323           PetscScalar     vals[LANDAU_MAX_Q_FACE*LANDAU_MAX_Q_FACE];
324           PetscInt        q,d,nc,idx = Idxs[g];
325           if (idx >= 0) {
326             nc = 1;
327           } else {
328             nc = nfaces;
329           }
330           for (q = 0; q < nr; q++) {
331             for (d = 0; d < nc; d++) {
332               vals[q*nc + d] = s_scale(f,q)*s_scale(g,d)*s_fieldMats(f,g);
333             }
334           }
335           MatSetValuesDevice(d_mat,nr,&s_idx(f,0),nc,&s_idx(g,0),vals,ADD_VALUES);
336         });
337     });
338   return 0;
339 }
340 
341 PetscErrorCode LandauKokkosJacobian(DM plex[], const PetscInt Nq, const PetscInt num_grids, const PetscInt a_numCells[], PetscReal a_Eq_m[], PetscScalar a_elem_closure[],
342                                     const PetscInt N, const PetscScalar a_xarray[], const LandauStaticData *SData_d, const PetscInt num_sub_blocks, const PetscReal shift,
343                                     const PetscLogEvent events[], const PetscInt a_mat_offset[], const PetscInt a_species_offset[], Mat subJ[], Mat JacP)
344 {
345   using scr_mem_t = Kokkos::DefaultExecutionSpace::scratch_memory_space;
346   using fieldMats_scr_t = Kokkos::View<PetscScalar**, Kokkos::LayoutRight, scr_mem_t>;
347   using idx_scr_t = Kokkos::View<PetscInt**, Kokkos::LayoutRight, scr_mem_t>;
348   using scale_scr_t = Kokkos::View<PetscReal**, Kokkos::LayoutRight, scr_mem_t>;
349   using g2_scr_t = Kokkos::View<PetscReal***, Kokkos::LayoutRight, scr_mem_t>;
350   using g3_scr_t = Kokkos::View<PetscReal****, Kokkos::LayoutRight, scr_mem_t>;
351   PetscErrorCode    ierr;
352   PetscInt          Nb=Nq,dim,num_cells_max,num_cells_tot,Nf_max,nip_global;
353   int               nfaces=0;
354   LandauCtx         *ctx;
355   PetscReal         *d_Eq_m=NULL;
356   PetscScalar       *d_vertex_f=NULL;
357   P4estVertexMaps   *maps[LANDAU_MAX_GRIDS]; // this gets captured
358   PetscSplitCSRDataStructure d_mat;
359   PetscContainer    container = NULL;
360   const int         conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp==0) ? Nq : 1;
361   auto              d_alpha_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->alpha); //static data
362   const PetscReal   *d_alpha = d_alpha_k->data();
363   const PetscInt    Nftot = d_alpha_k->size(); // total number of species
364   auto              d_beta_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->beta);
365   const PetscReal   *d_beta = d_beta_k->data();
366   auto              d_invMass_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invMass);
367   const PetscReal   *d_invMass = d_invMass_k->data();
368   auto              d_B_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->B);
369   const PetscReal   *d_BB = d_B_k->data();
370   auto              d_D_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->D);
371   const PetscReal   *d_DD = d_D_k->data();
372   auto              d_invJ_k = *static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invJ);     // use Kokkos vector in kernels
373   auto              d_x_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->x); //static data
374   const PetscReal   *d_x = d_x_k->data();
375   auto              d_y_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->y); //static data
376   const PetscReal   *d_y = d_y_k->data();
377   auto              d_z_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->z); //static data
378   const PetscReal   *d_z = (LANDAU_DIM==3) ? d_z_k->data() : NULL;
379   auto              d_w_k = *static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->w); //static data
380   const PetscReal   *d_w = d_w_k.data();
381   auto              d_Eq_m_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->Eq_m); // static storage, dynamci data - E(t), copy later
382   // grid offset d_numCells
383   auto              d_numCells_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->NCells);
384   const PetscInt   *d_numCells = d_numCells_k->data();
385   auto              d_species_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->species_offset);
386   const PetscInt   *d_species_offset = d_species_offset_k->data();
387   auto              d_mat_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->mat_offset);
388   const PetscInt   *d_mat_offset = d_mat_offset_k->data();
389   auto              d_ip_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ip_offset);
390   const PetscInt   *d_ip_offset = d_ip_offset_k->data();
391   auto              d_ipf_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ipf_offset);
392   const PetscInt   *d_ipf_offset = d_ipf_offset_k->data();
393   auto              d_elem_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->elem_offset);
394   const PetscInt   *d_elem_offset = d_elem_offset_k->data();
395 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
396   Kokkos::View<PetscReal**, Kokkos::LayoutLeft > d_fdf_k = *static_cast<Kokkos::View<PetscReal**, Kokkos::LayoutLeft >*>(SData_d->ipfdf_data);
397 #else
398   Kokkos::View<PetscReal**, Kokkos::LayoutRight > d_fdf_k = *static_cast<Kokkos::View<PetscReal**, Kokkos::LayoutRight >*>(SData_d->ipfdf_data);
399 #endif
400   PetscFunctionBegin;
401   ierr = PetscLogEventBegin(events[3],0,0,0,0);CHKERRQ(ierr);
402   ierr = DMGetApplicationContext(plex[0], &ctx);CHKERRQ(ierr);
403   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
404   ierr = DMGetDimension(plex[0], &dim);CHKERRQ(ierr);
405   if (LANDAU_DIM != dim) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM);
406   if (ctx->gpu_assembly) {
407     ierr = PetscObjectQuery((PetscObject) JacP, "assembly_maps", (PetscObject *) &container);CHKERRQ(ierr);
408     if (container) { // not here first call
409       P4estVertexMaps   *h_maps=NULL;
410       ierr = PetscContainerGetPointer(container, (void **) &h_maps);CHKERRQ(ierr);
411       for (PetscInt grid=0 ; grid<num_grids ; grid++) {
412         if (h_maps[grid].d_self) {
413           maps[grid] = h_maps[grid].d_self;
414           nfaces = h_maps[grid].num_face; // nface=0 for CPU assembly
415         } else {
416           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container");
417         }
418       }
419       // this does the setup the first time called
420       ierr = MatKokkosGetDeviceMatWrite(JacP,&d_mat);CHKERRQ(ierr);
421     } else { // kernel output - first call assembled on device
422       for (PetscInt grid=0 ; grid<num_grids ; grid++) maps[grid] = NULL;
423       nfaces = 0;
424     }
425   } else {
426     for (PetscInt grid=0 ; grid<num_grids ; grid++) maps[grid] = NULL;
427     nfaces = 0;
428   }
429   nip_global = num_cells_tot = Nf_max = num_cells_max = 0;
430   for (PetscInt grid=0 ; grid<num_grids ; grid++) {
431     int Nfloc = a_species_offset[grid+1] - a_species_offset[grid];
432     if (Nfloc > Nf_max) Nf_max = Nfloc;
433     if (a_numCells[grid] > num_cells_max) num_cells_max = a_numCells[grid];
434     num_cells_tot += a_numCells[grid];
435     nip_global += Nq*a_numCells[grid];
436   }
437   const PetscInt totDim_max = Nf_max*Nq, elem_mat_size_max = totDim_max*totDim_max;
438   const PetscInt elem_mat_num_cells_max_grid = container ? 0 : num_cells_max;
439   Kokkos::View<PetscScalar***, Kokkos::LayoutRight> d_elem_mats("element matrices", num_grids, elem_mat_num_cells_max_grid, elem_mat_size_max); // first call have large set of global element matrices
440   const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> >  h_Eq_m_k (a_Eq_m, Nftot);
441   if (a_elem_closure || a_xarray) {
442     Kokkos::deep_copy (*d_Eq_m_k, h_Eq_m_k);
443     d_Eq_m = d_Eq_m_k->data();
444   } else d_Eq_m = NULL;
445   ierr = PetscKokkosInitializeCheck();CHKERRQ(ierr);
446   ierr = PetscLogEventEnd(events[3],0,0,0,0);CHKERRQ(ierr);
447   if (a_elem_closure || a_xarray) { // Jacobian, create f & df
448     Kokkos::View<PetscScalar*, Kokkos::LayoutLeft> *d_vertex_f_k = NULL;
449     ierr = PetscLogEventBegin(events[1],0,0,0,0);CHKERRQ(ierr);
450     if (a_elem_closure) {
451       PetscInt closure_sz = 0; // argh, don't have this on the host!!!
452       for (PetscInt grid=0 ; grid<num_grids ; grid++) {
453         PetscInt nfloc = a_species_offset[grid+1] - a_species_offset[grid];
454         closure_sz     += Nq*nfloc*a_numCells[grid];
455       }
456       d_vertex_f_k = new Kokkos::View<PetscScalar*, Kokkos::LayoutLeft> ("closure",closure_sz);
457       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
458       Kokkos::deep_copy (*d_vertex_f_k, h_closure_k);
459       d_vertex_f  = d_vertex_f_k->data();
460     } else {
461       d_vertex_f = (PetscScalar*)a_xarray;
462     }
463     ierr = PetscLogEventEnd(events[1],0,0,0,0);CHKERRQ(ierr);
464     ierr = PetscLogEventBegin(events[8],0,0,0,0);CHKERRQ(ierr);
465     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
466     Kokkos::parallel_for("f, df", Kokkos::TeamPolicy<>(num_cells_tot, team_size, /* Kokkos::AUTO */ 16), KOKKOS_LAMBDA (const team_member team) {
467         // find my grid
468         PetscInt grid = 0, g_cell = team.league_rank();
469         while (g_cell >= d_elem_offset[grid+1]) grid++; // yuck search for grid
470         {
471           const PetscInt     moffset = d_mat_offset[grid], nip_loc = d_numCells[grid]*Nq, Nfloc = d_species_offset[grid+1] - d_species_offset[grid], elem = g_cell -  d_elem_offset[grid];
472           const PetscInt     IP_idx = d_ip_offset[grid], IPf_idx = d_ipf_offset[grid];
473           const PetscScalar  *coef;
474           PetscScalar        coef_buff[LANDAU_MAX_SPECIES*LANDAU_MAX_NQ];
475           // un pack IPData
476           if (!maps[grid]) {
477             coef = &d_vertex_f[elem*Nb*Nfloc + IPf_idx]; // closure and IP indexing are the same
478           } else {
479             coef = coef_buff;
480             for (int f = 0; f < Nfloc; ++f) {
481               LandauIdx *const Idxs = &maps[grid]->gIdx[elem][f][0];
482               for (int b = 0; b < Nb; ++b) {
483                 PetscInt idx = Idxs[b];
484                 if (idx >= 0) {
485                   coef_buff[f*Nb+b] = d_vertex_f[idx+moffset]; // xarray data, not IP
486                 } else {
487                   idx = -idx - 1;
488                   coef_buff[f*Nb+b] = 0;
489                   for (int q = 0; q < maps[grid]->num_face; q++) {
490                     PetscInt    id = maps[grid]->c_maps[idx][q].gid;
491                     PetscScalar scale = maps[grid]->c_maps[idx][q].scale;
492                     coef_buff[f*Nb+b] += scale*d_vertex_f[id+moffset];
493                   }
494                 }
495               }
496             }
497           }
498           Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nq), [=] (int myQi) {
499               const PetscInt          ipidx = IP_idx + myQi + elem * Nq, ipidx_g = myQi + elem * Nq;
500               const PetscReal *const  invJj = &d_invJ_k(ipidx*dim*dim);
501               const PetscReal         *Bq = &d_BB[myQi*Nb], *Dq = &d_DD[myQi*Nb*dim];
502               Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,(int)Nfloc), [=] (int f) {
503                   PetscInt   b, e, d;
504                   PetscReal  refSpaceDer[LANDAU_DIM];
505                   const PetscInt idx = IPf_idx + f*nip_loc + ipidx_g;
506                   d_fdf_k(0,idx) = 0.0;
507                   for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0;
508                   for (b = 0; b < Nb; ++b) {
509                     const PetscInt    cidx = b;
510                     d_fdf_k(0,idx) += Bq[cidx]*PetscRealPart(coef[f*Nb+cidx]);
511                     for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx*dim+d]*PetscRealPart(coef[f*Nb+cidx]);
512                   }
513                   for (d = 0; d < dim; ++d) {
514                     for (e = 0, d_fdf_k(d+1,idx) = 0.0; e < dim; ++e) {
515                       d_fdf_k(d+1,idx) += invJj[e*dim+d]*refSpaceDer[e];
516                     }
517                   }
518                 }); // f
519             }); // myQi
520         }
521       }); // global elems - fdf
522     // #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
523     //     ierr = PetscLogGpuFlops(nip_loc*(PetscLogDouble)(2*Nb*(1+dim)));CHKERRQ(ierr);
524     // #else
525     //     ierr = PetscLogFlops(nip_loc*(PetscLogDouble)(2*Nb*(1+dim)));CHKERRQ(ierr);
526     // #endif
527     Kokkos::fence();
528     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr); // is this a fence?
529     ierr = PetscLogEventEnd(events[8],0,0,0,0);CHKERRQ(ierr);
530     if (d_vertex_f_k) delete d_vertex_f_k;
531     // Jacobian
532     ierr = PetscLogEventBegin(events[4],0,0,0,0);CHKERRQ(ierr);
533     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
534     // #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
535     //           ierr = PetscLogGpuFlops(nip_loc*(PetscLogDouble)((nip_loc*(11*Nfloc+ 4*dim*dim) + 6*Nfloc*dim*dim*dim + 10*Nfloc*dim*dim + 4*Nfloc*dim + Nb*Nfloc*Nb*Nq*dim*dim*5)));CHKERRQ(ierr);
536     // #else
537     //           ierr = PetscLogFlops(nip_loc*(PetscLogDouble)((nip_loc*(11*Nfloc+ 4*dim*dim) + 6*Nfloc*dim*dim*dim + 10*Nfloc*dim*dim + 4*Nfloc*dim + Nb*Nfloc*Nb*Nq*dim*dim*5)));CHKERRQ(ierr);
538     // #endif
539     const int scr_bytes = 2*(g2_scr_t::shmem_size(dim,Nf_max,Nq) + g3_scr_t::shmem_size(dim,dim,Nf_max,Nq))+fieldMats_scr_t::shmem_size(Nb,Nb)+idx_scr_t::shmem_size(Nb,nfaces)+scale_scr_t::shmem_size(Nb,nfaces);
540     ierr = PetscInfo6(plex[0], "Jacobian shared memory size: %d bytes in level %d num_cells_tot=%D team size=%D #face=%D Nf_max=%D\n",scr_bytes,KOKKOS_SHARED_LEVEL,num_cells_tot,team_size,nfaces,Nf_max);CHKERRQ(ierr);
541     Kokkos::parallel_for("Jacobian", Kokkos::TeamPolicy<>(num_cells_tot, team_size, /* Kokkos::AUTO */ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes)), KOKKOS_LAMBDA (const team_member team) {
542         // find my grid
543         PetscInt grid = 0, g_cell = team.league_rank();
544         while (g_cell >= d_elem_offset[grid+1]) grid++; // yuck search for grid
545         {
546           const PetscInt  moffset = d_mat_offset[grid], Nfloc = d_species_offset[grid+1]-d_species_offset[grid], elem = g_cell-d_elem_offset[grid], totDim = Nfloc*Nq;
547           const PetscInt  IP_idx = d_ip_offset[grid];
548           const PetscInt  f_off = d_species_offset[grid];
549           g2_scr_t        g2(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,Nfloc,Nq);
550           g3_scr_t        g3(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,dim,Nfloc,Nq);
551           g2_scr_t        gg2(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,Nfloc,Nq);
552           g3_scr_t        gg3(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,dim,Nfloc,Nq);
553           // get g2[] & g3[]
554           Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nq), [=] (int myQi) {
555               using Kokkos::parallel_reduce;
556               const PetscInt                    jpidx_g = myQi + elem * Nq, jpidx = IP_idx + jpidx_g;
557               const PetscReal* const            invJj = &d_invJ_k(jpidx*dim*dim);
558               const PetscReal                   vj[3] = {d_x[jpidx], d_y[jpidx], d_z ? d_z[jpidx] : 0}, wj = d_w[jpidx];
559               landau_inner_red::TensorValueType gg_temp; // reduce on part of gg2 and g33 for IP jpidx_g
560               Kokkos::parallel_reduce(Kokkos::ThreadVectorRange (team, (int)nip_global), [=] (const int& ipidx, landau_inner_red::TensorValueType & ggg) {
561                   const PetscReal wi = d_w[ipidx], x = d_x[ipidx], y = d_y[ipidx];
562                   PetscReal       temp1[3] = {0, 0, 0}, temp2 = 0;
563                   PetscInt        fieldA,d2,d3,f_off_r,grid_r,ipidx_g,nip_loc_r,Nfloc_r,IPf_idx_r;
564 #if LANDAU_DIM==2
565                   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.;
566                   LandauTensor2D(vj, x, y, Ud, Uk, mask);
567 #else
568                   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.;
569                   LandauTensor3D(vj, x, y, z, U, mask);
570 #endif
571                   grid_r = 0;
572                   while (ipidx >= d_ip_offset[grid_r+1]) grid_r++; // yuck search for grid
573                   f_off_r = d_species_offset[grid_r];
574                   ipidx_g = ipidx - d_ip_offset[grid_r];
575                   nip_loc_r = d_numCells[grid_r]*Nq;
576                   Nfloc_r = d_species_offset[grid_r+1] - d_species_offset[grid_r];
577                   IPf_idx_r = d_ipf_offset[grid_r];
578                   for (fieldA = 0; fieldA < Nfloc_r; ++fieldA) {
579                     const PetscInt idx = IPf_idx_r + fieldA*nip_loc_r + ipidx_g;
580                     temp1[0] += d_fdf_k(1,idx)*d_beta[fieldA+f_off_r]*d_invMass[fieldA+f_off_r];
581                     temp1[1] += d_fdf_k(2,idx)*d_beta[fieldA+f_off_r]*d_invMass[fieldA+f_off_r];
582 #if LANDAU_DIM==3
583                     temp1[2] += d_fdf_k(3,idx)*d_beta[fieldA+f_off_r]*d_invMass[fieldA+f_off_r];
584 #endif
585                     temp2    += d_fdf_k(0,idx)*d_beta[fieldA+f_off_r];
586                   }
587                   temp1[0] *= wi;
588                   temp1[1] *= wi;
589 #if LANDAU_DIM==3
590                   temp1[2] *= wi;
591 #endif
592                   temp2    *= wi;
593 #if LANDAU_DIM==2
594                   for (d2 = 0; d2 < 2; d2++) {
595                     for (d3 = 0; d3 < 2; ++d3) {
596                       /* K = U * grad(f): g2=e: i,A */
597                       ggg.gg2[d2] += Uk[d2][d3]*temp1[d3];
598                       /* D = -U * (I \kron (fx)): g3=f: i,j,A */
599                       ggg.gg3[d2][d3] += Ud[d2][d3]*temp2;
600                     }
601                   }
602 #else
603                   for (d2 = 0; d2 < 3; ++d2) {
604                     for (d3 = 0; d3 < 3; ++d3) {
605                       /* K = U * grad(f): g2 = e: i,A */
606                       ggg.gg2[d2] += U[d2][d3]*temp1[d3];
607                       /* D = -U * (I \kron (fx)): g3 = f: i,j,A */
608                       ggg.gg3[d2][d3] += U[d2][d3]*temp2;
609                     }
610                   }
611 #endif
612                 }, Kokkos::Sum<landau_inner_red::TensorValueType>(gg_temp));
613               // add alpha and put in gg2/3
614               Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)Nfloc), [&] (const int& fieldA) {
615                   PetscInt d2,d3;
616                   for (d2 = 0; d2 < dim; d2++) {
617                     gg2(d2,fieldA,myQi) = gg_temp.gg2[d2]*d_alpha[fieldA+f_off];
618                     for (d3 = 0; d3 < dim; d3++) {
619                       gg3(d2,d3,fieldA,myQi) = -gg_temp.gg3[d2][d3]*d_alpha[fieldA+f_off]*d_invMass[fieldA+f_off];
620                     }
621                   }
622                 });
623               /* add electric field term once per IP */
624               Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)Nfloc), [&] (const int& fieldA) {
625                   gg2(dim-1,fieldA,myQi) += d_Eq_m[fieldA+f_off];
626                 });
627               Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)Nfloc), [=] (const int& fieldA) {
628                   int d,d2,d3,dp;
629                   /* Jacobian transform - g2, g3 - per thread (2D) */
630                   for (d = 0; d < dim; ++d) {
631                     g2(d,fieldA,myQi) = 0;
632                     for (d2 = 0; d2 < dim; ++d2) {
633                       g2(d,fieldA,myQi) += invJj[d*dim+d2]*gg2(d2,fieldA,myQi);
634                       g3(d,d2,fieldA,myQi) = 0;
635                       for (d3 = 0; d3 < dim; ++d3) {
636                         for (dp = 0; dp < dim; ++dp) {
637                           g3(d,d2,fieldA,myQi) += invJj[d*dim + d3]*gg3(d3,dp,fieldA,myQi)*invJj[d2*dim + dp];
638                         }
639                       }
640                       g3(d,d2,fieldA,myQi) *= wj;
641                     }
642                     g2(d,fieldA,myQi) *= wj;
643                   }
644                 });
645             }); // Nq
646           team.team_barrier();
647           { /* assemble */
648             fieldMats_scr_t s_fieldMats(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,Nb); // Only used for GPU assembly (ie, not first pass)
649             idx_scr_t       s_idx(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,nfaces);
650             scale_scr_t     s_scale(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,nfaces);
651             for (PetscInt fieldA = 0; fieldA < Nfloc; fieldA++) {
652               /* assemble */
653               Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) {
654                   Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) {
655                       PetscScalar t = 0;
656                       for (int qj = 0 ; qj < Nq ; qj++) { // look at others integration points
657                         const PetscReal *BJq = &d_BB[qj*Nb], *DIq = &d_DD[qj*Nb*dim];
658                         for (int d = 0; d < dim; ++d) {
659                           t += DIq[f*dim+d]*g2(d,fieldA,qj)*BJq[g];
660                           for (int d2 = 0; d2 < dim; ++d2) {
661                             t += DIq[f*dim + d]*g3(d,d2,fieldA,qj)*DIq[g*dim + d2];
662                           }
663                         }
664                       }
665                       if (elem_mat_num_cells_max_grid) { // CPU assembly
666                         const PetscInt fOff = (fieldA*Nb + f)*totDim + fieldA*Nb + g;
667                         d_elem_mats(grid,elem,fOff) = t;
668                       } else s_fieldMats(f,g) = t;
669                     });
670                 });
671               if (!elem_mat_num_cells_max_grid) { // GPU assembly
672                 landau_mat_assemble (d_mat, team, s_fieldMats, s_idx, s_scale, Nb, Nq, nfaces, moffset, elem, fieldA, maps[grid]);
673               }
674             }
675           }
676         }
677       });
678     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
679     ierr = PetscLogEventEnd(events[4],0,0,0,0);CHKERRQ(ierr);
680   } else { // mass
681     ierr = PetscLogEventBegin(events[4],0,0,0,0);CHKERRQ(ierr);
682     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
683     // #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
684     //           ierr = PetscLogGpuFlops(nip_loc*(PetscLogDouble)(Nb*Nfloc*Nb*Nq*4));CHKERRQ(ierr);
685     //           if (ctx->deviceType == LANDAU_CPU) PetscInfo(plex[grid], "Warning: Landau selected CPU but no support for Kokkos using CPU\n");
686     // #else
687     //           ierr = PetscLogFlops(nip_loc*(PetscLogDouble)(Nb*Nfloc*Nb*Nq*4));CHKERRQ(ierr);
688     // #endif
689     int scr_bytes = fieldMats_scr_t::shmem_size(Nq,Nq) + idx_scr_t::shmem_size(Nb,nfaces) + scale_scr_t::shmem_size(Nb,nfaces);
690     ierr = PetscInfo6(plex[0], "Mass shared memory size: %d bytes in level %d conc=%D team size=%D #face=%D Nb=%D\n",scr_bytes,KOKKOS_SHARED_LEVEL,conc,team_size,nfaces,Nb);CHKERRQ(ierr);
691     Kokkos::parallel_for("Mass", Kokkos::TeamPolicy<>(num_cells_tot, team_size, /* Kokkos::AUTO */ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes)), KOKKOS_LAMBDA (const team_member team) {
692         fieldMats_scr_t s_fieldMats(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,Nb);
693         idx_scr_t       s_idx(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,nfaces);
694         scale_scr_t     s_scale(team.team_scratch(KOKKOS_SHARED_LEVEL),Nb,nfaces);
695         // find my grid
696         PetscInt grid = 0, g_cell = team.league_rank();
697         while (g_cell >= d_elem_offset[grid+1]) grid++; // yuck search for grid
698         {
699           const PetscInt moffset = d_mat_offset[grid], Nfloc = d_species_offset[grid+1]-d_species_offset[grid], totDim = Nfloc*Nq, elem = g_cell-d_elem_offset[grid];
700           const PetscInt IP_idx = d_ip_offset[grid];
701           for (int fieldA = 0; fieldA < Nfloc; fieldA++) {
702             /* assemble */
703             Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) {
704                 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) {
705                     PetscScalar    t = 0;
706                     for (int qj = 0 ; qj < Nq ; qj++) { // look at others integration points
707                       const PetscReal *BJq = &d_BB[qj*Nb];
708                       const PetscInt jpidx = IP_idx + qj + elem * Nq;
709                       if (dim==2) {
710                         t += BJq[f] * d_w_k(jpidx) * shift * BJq[g] * 2. * PETSC_PI;
711                       } else {
712                         t += BJq[f] * d_w_k(jpidx) * shift * BJq[g];
713                       }
714                     }
715                     if (elem_mat_num_cells_max_grid) {
716                       const PetscInt fOff = (fieldA*Nb + f)*totDim + fieldA*Nb + g;
717                       d_elem_mats(grid,elem,fOff) = t;
718                     } else s_fieldMats(f,g) = t;
719                   });
720               });
721             if (!elem_mat_num_cells_max_grid) { // device assembly
722               landau_mat_assemble (d_mat, team, s_fieldMats, s_idx, s_scale, Nb, Nq,nfaces, moffset, elem, fieldA, maps[grid]);
723             } // else not using GPU assembly
724           }
725         }
726       });
727     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
728     ierr = PetscLogEventEnd(events[4],0,0,0,0);CHKERRQ(ierr);
729   }
730   Kokkos::fence();
731   if (elem_mat_num_cells_max_grid) { // CPU assembly
732     Kokkos::View<PetscScalar***, Kokkos::LayoutRight>::HostMirror h_elem_mats = Kokkos::create_mirror_view(d_elem_mats);
733     Kokkos::deep_copy (h_elem_mats, d_elem_mats);
734     for (PetscInt grid=0 ; grid<num_grids ; grid++) {
735       const PetscInt  Nfloc = a_species_offset[grid+1]-a_species_offset[grid], totDim = Nfloc*Nq;
736       PetscSection    section, globalSection;
737       PetscInt        cStart,cEnd;
738       ierr = PetscLogEventBegin(events[5],0,0,0,0);CHKERRQ(ierr);
739       ierr = DMPlexGetHeightStratum(plex[grid],0,&cStart,&cEnd);CHKERRQ(ierr);
740       ierr = DMGetLocalSection(plex[grid], &section);CHKERRQ(ierr);
741       ierr = DMGetGlobalSection(plex[grid], &globalSection);CHKERRQ(ierr);
742       ierr = PetscLogEventEnd(events[5],0,0,0,0);CHKERRQ(ierr);
743       ierr = PetscLogEventBegin(events[6],0,0,0,0);CHKERRQ(ierr);
744       for (PetscInt ej = cStart ; ej < cEnd; ++ej) {
745         const PetscScalar *elMat = &h_elem_mats(grid,ej-cStart,0);
746         ierr = DMPlexMatSetClosure(plex[grid], section, globalSection, subJ[grid], ej, elMat, ADD_VALUES);CHKERRQ(ierr);
747         if (grid==0 && ej==-1) {
748           int d,f;
749           PetscPrintf(PETSC_COMM_SELF,"Kokkos Element matrix %d/%d\n",1,(int)a_numCells[grid]);
750           for (d = 0; d < totDim; ++d){
751             for (f = 0; f < totDim; ++f) PetscPrintf(PETSC_COMM_SELF," %12.5e",  PetscRealPart(elMat[d*totDim + f]));
752             PetscPrintf(PETSC_COMM_SELF,"\n");
753           }
754           exit(14);
755         }
756       }
757       // transition to use of maps for a Kokkos VecGetClosure
758       if (ctx->gpu_assembly) {
759         if (!(a_elem_closure || a_xarray)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "transition in Mass");
760       }
761       if (!container) {   // move nest matrix to global JacP
762         PetscInt          moffset = a_mat_offset[grid], nloc, nzl, colbuf[1024], row;
763         const PetscInt    *cols;
764         const PetscScalar *vals;
765         Mat               B = subJ[grid];
766         ierr = MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767         ierr = MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
768         ierr = MatGetSize(B, &nloc, NULL);CHKERRQ(ierr);
769         if (nloc != a_mat_offset[grid+1] - moffset) SETERRQ2(PetscObjectComm((PetscObject) B), PETSC_ERR_PLIB, "nloc %D != mat_offset[grid+1] - moffset = %D",nloc, a_mat_offset[grid+1] - moffset);
770         for (int i=0 ; i<nloc ; i++) {
771           ierr = MatGetRow(B,i,&nzl,&cols,&vals);CHKERRQ(ierr);
772           if (nzl>1024) SETERRQ1(PetscObjectComm((PetscObject) B), PETSC_ERR_PLIB, "Row too big: %D",nzl);
773           for (int j=0; j<nzl; j++) colbuf[j] = cols[j] + moffset;
774           row = i + moffset;
775           ierr = MatSetValues(JacP,1,&row,nzl,colbuf,vals,ADD_VALUES);CHKERRQ(ierr);
776           ierr = MatRestoreRow(B,i,&nzl,&cols,&vals);CHKERRQ(ierr);
777         }
778         ierr = MatDestroy(&subJ[grid]);CHKERRQ(ierr);
779       }
780       ierr = PetscLogEventEnd(events[6],0,0,0,0);CHKERRQ(ierr);
781     } // grids
782   }
783   PetscFunctionReturn(0);
784 }
785 } // extern "C"
786