xref: /petsc/src/ts/utils/dmplexlandau/kokkos/landau.kokkos.cxx (revision f155c23239e6d3f7c7ec79ff00b4f28519d0ce99)
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   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   PetscCheckFalse(LANDAU_DIM != dim,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 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
225     auto ipfdf_data =  new Kokkos::View<PetscReal***, Kokkos::LayoutLeft > ("fdf", batch_sz, dim+1, IPf_sz);
226 #else
227     auto ipfdf_data =  new Kokkos::View<PetscReal***, Kokkos::LayoutRight > ("fdf",batch_sz,  dim+1, IPf_sz);
228 #endif
229     SData_d->ipfdf_data = static_cast<void*>(ipfdf_data);
230     auto Eq_m = new Kokkos::View<PetscReal*, Kokkos::LayoutLeft> ("Eq_m",Nftot); // allocate but do not set, same for whole batch
231     SData_d->Eq_m = static_cast<void*>(Eq_m);
232     const Kokkos::View<LandauIdx*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_coo_elem_offsets (SData_d->coo_elem_offsets, SData_d->n_coo_cellsTot+1);
233     auto coo_elem_offsets = new Kokkos::View<LandauIdx*, Kokkos::LayoutLeft> ("coo_elem_offsets",SData_d->n_coo_cellsTot+1);
234     const Kokkos::View<LandauIdx*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_coo_elem_fullNb (SData_d->coo_elem_fullNb, SData_d->n_coo_cellsTot);
235     auto coo_elem_fullNb = new Kokkos::View<LandauIdx*, Kokkos::LayoutLeft> ("coo_elem_offsets",SData_d->n_coo_cellsTot);
236     const Kokkos::View<LandauIdx**, Kokkos::LayoutRight, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_coo_elem_point_offsets ((LandauIdx*)SData_d->coo_elem_point_offsets, SData_d->n_coo_cellsTot, LANDAU_MAX_NQ+1);
237     auto coo_elem_point_offsets = new Kokkos::View<LandauIdx**, Kokkos::LayoutRight> ("coo_elem_point_offsets",SData_d->n_coo_cellsTot,LANDAU_MAX_NQ+1);
238     // assign & copy
239     Kokkos::deep_copy (*coo_elem_offsets,  h_coo_elem_offsets);
240     Kokkos::deep_copy (*coo_elem_fullNb,  h_coo_elem_fullNb);
241     Kokkos::deep_copy (*coo_elem_point_offsets,  h_coo_elem_point_offsets);
242     // need to free this now and use pointer space
243     ierr = PetscFree3(SData_d->coo_elem_offsets,SData_d->coo_elem_fullNb,SData_d->coo_elem_point_offsets);CHKERRQ(ierr);
244     SData_d->coo_elem_offsets       = reinterpret_cast<LandauIdx*>(coo_elem_offsets);
245     SData_d->coo_elem_fullNb        = reinterpret_cast<LandauIdx*>(coo_elem_fullNb);
246     SData_d->coo_elem_point_offsets = reinterpret_cast<LandauIdx (*)[LANDAU_MAX_NQ+1]>(coo_elem_point_offsets);
247   }
248   SData_d->maps = NULL; // not used
249   PetscFunctionReturn(0);
250 }
251 
252 PetscErrorCode LandauKokkosStaticDataClear(LandauStaticData *SData_d)
253 {
254   PetscFunctionBegin;
255   if (SData_d->alpha) {
256     auto alpha = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->alpha);
257     delete alpha;
258     SData_d->alpha = NULL;
259     auto beta = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->beta);
260     delete beta;
261     auto invMass = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invMass);
262     delete invMass;
263     auto B = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->B);
264     delete B;
265     auto D = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->D);
266     delete D;
267     auto invJ = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invJ);
268     delete invJ;
269     auto x = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->x);
270     delete x;
271     auto y = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->y);
272     delete y;
273     if (SData_d->z) {
274       auto z = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->z);
275       delete z;
276     }
277 #if defined(LANDAU_LAYOUT_LEFT) // preallocate dynamic data, no copy
278     auto z = static_cast<Kokkos::View<PetscReal***, Kokkos::LayoutLeft>*>(SData_d->ipfdf_data);
279 #else
280     auto z = static_cast<Kokkos::View<PetscReal***, Kokkos::LayoutRight>*>(SData_d->ipfdf_data);
281 #endif
282     delete z;
283     auto w = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->w);
284     delete w;
285     auto Eq_m = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->Eq_m);
286     delete Eq_m;
287     // offset
288     auto nc = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->NCells);
289     delete nc;
290     auto soff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->species_offset);
291     delete soff;
292     auto moff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->mat_offset);
293     delete moff;
294     auto ipoff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ip_offset);
295     delete ipoff;
296     auto eoff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->elem_offset);
297     delete eoff;
298     auto ipfoff = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ipf_offset);
299     delete ipfoff;
300     auto coo_elem_offsets = static_cast<Kokkos::View<LandauIdx*, Kokkos::LayoutLeft>*>((void*)SData_d->coo_elem_offsets);
301     delete coo_elem_offsets;
302     auto coo_elem_fullNb = static_cast<Kokkos::View<LandauIdx*, Kokkos::LayoutLeft>*>((void*)SData_d->coo_elem_fullNb);
303     delete coo_elem_fullNb;
304     auto coo_elem_point_offsets = static_cast<Kokkos::View<LandauIdx**, Kokkos::LayoutRight>*>((void*)SData_d->coo_elem_point_offsets);
305     delete coo_elem_point_offsets;
306     SData_d->coo_elem_offsets = NULL;
307     SData_d->coo_elem_point_offsets = NULL;
308     SData_d->coo_elem_fullNb = NULL;
309   }
310   PetscFunctionReturn(0);
311 }
312 
313 #define KOKKOS_SHARED_LEVEL 1
314 KOKKOS_INLINE_FUNCTION
315 PetscErrorCode landau_mat_assemble(PetscSplitCSRDataStructure d_mat, PetscScalar *coo_vals, const PetscScalar Aij, const PetscInt f, const PetscInt g, const PetscInt Nb,
316                                    PetscInt moffset, const PetscInt elem, const PetscInt fieldA, const P4estVertexMaps *d_maps,
317                                    const LandauIdx coo_elem_offsets[], const LandauIdx coo_elem_fullNb[], Kokkos::View<LandauIdx**, Kokkos::LayoutRight> &coo_elem_point_offsets_k,
318                                    const PetscInt glb_elem_idx, const PetscInt bid_coo_sz_batch)
319 {
320   PetscInt               idx,q,nr,nc,d;
321   const LandauIdx *const Idxs = &d_maps->gIdx[elem][fieldA][0];
322   PetscScalar            row_scale[LANDAU_MAX_Q_FACE]={0},col_scale[LANDAU_MAX_Q_FACE]={0};
323   if (coo_vals) { // mirror (i,j) in CreateStaticGPUData
324     const int fullNb = coo_elem_fullNb[glb_elem_idx],fullNb2=fullNb*fullNb;
325     idx = Idxs[f];
326     if (idx >= 0) {
327       nr = 1;
328       row_scale[0] = 1.;
329     } else {
330       idx = -idx - 1;
331       for (q = 0, nr = 0; q < d_maps->num_face; q++, nr++) {
332         if (d_maps->c_maps[idx][q].gid < 0) break;
333         row_scale[q] = d_maps->c_maps[idx][q].scale;
334       }
335     }
336     idx = Idxs[g];
337     if (idx >= 0) {
338       nc = 1;
339       col_scale[0] = 1.;
340     } else {
341       idx = -idx - 1;
342       nc = d_maps->num_face;
343       for (q = 0, nc = 0; q < d_maps->num_face; q++, nc++) {
344         if (d_maps->c_maps[idx][q].gid < 0) break;
345         col_scale[q] = d_maps->c_maps[idx][q].scale;
346       }
347     }
348     const int idx0 = bid_coo_sz_batch + coo_elem_offsets[glb_elem_idx] + fieldA*fullNb2 + fullNb * coo_elem_point_offsets_k(glb_elem_idx,f) + nr * coo_elem_point_offsets_k(glb_elem_idx,g);
349     for (int q = 0, idx2 = idx0; q < nr; q++) {
350       for (int d = 0; d < nc; d++, idx2++) {
351         coo_vals[idx2] = row_scale[q]*col_scale[d]*Aij;
352       }
353     }
354   } else {
355     PetscScalar            vals[LANDAU_MAX_Q_FACE*LANDAU_MAX_Q_FACE]={0};
356     PetscInt               rows[LANDAU_MAX_Q_FACE],cols[LANDAU_MAX_Q_FACE];
357     idx = Idxs[f];
358     if (idx >= 0) {
359       nr = 1;
360       rows[0] = idx;
361       row_scale[0] = 1.;
362     } else {
363       idx = -idx - 1;
364       for (q = 0, nr = 0; q < d_maps->num_face; q++, nr++) {
365         if (d_maps->c_maps[idx][q].gid < 0) break;
366         rows[q]     = d_maps->c_maps[idx][q].gid;
367         row_scale[q] = d_maps->c_maps[idx][q].scale;
368       }
369     }
370     idx = Idxs[g];
371     if (idx >= 0) {
372       nc = 1;
373       cols[0] = idx;
374       col_scale[0] = 1.;
375     } else {
376       idx = -idx - 1;
377       for (q = 0, nc = 0; q < d_maps->num_face; q++, nc++) {
378         if (d_maps->c_maps[idx][q].gid < 0) break;
379         cols[q]     = d_maps->c_maps[idx][q].gid;
380         col_scale[q] = d_maps->c_maps[idx][q].scale;
381       }
382     }
383 
384     for (q = 0; q < nr; q++) rows[q] = rows[q] + moffset;
385     for (q = 0; q < nc; q++) cols[q] = cols[q] + moffset;
386     for (q = 0; q < nr; q++) {
387       for (d = 0; d < nc; d++) {
388         vals[q*nc + d] = row_scale[q]*col_scale[d]*Aij;
389       }
390     }
391     MatSetValuesDevice(d_mat,nr,rows,nc,cols,vals,ADD_VALUES);
392   }
393   return 0;
394 }
395 
396 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[],
397                                     const PetscScalar a_xarray[], const LandauStaticData *SData_d, const PetscReal shift,
398                                     const PetscLogEvent events[], const PetscInt a_mat_offset[], const PetscInt a_species_offset[], Mat subJ[], Mat JacP)
399 {
400   using scr_mem_t = Kokkos::DefaultExecutionSpace::scratch_memory_space;
401   using real2_scr_t = Kokkos::View<PetscScalar**, Kokkos::LayoutRight, scr_mem_t>;
402   using g2_scr_t = Kokkos::View<PetscReal***, Kokkos::LayoutRight, scr_mem_t>;
403   using g3_scr_t = Kokkos::View<PetscReal****, Kokkos::LayoutRight, scr_mem_t>;
404   PetscErrorCode    ierr;
405   PetscInt          Nb=Nq,dim,num_cells_max,Nf_max,num_cells_batch;
406   int               nfaces=0;
407   LandauCtx         *ctx;
408   PetscReal         *d_Eq_m=NULL;
409   PetscScalar       *d_vertex_f=NULL;
410   P4estVertexMaps   *maps[LANDAU_MAX_GRIDS]; // this gets captured
411   PetscSplitCSRDataStructure d_mat;
412   PetscContainer    container;
413   const int         conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp==0) ? Nq : 1;
414   const PetscInt    coo_sz_batch = SData_d->coo_size/batch_sz; // capture
415   auto              d_alpha_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->alpha); //static data
416   const PetscReal   *d_alpha = d_alpha_k->data();
417   const PetscInt    Nftot = d_alpha_k->size(); // total number of species
418   auto              d_beta_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->beta);
419   const PetscReal   *d_beta = d_beta_k->data();
420   auto              d_invMass_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invMass);
421   const PetscReal   *d_invMass = d_invMass_k->data();
422   auto              d_B_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->B);
423   const PetscReal   *d_BB = d_B_k->data();
424   auto              d_D_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->D);
425   const PetscReal   *d_DD = d_D_k->data();
426   auto              d_invJ_k = *static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->invJ);     // use Kokkos vector in kernels
427   auto              d_x_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->x); //static data
428   const PetscReal   *d_x = d_x_k->data();
429   auto              d_y_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->y); //static data
430   const PetscReal   *d_y = d_y_k->data();
431   auto              d_z_k = static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->z); //static data
432   const PetscReal   *d_z = (LANDAU_DIM==3) ? d_z_k->data() : NULL;
433   auto              d_w_k = *static_cast<Kokkos::View<PetscReal*, Kokkos::LayoutLeft>*>(SData_d->w); //static data
434   const PetscReal   *d_w = d_w_k.data();
435   // grid offsets - single vertex grid data
436   auto              d_numCells_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->NCells);
437   const PetscInt   *d_numCells = d_numCells_k->data();
438   auto              d_species_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->species_offset);
439   const PetscInt   *d_species_offset = d_species_offset_k->data();
440   auto              d_mat_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->mat_offset);
441   const PetscInt   *d_mat_offset = d_mat_offset_k->data();
442   auto              d_ip_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ip_offset);
443   const PetscInt   *d_ip_offset = d_ip_offset_k->data();
444   auto              d_ipf_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->ipf_offset);
445   const PetscInt   *d_ipf_offset = d_ipf_offset_k->data();
446   auto              d_elem_offset_k = static_cast<Kokkos::View<PetscInt*, Kokkos::LayoutLeft>*>(SData_d->elem_offset);
447   const PetscInt   *d_elem_offset = d_elem_offset_k->data();
448 #if defined(LANDAU_LAYOUT_LEFT)         // preallocate dynamic data, including batched vertices
449   Kokkos::View<PetscReal***, Kokkos::LayoutLeft > d_fdf_k = *static_cast<Kokkos::View<PetscReal***, Kokkos::LayoutLeft >*>(SData_d->ipfdf_data);
450 #else
451   Kokkos::View<PetscReal***, Kokkos::LayoutRight > d_fdf_k = *static_cast<Kokkos::View<PetscReal***, Kokkos::LayoutRight >*>(SData_d->ipfdf_data);
452 #endif
453   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
454   // COO
455   auto d_coo_elem_offsets_k = static_cast<Kokkos::View<LandauIdx*, Kokkos::LayoutLeft>*>((void*)SData_d->coo_elem_offsets);
456   auto d_coo_elem_fullNb_k = static_cast<Kokkos::View<LandauIdx*, Kokkos::LayoutLeft>*>((void*)SData_d->coo_elem_fullNb);
457   auto d_coo_elem_point_offsets_k = static_cast<Kokkos::View<LandauIdx**, Kokkos::LayoutRight>*>((void*)SData_d->coo_elem_point_offsets);
458   Kokkos::View<PetscScalar*, Kokkos::LayoutRight,Kokkos::DefaultExecutionSpace> d_coo_vals_k("coo_vals", SData_d->coo_size); // device data (default space)
459   PetscScalar* d_coo_vals = (SData_d->coo_size==0) ? NULL : d_coo_vals_k.data();
460 
461   PetscFunctionBegin;
462   ierr = PetscLogEventBegin(events[3],0,0,0,0);CHKERRQ(ierr);
463   ierr = DMGetApplicationContext(plex[0], &ctx);CHKERRQ(ierr);
464   PetscCheckFalse(!ctx,PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
465   ierr = DMGetDimension(plex[0], &dim);CHKERRQ(ierr);
466   PetscCheckFalse(LANDAU_DIM != dim,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM);
467   if (ctx->gpu_assembly) {
468     ierr = PetscObjectQuery((PetscObject) JacP, "assembly_maps", (PetscObject *) &container);CHKERRQ(ierr);
469     if (container) {
470       P4estVertexMaps   *h_maps=NULL;
471       ierr = PetscContainerGetPointer(container, (void **) &h_maps);CHKERRQ(ierr);
472       for (PetscInt grid=0 ; grid<num_grids ; grid++) {
473         if (h_maps[grid].d_self) {
474           maps[grid] = h_maps[grid].d_self;
475           nfaces = h_maps[grid].num_face; // nface=0 for CPU assembly
476         } else {
477           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container");
478         }
479       }
480       if (!d_coo_vals) {
481         // this does the setup the first time called
482         ierr = MatKokkosGetDeviceMatWrite(JacP,&d_mat);CHKERRQ(ierr);
483       } else {
484         d_mat = NULL;
485       }
486     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "GPU assembly but no metadata in container");
487   } else {
488     for (PetscInt grid=0 ; grid<num_grids ; grid++) maps[grid] = NULL;
489     nfaces = 0;
490     d_mat = NULL;
491     container = NULL;
492   }
493   num_cells_batch = Nf_max = num_cells_max = 0;
494   for (PetscInt grid=0 ; grid<num_grids ; grid++) {
495     int Nfloc = a_species_offset[grid+1] - a_species_offset[grid];
496     if (Nfloc > Nf_max) Nf_max = Nfloc;
497     if (a_numCells[grid] > num_cells_max) num_cells_max = a_numCells[grid];
498     num_cells_batch += a_numCells[grid]; // we don't have a host element offset here (but in ctx)
499   }
500   const PetscInt totDim_max = Nf_max*Nq, elem_mat_size_max = totDim_max*totDim_max;
501   const PetscInt elem_mat_num_cells_max_grid = container ? 0 : num_cells_max;
502   Kokkos::View<PetscScalar****, Kokkos::LayoutRight> d_elem_mats("element matrices", batch_sz, num_grids, elem_mat_num_cells_max_grid, elem_mat_size_max); // first call have large set of global (Jac) element matrices
503   const Kokkos::View<PetscReal*, Kokkos::LayoutLeft, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> >  h_Eq_m_k (a_Eq_m, Nftot);
504   if (a_elem_closure || a_xarray) {
505     Kokkos::deep_copy (*d_Eq_m_k, h_Eq_m_k);
506     d_Eq_m = d_Eq_m_k->data();
507   } else d_Eq_m = NULL;
508   ierr = PetscKokkosInitializeCheck();CHKERRQ(ierr);
509   ierr = PetscLogEventEnd(events[3],0,0,0,0);CHKERRQ(ierr);
510   if (a_elem_closure || a_xarray) { // Jacobian, create f & df
511     Kokkos::View<PetscScalar*, Kokkos::LayoutLeft> *d_vertex_f_k = NULL;
512     ierr = PetscLogEventBegin(events[1],0,0,0,0);CHKERRQ(ierr);
513     if (a_elem_closure) {
514       PetscInt closure_sz = 0; // argh, don't have this on the host!!!
515       for (PetscInt grid=0 ; grid<num_grids ; grid++) {
516         PetscInt nfloc = a_species_offset[grid+1] - a_species_offset[grid];
517         closure_sz     += Nq*nfloc*a_numCells[grid];
518       }
519       closure_sz *= batch_sz;
520       d_vertex_f_k = new Kokkos::View<PetscScalar*, Kokkos::LayoutLeft> ("closure",closure_sz);
521       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
522       Kokkos::deep_copy (*d_vertex_f_k, h_closure_k);
523       d_vertex_f  = d_vertex_f_k->data();
524     } else {
525       d_vertex_f = (PetscScalar*)a_xarray;
526     }
527     ierr = PetscLogEventEnd(events[1],0,0,0,0);CHKERRQ(ierr);
528     ierr = PetscLogEventBegin(events[8],0,0,0,0);CHKERRQ(ierr);
529     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
530 
531     const int scr_bytes_fdf = real2_scr_t::shmem_size(Nf_max,Nb);
532     ierr = PetscInfo(plex[0], "Jacobian shared memory size: %d bytes in level %d num cells total=%D team size=%D #face=%D Nf_max=%D\n",scr_bytes_fdf,KOKKOS_SHARED_LEVEL,num_cells_batch*batch_sz,team_size,nfaces,Nf_max);CHKERRQ(ierr);
533     Kokkos::parallel_for("f, df", Kokkos::TeamPolicy<>(num_cells_batch*batch_sz, team_size, /* Kokkos::AUTO */ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes_fdf)), KOKKOS_LAMBDA (const team_member team) {
534         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];
535         // find my grid
536         PetscInt grid = 0;
537         while (b_elem_idx >= d_elem_offset[grid+1]) grid++;
538         {
539           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];
540           const PetscInt     moffset = LAND_MOFFSET(b_id,grid,batch_sz,num_grids,d_mat_offset);
541           {
542             real2_scr_t    s_coef_k(team.team_scratch(KOKKOS_SHARED_LEVEL),Nf_max,Nb);
543             PetscScalar     *coef;
544             const PetscReal *invJe = &d_invJ_k((d_ip_offset[grid] + loc_elem*Nq)*dim*dim);
545             // un pack IPData
546             if (!maps[grid]) {
547               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
548             } else {
549               coef = s_coef_k.data();
550               Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,(int)loc_Nf), [=] (int f) {
551                   //for (int f = 0; f < loc_Nf; ++f) {
552                   Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,(int)Nb), [=] (int b) {
553                       //for (int b = 0; b < Nb; ++b) {
554                       LandauIdx idx = maps[grid]->gIdx[loc_elem][f][b];
555                       if (idx >= 0) {
556                         coef[f*Nb+b] = d_vertex_f[idx+moffset]; // xarray data, not IP, need raw array to deal with CPU assemble case (not used)
557                       } else {
558                         idx = -idx - 1;
559                         coef[f*Nb+b] = 0;
560                         for (int q = 0; q < maps[grid]->num_face; q++) {
561                           PetscInt    id = maps[grid]->c_maps[idx][q].gid;
562                           PetscScalar scale = maps[grid]->c_maps[idx][q].scale;
563                           if (id >= 0) coef[f*Nb+b] += scale*d_vertex_f[id+moffset];
564                         }
565                       }
566                     });
567                 });
568             }
569             team.team_barrier();
570             Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nq), [=] (int myQi) {
571                 const PetscReal *const  invJ = &invJe[myQi*dim*dim]; // b_elem_idx: batch element index
572                 const PetscReal         *Bq = &d_BB[myQi*Nb], *Dq = &d_DD[myQi*Nb*dim];
573                 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,(int)loc_Nf), [=] (int f) {
574                     PetscInt   b, e, d;
575                     PetscReal  refSpaceDer[LANDAU_DIM];
576                     const PetscInt idx = d_ipf_offset[grid] + f*loc_nip + loc_elem*Nq + myQi;
577                     d_fdf_k(b_id,0,idx) = 0.0;
578                     for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0;
579                     for (b = 0; b < Nb; ++b) {
580                       d_fdf_k(b_id,0,idx) += Bq[b]*PetscRealPart(coef[f*Nb+b]);
581                       for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[b*dim+d]*PetscRealPart(coef[f*Nb+b]);
582                     }
583                     for (d = 0; d < dim; ++d) {
584                       for (e = 0, d_fdf_k(b_id,d+1,idx) = 0.0; e < dim; ++e) {
585                         d_fdf_k(b_id,d+1,idx) += invJ[e*dim+d]*refSpaceDer[e];
586                       }
587                     }
588                   }); // f
589               }); // myQi
590           }
591           team.team_barrier();
592         } // 'grid' scope
593       }); // global elems - fdf
594     Kokkos::fence();
595     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr); // is this a fence?
596     ierr = PetscLogEventEnd(events[8],0,0,0,0);CHKERRQ(ierr);
597     // Jacobian
598     auto jac_lambda = KOKKOS_LAMBDA (const team_member team) {
599       const PetscInt b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank()%b_Nelem, b_id = team.league_rank()/b_Nelem;
600       // find my grid
601       PetscInt grid = 0;
602       while (b_elem_idx >= d_elem_offset[grid+1]) grid++;
603       {
604         const PetscInt  loc_Nf = d_species_offset[grid+1]-d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid];
605         const PetscInt  moffset = LAND_MOFFSET(b_id,grid,batch_sz,num_grids,d_mat_offset);
606         const PetscInt  f_off = d_species_offset[grid], totDim = loc_Nf*Nq;
607         g2_scr_t        g2(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,loc_Nf,Nq);
608         g3_scr_t        g3(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,dim,loc_Nf,Nq);
609         g2_scr_t        gg2(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,loc_Nf,Nq);
610         g3_scr_t        gg3(team.team_scratch(KOKKOS_SHARED_LEVEL),dim,dim,loc_Nf,Nq);
611         // get g2[] & g3[]
612         Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nq), [=] (int myQi) {
613             using Kokkos::parallel_reduce;
614             const PetscInt  jpidx_glb = d_ip_offset[grid] + loc_elem * Nq + myQi;
615             const PetscReal *invJ = &d_invJ_k(jpidx_glb*dim*dim);
616             const PetscReal vj[3] = {d_x[jpidx_glb], d_y[jpidx_glb], d_z ? d_z[jpidx_glb] : 0}, wj = d_w[jpidx_glb];
617             landau_inner_red::TensorValueType gg_temp; // reduce on part of gg2 and g33 for IP jpidx_g
618             Kokkos::parallel_reduce(Kokkos::ThreadVectorRange (team, (int)d_ip_offset[num_grids]), [=] (const int& ipidx, landau_inner_red::TensorValueType & ggg) {
619                 const PetscReal wi = d_w[ipidx], x = d_x[ipidx], y = d_y[ipidx];
620                 PetscReal       temp1[3] = {0, 0, 0}, temp2 = 0;
621                 PetscInt        fieldA,d2,d3,f_off_r,grid_r,ipidx_g,nip_loc_r,loc_Nf_r;
622 #if LANDAU_DIM==2
623                 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.;
624                 LandauTensor2D(vj, x, y, Ud, Uk, mask);
625 #else
626                 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.;
627                 LandauTensor3D(vj, x, y, z, U, mask);
628 #endif
629                 grid_r = 0;
630                 while (ipidx >= d_ip_offset[grid_r+1]) grid_r++; // yuck search for grid
631                 f_off_r = d_species_offset[grid_r];
632                 ipidx_g = ipidx - d_ip_offset[grid_r];
633                 nip_loc_r = d_numCells[grid_r]*Nq;
634                 loc_Nf_r = d_species_offset[grid_r+1] - d_species_offset[grid_r];
635                 for (fieldA = 0; fieldA < loc_Nf_r; ++fieldA) {
636                   const PetscInt idx = d_ipf_offset[grid_r] + fieldA*nip_loc_r + ipidx_g;
637                   temp1[0] += d_fdf_k(b_id,1,idx)*d_beta[fieldA+f_off_r]*d_invMass[fieldA+f_off_r];
638                   temp1[1] += d_fdf_k(b_id,2,idx)*d_beta[fieldA+f_off_r]*d_invMass[fieldA+f_off_r];
639 #if LANDAU_DIM==3
640                   temp1[2] += d_fdf_k(b_id,3,idx)*d_beta[fieldA+f_off_r]*d_invMass[fieldA+f_off_r];
641 #endif
642                   temp2    += d_fdf_k(b_id,0,idx)*d_beta[fieldA+f_off_r];
643                 }
644                 temp1[0] *= wi;
645                 temp1[1] *= wi;
646 #if LANDAU_DIM==3
647                 temp1[2] *= wi;
648 #endif
649                 temp2    *= wi;
650 #if LANDAU_DIM==2
651                 for (d2 = 0; d2 < 2; d2++) {
652                   for (d3 = 0; d3 < 2; ++d3) {
653                     /* K = U * grad(f): g2=e: i,A */
654                     ggg.gg2[d2] += Uk[d2][d3]*temp1[d3];
655                     /* D = -U * (I \kron (fx)): g3=f: i,j,A */
656                     ggg.gg3[d2][d3] += Ud[d2][d3]*temp2;
657                   }
658                 }
659 #else
660                 for (d2 = 0; d2 < 3; ++d2) {
661                   for (d3 = 0; d3 < 3; ++d3) {
662                     /* K = U * grad(f): g2 = e: i,A */
663                     ggg.gg2[d2] += U[d2][d3]*temp1[d3];
664                     /* D = -U * (I \kron (fx)): g3 = f: i,j,A */
665                     ggg.gg3[d2][d3] += U[d2][d3]*temp2;
666                   }
667                 }
668 #endif
669               }, Kokkos::Sum<landau_inner_red::TensorValueType>(gg_temp));
670             // add alpha and put in gg2/3
671             Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)loc_Nf), [&] (const int& fieldA) {
672                 PetscInt d2,d3;
673                 for (d2 = 0; d2 < dim; d2++) {
674                   gg2(d2,fieldA,myQi) = gg_temp.gg2[d2]*d_alpha[fieldA+f_off];
675                   for (d3 = 0; d3 < dim; d3++) {
676                     gg3(d2,d3,fieldA,myQi) = -gg_temp.gg3[d2][d3]*d_alpha[fieldA+f_off]*d_invMass[fieldA+f_off];
677                   }
678                 }
679               });
680             /* add electric field term once per IP */
681             Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)loc_Nf), [&] (const int& fieldA) {
682                 gg2(dim-1,fieldA,myQi) += d_Eq_m[fieldA+f_off];
683               });
684             Kokkos::parallel_for(Kokkos::ThreadVectorRange (team, (int)loc_Nf), [=] (const int& fieldA) {
685                 int d,d2,d3,dp;
686                 /* Jacobian transform - g2, g3 - per thread (2D) */
687                 for (d = 0; d < dim; ++d) {
688                   g2(d,fieldA,myQi) = 0;
689                   for (d2 = 0; d2 < dim; ++d2) {
690                     g2(d,fieldA,myQi) += invJ[d*dim+d2]*gg2(d2,fieldA,myQi);
691                     g3(d,d2,fieldA,myQi) = 0;
692                     for (d3 = 0; d3 < dim; ++d3) {
693                       for (dp = 0; dp < dim; ++dp) {
694                         g3(d,d2,fieldA,myQi) += invJ[d*dim + d3]*gg3(d3,dp,fieldA,myQi)*invJ[d2*dim + dp];
695                       }
696                     }
697                     g3(d,d2,fieldA,myQi) *= wj;
698                   }
699                   g2(d,fieldA,myQi) *= wj;
700                 }
701               });
702           }); // Nq
703         team.team_barrier();
704         { /* assemble */
705           for (PetscInt fieldA = 0; fieldA < loc_Nf; fieldA++) {
706             /* assemble */
707             Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) {
708                 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) {
709                     PetscScalar t = 0;
710                     for (int qj = 0 ; qj < Nq ; qj++) { // look at others integration points
711                       const PetscReal *BJq = &d_BB[qj*Nb], *DIq = &d_DD[qj*Nb*dim];
712                       for (int d = 0; d < dim; ++d) {
713                         t += DIq[f*dim+d]*g2(d,fieldA,qj)*BJq[g];
714                         for (int d2 = 0; d2 < dim; ++d2) {
715                           t += DIq[f*dim + d]*g3(d,d2,fieldA,qj)*DIq[g*dim + d2];
716                         }
717                       }
718                     }
719                     if (elem_mat_num_cells_max_grid) { // CPU assembly
720                       const PetscInt fOff = (fieldA*Nb + f)*totDim + fieldA*Nb + g;
721                       d_elem_mats(b_id,grid,loc_elem,fOff) = t;
722                     } else {
723                       landau_mat_assemble (d_mat, d_coo_vals, t, f, g, Nb, moffset, loc_elem, fieldA, maps[grid], d_coo_elem_offsets_k->data(), d_coo_elem_fullNb_k->data(), *d_coo_elem_point_offsets_k, b_elem_idx, b_id*coo_sz_batch);
724                     }
725                   });
726               });
727           }
728         }
729       } // scope with 'grid'
730     };
731     ierr = PetscLogEventBegin(events[4],0,0,0,0);CHKERRQ(ierr);
732     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
733     const int scr_bytes = 2*(g2_scr_t::shmem_size(dim,Nf_max,Nq) + g3_scr_t::shmem_size(dim,dim,Nf_max,Nq));
734     Kokkos::parallel_for("Jacobian", Kokkos::TeamPolicy<Kokkos::LaunchBounds<512,1>>(num_cells_batch*batch_sz, team_size, /* Kokkos::AUTO */ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes)), jac_lambda);
735     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
736     ierr = PetscLogEventEnd(events[4],0,0,0,0);CHKERRQ(ierr);
737     if (d_vertex_f_k) delete d_vertex_f_k;
738   } else { // mass - 2*Nb is guess at max size (2D Q3 is 31 =< 2*Nb = 32)
739     using fieldMats_scr_t = Kokkos::View<PetscScalar*, Kokkos::LayoutRight, scr_mem_t>;
740     PetscInt loc_ass_sz = 1; // captured
741     ierr = PetscLogEventBegin(events[16],0,0,0,0);CHKERRQ(ierr);
742     ierr = PetscLogGpuTimeBegin();CHKERRQ(ierr);
743     if (loc_ass_sz) loc_ass_sz = ctx->SData_d.coo_max_fullnb*ctx->SData_d.coo_max_fullnb;
744     int scr_bytes = fieldMats_scr_t::shmem_size(loc_ass_sz); // + idx_scr_t::shmem_size(Nb,nfaces) + scale_scr_t::shmem_size(Nb,nfaces);
745     ierr = PetscInfo(plex[0], "Mass shared memory size: %d bytes in level %d conc=%D team size=%D #face=%D Nb=%D, %s assembly\n",scr_bytes,KOKKOS_SHARED_LEVEL,conc,team_size,nfaces,Nb, d_coo_vals ? (loc_ass_sz==0 ? "COO" : "optimized COO") : "CSR");CHKERRQ(ierr);
746     Kokkos::parallel_for("Mass", Kokkos::TeamPolicy<>(num_cells_batch*batch_sz, team_size, /* Kokkos::AUTO */ 16).set_scratch_size(KOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes)), KOKKOS_LAMBDA (const team_member team) {
747         fieldMats_scr_t s_fieldMats(team.team_scratch(KOKKOS_SHARED_LEVEL),loc_ass_sz); // Only used for GPU assembly (ie, not first pass)
748         const PetscInt  b_Nelem = d_elem_offset[num_grids], b_elem_idx = team.league_rank()%b_Nelem, b_id = team.league_rank()/b_Nelem;
749         // find my grid
750         PetscInt grid = 0;
751         while (b_elem_idx >= d_elem_offset[grid+1]) grid++;
752         {
753           const PetscInt loc_Nf = d_species_offset[grid+1]-d_species_offset[grid], loc_elem = b_elem_idx - d_elem_offset[grid], totDim = loc_Nf*Nq, jpidx_0 = d_ip_offset[grid] + loc_elem * Nq;
754           const PetscInt moffset = LAND_MOFFSET(b_id,grid,batch_sz,num_grids,d_mat_offset);
755           for (int fieldA = 0; fieldA < loc_Nf; fieldA++) {
756             /* assemble */
757             Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,Nb), [=] (int f) {
758                 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,Nb), [=] (int g) {
759                     PetscScalar    t = 0;
760                     for (int qj = 0 ; qj < Nq ; qj++) { // look at others integration points
761                       const PetscReal *BJq = &d_BB[qj*Nb];
762                       const PetscInt jpidx_glb = jpidx_0 + qj;
763                       if (dim==2) {
764                         t += BJq[f] * d_w_k(jpidx_glb) * shift * BJq[g] * 2. * PETSC_PI;
765                       } else {
766                         t += BJq[f] * d_w_k(jpidx_glb) * shift * BJq[g];
767                       }
768                     }
769                     if (elem_mat_num_cells_max_grid) {
770                       const PetscInt fOff = (fieldA*Nb + f)*totDim + fieldA*Nb + g;
771                       d_elem_mats(b_id,grid,loc_elem,fOff) = t;
772                     } else {
773                       if (d_coo_vals && loc_ass_sz) {
774                         PetscInt               idx,q,nr,nc;
775                         const LandauIdx *const Idxs = &maps[grid]->gIdx[loc_elem][fieldA][0];
776                         PetscScalar            row_scale[LANDAU_MAX_Q_FACE]={0},col_scale[LANDAU_MAX_Q_FACE]={0};
777                         const int              fullNb = (*d_coo_elem_fullNb_k)(b_elem_idx), num_face = maps[grid]->num_face;
778                         const PetscScalar      Aij = t;
779                         idx = Idxs[f];
780                         if (idx >= 0) {
781                           nr = 1;
782                           row_scale[0] = 1.;
783                         } else {
784                           idx = -idx - 1;
785                           for (q = 0, nr = 0; q < num_face; q++, nr++) {
786                             if (maps[grid]->c_maps[idx][q].gid < 0) break;
787                             row_scale[q] = maps[grid]->c_maps[idx][q].scale;
788                           }
789                         }
790                         idx = Idxs[g];
791                         if (idx >= 0) {
792                           nc = 1;
793                           col_scale[0] = 1.;
794                         } else {
795                           idx = -idx - 1;
796                           for (q = 0, nc = 0; q < num_face; q++, nc++) {
797                             if (maps[grid]->c_maps[idx][q].gid < 0) break;
798                             col_scale[q] = maps[grid]->c_maps[idx][q].scale;
799                           }
800                         }
801                         const int idx0 = fullNb * (*d_coo_elem_point_offsets_k)(b_elem_idx,f) + nr * (*d_coo_elem_point_offsets_k)(b_elem_idx,g);
802                         for (int q = 0, idx2 = idx0; q < nr; q++) {
803                           for (int d = 0; d < nc; d++, idx2++) {
804                             s_fieldMats(idx2) = row_scale[q]*col_scale[d]*Aij;
805                           }
806                         }
807                       } else {
808                         landau_mat_assemble (d_mat, d_coo_vals, t, f, g, Nb, moffset, loc_elem, fieldA, maps[grid], d_coo_elem_offsets_k->data(), d_coo_elem_fullNb_k->data(), *d_coo_elem_point_offsets_k, b_elem_idx, b_id*coo_sz_batch);
809                       }
810                     }
811                   });
812               });
813             if (!elem_mat_num_cells_max_grid) {
814               if (d_coo_vals && loc_ass_sz) {
815                 const int fullNb = (*d_coo_elem_fullNb_k)(b_elem_idx), fullNb2=fullNb*fullNb;
816                 const int idx0 =  b_id*coo_sz_batch + (*d_coo_elem_offsets_k)(b_elem_idx) + fieldA*fullNb2;
817                 team.team_barrier();
818 #if 0
819                 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,0,fullNb), [=] (int gg) {
820                   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,fullNb), [=] (int ff) {
821                         const int idx = fullNb * ff + gg;
822                         d_coo_vals[idx0 + idx] = s_fieldMats(idx);
823                       });
824                   });
825 #else
826                 Kokkos::parallel_for(Kokkos::TeamVectorRange(team,0,fullNb2), [=] (int idx) { d_coo_vals[idx0 + idx] = s_fieldMats(idx); });
827 #endif
828               }
829             }
830           } // field
831         } // grid
832       });
833     ierr = PetscLogGpuTimeEnd();CHKERRQ(ierr);
834     ierr = PetscLogEventEnd(events[16],0,0,0,0);CHKERRQ(ierr);
835   }
836   Kokkos::fence();
837   if (d_coo_vals) {
838     ierr = MatSetValuesCOO(JacP,d_coo_vals,ADD_VALUES);CHKERRQ(ierr);
839   } else if (elem_mat_num_cells_max_grid) { // CPU assembly
840     Kokkos::View<PetscScalar****, Kokkos::LayoutRight>::HostMirror h_elem_mats = Kokkos::create_mirror_view(d_elem_mats);
841     Kokkos::deep_copy (h_elem_mats, d_elem_mats);
842      PetscCheck(!container,PETSC_COMM_SELF, PETSC_ERR_PLIB, "?????");
843     for (PetscInt b_id = 0 ; b_id < batch_sz ; b_id++) { // OpenMP (once)
844       for (PetscInt grid=0 ; grid<num_grids ; grid++) {
845         PetscSection      section, globalSection;
846         PetscInt          cStart,cEnd;
847         Mat               B = subJ[ LAND_PACK_IDX(b_id,grid) ];
848         PetscInt          moffset = LAND_MOFFSET(b_id,grid,batch_sz,num_grids,a_mat_offset), nloc, nzl, colbuf[1024], row;
849         const PetscInt    *cols;
850         const PetscScalar *vals;
851         ierr = PetscLogEventBegin(events[5],0,0,0,0);CHKERRQ(ierr);
852         ierr = DMPlexGetHeightStratum(plex[grid],0,&cStart,&cEnd);CHKERRQ(ierr);
853         ierr = DMGetLocalSection(plex[grid], &section);CHKERRQ(ierr);
854         ierr = DMGetGlobalSection(plex[grid], &globalSection);CHKERRQ(ierr);
855         ierr = PetscLogEventEnd(events[5],0,0,0,0);CHKERRQ(ierr);
856         ierr = PetscLogEventBegin(events[6],0,0,0,0);CHKERRQ(ierr);
857         for (PetscInt ej = cStart ; ej < cEnd; ++ej) {
858           const PetscScalar *elMat = &h_elem_mats(b_id,grid,ej-cStart,0);
859           ierr = DMPlexMatSetClosure(plex[grid], section, globalSection, B, ej, elMat, ADD_VALUES);CHKERRQ(ierr);
860           if (grid==0 && ej==-1) {
861             const PetscInt  loc_Nf = a_species_offset[grid+1]-a_species_offset[grid], totDim = loc_Nf*Nq;
862             int d,f;
863             PetscPrintf(PETSC_COMM_SELF,"Kokkos Element matrix %d/%d\n",1,(int)a_numCells[grid]);
864             for (d = 0; d < totDim; ++d){
865               for (f = 0; f < totDim; ++f) PetscPrintf(PETSC_COMM_SELF," %12.5e",  PetscRealPart(elMat[d*totDim + f]));
866               PetscPrintf(PETSC_COMM_SELF,"\n");
867             }
868             exit(14);
869           }
870         }
871         // move nest matrix to global JacP
872         ierr = MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
873         ierr = MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
874         ierr = MatGetSize(B, &nloc, NULL);CHKERRQ(ierr);
875         for (int i=0 ; i<nloc ; i++) {
876           ierr = MatGetRow(B,i,&nzl,&cols,&vals);CHKERRQ(ierr);
877           PetscCheckFalse(nzl>1024,PetscObjectComm((PetscObject) B), PETSC_ERR_PLIB, "Row too big: %D",nzl);
878           for (int j=0; j<nzl; j++) colbuf[j] = cols[j] + moffset;
879           row = i + moffset;
880           ierr = MatSetValues(JacP,1,&row,nzl,colbuf,vals,ADD_VALUES);CHKERRQ(ierr);
881           ierr = MatRestoreRow(B,i,&nzl,&cols,&vals);CHKERRQ(ierr);
882         }
883         ierr = MatDestroy(&subJ[ LAND_PACK_IDX(b_id,grid) ]);CHKERRQ(ierr);
884         ierr = PetscLogEventEnd(events[6],0,0,0,0);CHKERRQ(ierr);
885       } // grids
886     }
887   }
888   PetscFunctionReturn(0);
889 }
890 } // extern "C"
891