xref: /petsc/src/ts/utils/dmplexlandau/plexland.c (revision 9d0448ce9997fbe84686ae47a64f2ca347e5bb62)
1 #include <../src/mat/impls/aij/seq/aij.h>
2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
3 #include <petsclandau.h>              /*I "petsclandau.h"   I*/
4 #include <petscts.h>
5 #include <petscdmforest.h>
6 #include <petscdmcomposite.h>
7 
8 /* Landau collision operator */
9 
10 /* relativistic terms */
11 #if defined(PETSC_USE_REAL_SINGLE)
12   #define SPEED_OF_LIGHT 2.99792458e8F
13   #define C_0(v0)        (SPEED_OF_LIGHT / v0) /* needed for relativistic tensor on all architectures */
14 #else
15   #define SPEED_OF_LIGHT 2.99792458e8
16   #define C_0(v0)        (SPEED_OF_LIGHT / v0) /* needed for relativistic tensor on all architectures */
17 #endif
18 
19 #include "land_tensors.h"
20 
21 #if defined(PETSC_HAVE_OPENMP)
22   #include <omp.h>
23 #endif
24 
25 static PetscErrorCode LandauGPUMapsDestroy(void *ptr)
26 {
27   P4estVertexMaps *maps = (P4estVertexMaps *)ptr;
28   PetscFunctionBegin;
29   // free device data
30   if (maps[0].deviceType != LANDAU_CPU) {
31 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
32     if (maps[0].deviceType == LANDAU_KOKKOS) {
33       PetscCall(LandauKokkosDestroyMatMaps(maps, maps[0].numgrids)); // implies Kokkos does
34     }                                                                // else could be CUDA
35 #elif defined(PETSC_HAVE_CUDA)
36     if (maps[0].deviceType == LANDAU_CUDA) {
37       PetscCall(LandauCUDADestroyMatMaps(maps, maps[0].numgrids));
38     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps->deviceType %d ?????", maps->deviceType);
39 #endif
40   }
41   // free host data
42   for (PetscInt grid = 0; grid < maps[0].numgrids; grid++) {
43     PetscCall(PetscFree(maps[grid].c_maps));
44     PetscCall(PetscFree(maps[grid].gIdx));
45   }
46   PetscCall(PetscFree(maps));
47 
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 static PetscErrorCode energy_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
51 {
52   PetscReal v2 = 0;
53   PetscFunctionBegin;
54   /* compute v^2 / 2 */
55   for (int i = 0; i < dim; ++i) v2 += x[i] * x[i];
56   /* evaluate the Maxwellian */
57   u[0] = v2 / 2;
58   PetscFunctionReturn(PETSC_SUCCESS);
59 }
60 
61 /* needs double */
62 static PetscErrorCode gamma_m1_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
63 {
64   PetscReal *c2_0_arr = ((PetscReal *)actx);
65   double     u2 = 0, c02 = (double)*c2_0_arr, xx;
66 
67   PetscFunctionBegin;
68   /* compute u^2 / 2 */
69   for (int i = 0; i < dim; ++i) u2 += x[i] * x[i];
70   /* gamma - 1 = g_eps, for conditioning and we only take derivatives */
71   xx = u2 / c02;
72 #if defined(PETSC_USE_DEBUG)
73   u[0] = PetscSqrtReal(1. + xx);
74 #else
75   u[0] = xx / (PetscSqrtReal(1. + xx) + 1.) - 1.; // better conditioned. -1 might help condition and only used for derivative
76 #endif
77   PetscFunctionReturn(PETSC_SUCCESS);
78 }
79 
80 /*
81  LandauFormJacobian_Internal - Evaluates Jacobian matrix.
82 
83  Input Parameters:
84  .  globX - input vector
85  .  actx - optional user-defined context
86  .  dim - dimension
87 
88  Output Parameter:
89  .  J0acP - Jacobian matrix filled, not created
90  */
91 static PetscErrorCode LandauFormJacobian_Internal(Vec a_X, Mat JacP, const PetscInt dim, PetscReal shift, void *a_ctx)
92 {
93   LandauCtx         *ctx = (LandauCtx *)a_ctx;
94   PetscInt           numCells[LANDAU_MAX_GRIDS], Nq, Nb;
95   PetscQuadrature    quad;
96   PetscReal          Eq_m[LANDAU_MAX_SPECIES]; // could be static data w/o quench (ex2)
97   PetscScalar       *cellClosure = NULL;
98   const PetscScalar *xdata       = NULL;
99   PetscDS            prob;
100   PetscContainer     container;
101   P4estVertexMaps   *maps;
102   Mat                subJ[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ];
103 
104   PetscFunctionBegin;
105   PetscValidHeaderSpecific(a_X, VEC_CLASSID, 1);
106   PetscValidHeaderSpecific(JacP, MAT_CLASSID, 2);
107   PetscValidPointer(ctx, 5);
108   /* check for matrix container for GPU assembly. Support CPU assembly for debugging */
109   PetscCheck(ctx->plex[0] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
110   PetscCall(PetscLogEventBegin(ctx->events[10], 0, 0, 0, 0));
111   PetscCall(DMGetDS(ctx->plex[0], &prob)); // same DS for all grids
112   PetscCall(PetscObjectQuery((PetscObject)JacP, "assembly_maps", (PetscObject *)&container));
113   if (container) {
114     PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "maps but no GPU assembly");
115     PetscCall(PetscContainerGetPointer(container, (void **)&maps));
116     PetscCheck(maps, ctx->comm, PETSC_ERR_ARG_WRONG, "empty GPU matrix container");
117     for (PetscInt i = 0; i < ctx->num_grids * ctx->batch_sz; i++) subJ[i] = NULL;
118   } else {
119     PetscCheck(!ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "No maps but GPU assembly");
120     for (PetscInt tid = 0; tid < ctx->batch_sz; tid++) {
121       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCreateMatrix(ctx->plex[grid], &subJ[LAND_PACK_IDX(tid, grid)]));
122     }
123     maps = NULL;
124   }
125   // get dynamic data (Eq is odd, for quench and Spitzer test) for CPU assembly and raw data for Jacobian GPU assembly. Get host numCells[], Nq (yuck)
126   PetscCall(PetscFEGetQuadrature(ctx->fe[0], &quad));
127   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
128   PetscCall(PetscFEGetDimension(ctx->fe[0], &Nb));
129   PetscCheck(Nb == Nq, ctx->comm, PETSC_ERR_ARG_WRONG, "Nb (%d) != Nq (%d) - not true for some simplex elements (use P2)", (int)Nb, (int)Nq);
130   PetscCheck(Nq <= LANDAU_MAX_NQ, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQ (%d)", Nq, LANDAU_MAX_NQ);
131   // get metadata for collecting dynamic data
132   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
133     PetscInt cStart, cEnd;
134     PetscCheck(ctx->plex[grid] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
135     PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
136     numCells[grid] = cEnd - cStart; // grids can have different topology
137   }
138   PetscCall(PetscLogEventEnd(ctx->events[10], 0, 0, 0, 0));
139   if (shift == 0) { /* create dynamic point data: f_alpha for closure of each cell (cellClosure[nbatch,ngrids,ncells[g],f[Nb,ns[g]]]) or xdata */
140     DM pack;
141     PetscCall(VecGetDM(a_X, &pack));
142     PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pack has no DM");
143     PetscCall(PetscLogEventBegin(ctx->events[1], 0, 0, 0, 0));
144     for (PetscInt fieldA = 0; fieldA < ctx->num_species; fieldA++) {
145       Eq_m[fieldA] = ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */
146       if (dim == 2) Eq_m[fieldA] *= 2 * PETSC_PI;                                                  /* add the 2pi term that is not in Landau */
147     }
148     if (!ctx->gpu_assembly) {
149       Vec         *locXArray, *globXArray;
150       PetscScalar *cellClosure_it;
151       PetscInt     cellClosure_sz = 0, nDMs, Nf[LANDAU_MAX_GRIDS];
152       PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS];
153       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
154         PetscCall(DMGetLocalSection(ctx->plex[grid], &section[grid]));
155         PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid]));
156         PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid]));
157       }
158       /* count cellClosure size */
159       PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
160       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) cellClosure_sz += Nb * Nf[grid] * numCells[grid];
161       PetscCall(PetscMalloc1(cellClosure_sz * ctx->batch_sz, &cellClosure));
162       cellClosure_it = cellClosure;
163       PetscCall(PetscMalloc(sizeof(*locXArray) * nDMs, &locXArray));
164       PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
165       PetscCall(DMCompositeGetLocalAccessArray(pack, a_X, nDMs, NULL, locXArray));
166       PetscCall(DMCompositeGetAccessArray(pack, a_X, nDMs, NULL, globXArray));
167       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // OpenMP (once)
168         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
169           Vec      locX = locXArray[LAND_PACK_IDX(b_id, grid)], globX = globXArray[LAND_PACK_IDX(b_id, grid)], locX2;
170           PetscInt cStart, cEnd, ei;
171           PetscCall(VecDuplicate(locX, &locX2));
172           PetscCall(DMGlobalToLocalBegin(ctx->plex[grid], globX, INSERT_VALUES, locX2));
173           PetscCall(DMGlobalToLocalEnd(ctx->plex[grid], globX, INSERT_VALUES, locX2));
174           PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
175           for (ei = cStart; ei < cEnd; ++ei) {
176             PetscScalar *coef = NULL;
177             PetscCall(DMPlexVecGetClosure(ctx->plex[grid], section[grid], locX2, ei, NULL, &coef));
178             PetscCall(PetscMemcpy(cellClosure_it, coef, Nb * Nf[grid] * sizeof(*cellClosure_it))); /* change if LandauIPReal != PetscScalar */
179             PetscCall(DMPlexVecRestoreClosure(ctx->plex[grid], section[grid], locX2, ei, NULL, &coef));
180             cellClosure_it += Nb * Nf[grid];
181           }
182           PetscCall(VecDestroy(&locX2));
183         }
184       }
185       PetscCheck(cellClosure_it - cellClosure == cellClosure_sz * ctx->batch_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "iteration wrong %" PetscCount_FMT " != cellClosure_sz = %" PetscInt_FMT, (PetscCount)(cellClosure_it - cellClosure),
186                  cellClosure_sz * ctx->batch_sz);
187       PetscCall(DMCompositeRestoreLocalAccessArray(pack, a_X, nDMs, NULL, locXArray));
188       PetscCall(DMCompositeRestoreAccessArray(pack, a_X, nDMs, NULL, globXArray));
189       PetscCall(PetscFree(locXArray));
190       PetscCall(PetscFree(globXArray));
191       xdata = NULL;
192     } else {
193       PetscMemType mtype;
194       if (ctx->jacobian_field_major_order) { // get data in batch ordering
195         PetscCall(VecScatterBegin(ctx->plex_batch, a_X, ctx->work_vec, INSERT_VALUES, SCATTER_FORWARD));
196         PetscCall(VecScatterEnd(ctx->plex_batch, a_X, ctx->work_vec, INSERT_VALUES, SCATTER_FORWARD));
197         PetscCall(VecGetArrayReadAndMemType(ctx->work_vec, &xdata, &mtype));
198       } else {
199         PetscCall(VecGetArrayReadAndMemType(a_X, &xdata, &mtype));
200       }
201       PetscCheck(mtype == PETSC_MEMTYPE_HOST || ctx->deviceType != LANDAU_CPU, ctx->comm, PETSC_ERR_ARG_WRONG, "CPU run with device data: use -mat_type aij");
202       cellClosure = NULL;
203     }
204     PetscCall(PetscLogEventEnd(ctx->events[1], 0, 0, 0, 0));
205   } else xdata = cellClosure = NULL;
206 
207   /* do it */
208   if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) {
209     if (ctx->deviceType == LANDAU_CUDA) {
210 #if defined(PETSC_HAVE_CUDA)
211       PetscCall(LandauCUDAJacobian(ctx->plex, Nq, ctx->batch_sz, ctx->num_grids, numCells, Eq_m, cellClosure, xdata, &ctx->SData_d, shift, ctx->events, ctx->mat_offset, ctx->species_offset, subJ, JacP));
212 #else
213       SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "cuda");
214 #endif
215     } else if (ctx->deviceType == LANDAU_KOKKOS) {
216 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
217       PetscCall(LandauKokkosJacobian(ctx->plex, Nq, ctx->batch_sz, ctx->num_grids, numCells, Eq_m, cellClosure, xdata, &ctx->SData_d, shift, ctx->events, ctx->mat_offset, ctx->species_offset, subJ, JacP));
218 #else
219       SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "kokkos");
220 #endif
221     }
222   } else {               /* CPU version */
223     PetscTabulation *Tf; // used for CPU and print info. Same on all grids and all species
224     PetscInt         ip_offset[LANDAU_MAX_GRIDS + 1], ipf_offset[LANDAU_MAX_GRIDS + 1], elem_offset[LANDAU_MAX_GRIDS + 1], IPf_sz_glb, IPf_sz_tot, num_grids = ctx->num_grids, Nf[LANDAU_MAX_GRIDS];
225     PetscReal       *ff, *dudx, *dudy, *dudz, *invJ_a = (PetscReal *)ctx->SData_d.invJ, *xx = (PetscReal *)ctx->SData_d.x, *yy = (PetscReal *)ctx->SData_d.y, *zz = (PetscReal *)ctx->SData_d.z, *ww = (PetscReal *)ctx->SData_d.w;
226     PetscReal       *nu_alpha = (PetscReal *)ctx->SData_d.alpha, *nu_beta = (PetscReal *)ctx->SData_d.beta, *invMass = (PetscReal *)ctx->SData_d.invMass;
227     PetscReal(*lambdas)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS] = (PetscReal(*)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS])ctx->SData_d.lambdas;
228     PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS];
229     PetscScalar *coo_vals = NULL;
230     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
231       PetscCall(DMGetLocalSection(ctx->plex[grid], &section[grid]));
232       PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid]));
233       PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid]));
234     }
235     /* count IPf size, etc */
236     PetscCall(PetscDSGetTabulation(prob, &Tf)); // Bf, &Df same for all grids
237     const PetscReal *const BB = Tf[0]->T[0], *const DD = Tf[0]->T[1];
238     ip_offset[0] = ipf_offset[0] = elem_offset[0] = 0;
239     for (PetscInt grid = 0; grid < num_grids; grid++) {
240       PetscInt nfloc        = ctx->species_offset[grid + 1] - ctx->species_offset[grid];
241       elem_offset[grid + 1] = elem_offset[grid] + numCells[grid];
242       ip_offset[grid + 1]   = ip_offset[grid] + numCells[grid] * Nq;
243       ipf_offset[grid + 1]  = ipf_offset[grid] + Nq * nfloc * numCells[grid];
244     }
245     IPf_sz_glb = ipf_offset[num_grids];
246     IPf_sz_tot = IPf_sz_glb * ctx->batch_sz;
247     // prep COO
248     if (ctx->coo_assembly) {
249       PetscCall(PetscMalloc1(ctx->SData_d.coo_size, &coo_vals)); // allocate every time?
250       PetscCall(PetscInfo(ctx->plex[0], "COO Allocate %" PetscInt_FMT " values\n", (PetscInt)ctx->SData_d.coo_size));
251     }
252     if (shift == 0.0) { /* compute dynamic data f and df and init data for Jacobian */
253 #if defined(PETSC_HAVE_THREADSAFETY)
254       double starttime, endtime;
255       starttime = MPI_Wtime();
256 #endif
257       PetscCall(PetscLogEventBegin(ctx->events[8], 0, 0, 0, 0));
258       PetscCall(PetscMalloc4(IPf_sz_tot, &ff, IPf_sz_tot, &dudx, IPf_sz_tot, &dudy, dim == 3 ? IPf_sz_tot : 0, &dudz));
259       // F df/dx
260       for (PetscInt tid = 0; tid < ctx->batch_sz * elem_offset[num_grids]; tid++) {                        // for each element
261         const PetscInt b_Nelem = elem_offset[num_grids], b_elem_idx = tid % b_Nelem, b_id = tid / b_Nelem; // b_id == OMP thd_id in batch
262         // find my grid:
263         PetscInt grid = 0;
264         while (b_elem_idx >= elem_offset[grid + 1]) grid++; // yuck search for grid
265         {
266           const PetscInt loc_nip = numCells[grid] * Nq, loc_Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], loc_elem = b_elem_idx - elem_offset[grid];
267           const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); //b_id*b_N + ctx->mat_offset[grid];
268           PetscScalar   *coef, coef_buff[LANDAU_MAX_SPECIES * LANDAU_MAX_NQ];
269           PetscReal     *invJe = &invJ_a[(ip_offset[grid] + loc_elem * Nq) * dim * dim]; // ingJ is static data on batch 0
270           PetscInt       b, f, q;
271           if (cellClosure) {
272             coef = &cellClosure[b_id * IPf_sz_glb + ipf_offset[grid] + loc_elem * Nb * loc_Nf]; // this is const
273           } else {
274             coef = coef_buff;
275             for (f = 0; f < loc_Nf; ++f) {
276               LandauIdx *const Idxs = &maps[grid].gIdx[loc_elem][f][0];
277               for (b = 0; b < Nb; ++b) {
278                 PetscInt idx = Idxs[b];
279                 if (idx >= 0) {
280                   coef[f * Nb + b] = xdata[idx + moffset];
281                 } else {
282                   idx              = -idx - 1;
283                   coef[f * Nb + b] = 0;
284                   for (q = 0; q < maps[grid].num_face; q++) {
285                     PetscInt    id    = maps[grid].c_maps[idx][q].gid;
286                     PetscScalar scale = maps[grid].c_maps[idx][q].scale;
287                     coef[f * Nb + b] += scale * xdata[id + moffset];
288                   }
289                 }
290               }
291             }
292           }
293           /* get f and df */
294           for (PetscInt qi = 0; qi < Nq; qi++) {
295             const PetscReal *invJ = &invJe[qi * dim * dim];
296             const PetscReal *Bq   = &BB[qi * Nb];
297             const PetscReal *Dq   = &DD[qi * Nb * dim];
298             PetscReal        u_x[LANDAU_DIM];
299             /* get f & df */
300             for (f = 0; f < loc_Nf; ++f) {
301               const PetscInt idx = b_id * IPf_sz_glb + ipf_offset[grid] + f * loc_nip + loc_elem * Nq + qi;
302               PetscInt       b, e;
303               PetscReal      refSpaceDer[LANDAU_DIM];
304               ff[idx] = 0.0;
305               for (int d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0;
306               for (b = 0; b < Nb; ++b) {
307                 const PetscInt cidx = b;
308                 ff[idx] += Bq[cidx] * PetscRealPart(coef[f * Nb + cidx]);
309                 for (int d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx * dim + d] * PetscRealPart(coef[f * Nb + cidx]);
310               }
311               for (int d = 0; d < LANDAU_DIM; ++d) {
312                 for (e = 0, u_x[d] = 0.0; e < LANDAU_DIM; ++e) u_x[d] += invJ[e * dim + d] * refSpaceDer[e];
313               }
314               dudx[idx] = u_x[0];
315               dudy[idx] = u_x[1];
316 #if LANDAU_DIM == 3
317               dudz[idx] = u_x[2];
318 #endif
319             }
320           } // q
321         }   // grid
322       }     // grid*batch
323       PetscCall(PetscLogEventEnd(ctx->events[8], 0, 0, 0, 0));
324 #if defined(PETSC_HAVE_THREADSAFETY)
325       endtime = MPI_Wtime();
326       if (ctx->stage) ctx->times[LANDAU_F_DF] += (endtime - starttime);
327 #endif
328     } // Jacobian setup
329     // assemble Jacobian (or mass)
330     for (PetscInt tid = 0; tid < ctx->batch_sz * elem_offset[num_grids]; tid++) { // for each element
331       const PetscInt b_Nelem      = elem_offset[num_grids];
332       const PetscInt glb_elem_idx = tid % b_Nelem, b_id = tid / b_Nelem;
333       PetscInt       grid = 0;
334 #if defined(PETSC_HAVE_THREADSAFETY)
335       double starttime, endtime;
336       starttime = MPI_Wtime();
337 #endif
338       while (glb_elem_idx >= elem_offset[grid + 1]) grid++;
339       {
340         const PetscInt   loc_Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], loc_elem = glb_elem_idx - elem_offset[grid];
341         const PetscInt   moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset), totDim = loc_Nf * Nq, elemMatSize = totDim * totDim;
342         PetscScalar     *elemMat;
343         const PetscReal *invJe = &invJ_a[(ip_offset[grid] + loc_elem * Nq) * dim * dim];
344         PetscCall(PetscMalloc1(elemMatSize, &elemMat));
345         PetscCall(PetscMemzero(elemMat, elemMatSize * sizeof(*elemMat)));
346         if (shift == 0.0) { // Jacobian
347           PetscCall(PetscLogEventBegin(ctx->events[4], 0, 0, 0, 0));
348         } else { // mass
349           PetscCall(PetscLogEventBegin(ctx->events[16], 0, 0, 0, 0));
350         }
351         for (PetscInt qj = 0; qj < Nq; ++qj) {
352           const PetscInt jpidx_glb = ip_offset[grid] + qj + loc_elem * Nq;
353           PetscReal      g0[LANDAU_MAX_SPECIES], g2[LANDAU_MAX_SPECIES][LANDAU_DIM], g3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM]; // could make a LANDAU_MAX_SPECIES_GRID ~ number of ions - 1
354           PetscInt       d, d2, dp, d3, IPf_idx;
355           if (shift == 0.0) { // Jacobian
356             const PetscReal *const invJj = &invJe[qj * dim * dim];
357             PetscReal              gg2[LANDAU_MAX_SPECIES][LANDAU_DIM], gg3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM], gg2_temp[LANDAU_DIM], gg3_temp[LANDAU_DIM][LANDAU_DIM];
358             const PetscReal        vj[3] = {xx[jpidx_glb], yy[jpidx_glb], zz ? zz[jpidx_glb] : 0}, wj = ww[jpidx_glb];
359             // create g2 & g3
360             for (d = 0; d < LANDAU_DIM; d++) { // clear accumulation data D & K
361               gg2_temp[d] = 0;
362               for (d2 = 0; d2 < LANDAU_DIM; d2++) gg3_temp[d][d2] = 0;
363             }
364             /* inner beta reduction */
365             IPf_idx = 0;
366             for (PetscInt grid_r = 0, f_off = 0, ipidx = 0; grid_r < ctx->num_grids; grid_r++, f_off = ctx->species_offset[grid_r]) { // IPf_idx += nip_loc_r*Nfloc_r
367               PetscInt nip_loc_r = numCells[grid_r] * Nq, Nfloc_r = Nf[grid_r];
368               for (PetscInt ei_r = 0, loc_fdf_idx = 0; ei_r < numCells[grid_r]; ++ei_r) {
369                 for (PetscInt qi = 0; qi < Nq; qi++, ipidx++, loc_fdf_idx++) {
370                   const PetscReal wi = ww[ipidx], x = xx[ipidx], y = yy[ipidx];
371                   PetscReal       temp1[3] = {0, 0, 0}, temp2 = 0;
372 #if LANDAU_DIM == 2
373                   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.;
374                   LandauTensor2D(vj, x, y, Ud, Uk, mask);
375 #else
376                   PetscReal U[3][3], z = zz[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.;
377                   if (ctx->use_relativistic_corrections) {
378                     LandauTensor3DRelativistic(vj, x, y, z, U, mask, C_0(ctx->v_0));
379                   } else {
380                     LandauTensor3D(vj, x, y, z, U, mask);
381                   }
382 #endif
383                   for (int f = 0; f < Nfloc_r; ++f) {
384                     const PetscInt idx = b_id * IPf_sz_glb + ipf_offset[grid_r] + f * nip_loc_r + ei_r * Nq + qi; // IPf_idx + f*nip_loc_r + loc_fdf_idx;
385                     temp1[0] += dudx[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r];
386                     temp1[1] += dudy[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r];
387 #if LANDAU_DIM == 3
388                     temp1[2] += dudz[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r];
389 #endif
390                     temp2 += ff[idx] * nu_beta[f + f_off] * (*lambdas)[grid][grid_r];
391                   }
392                   temp1[0] *= wi;
393                   temp1[1] *= wi;
394 #if LANDAU_DIM == 3
395                   temp1[2] *= wi;
396 #endif
397                   temp2 *= wi;
398 #if LANDAU_DIM == 2
399                   for (d2 = 0; d2 < 2; d2++) {
400                     for (d3 = 0; d3 < 2; ++d3) {
401                       /* K = U * grad(f): g2=e: i,A */
402                       gg2_temp[d2] += Uk[d2][d3] * temp1[d3];
403                       /* D = -U * (I \kron (fx)): g3=f: i,j,A */
404                       gg3_temp[d2][d3] += Ud[d2][d3] * temp2;
405                     }
406                   }
407 #else
408                   for (d2 = 0; d2 < 3; ++d2) {
409                     for (d3 = 0; d3 < 3; ++d3) {
410                       /* K = U * grad(f): g2 = e: i,A */
411                       gg2_temp[d2] += U[d2][d3] * temp1[d3];
412                       /* D = -U * (I \kron (fx)): g3 = f: i,j,A */
413                       gg3_temp[d2][d3] += U[d2][d3] * temp2;
414                     }
415                   }
416 #endif
417                 } // qi
418               }   // ei_r
419               IPf_idx += nip_loc_r * Nfloc_r;
420             } /* grid_r - IPs */
421             PetscCheck(IPf_idx == IPf_sz_glb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IPf_idx != IPf_sz %" PetscInt_FMT " %" PetscInt_FMT, IPf_idx, IPf_sz_glb);
422             // add alpha and put in gg2/3
423             for (PetscInt fieldA = 0, f_off = ctx->species_offset[grid]; fieldA < loc_Nf; ++fieldA) {
424               for (d2 = 0; d2 < LANDAU_DIM; d2++) {
425                 gg2[fieldA][d2] = gg2_temp[d2] * nu_alpha[fieldA + f_off];
426                 for (d3 = 0; d3 < LANDAU_DIM; d3++) gg3[fieldA][d2][d3] = -gg3_temp[d2][d3] * nu_alpha[fieldA + f_off] * invMass[fieldA + f_off];
427               }
428             }
429             /* add electric field term once per IP */
430             for (PetscInt fieldA = 0, f_off = ctx->species_offset[grid]; fieldA < loc_Nf; ++fieldA) gg2[fieldA][LANDAU_DIM - 1] += Eq_m[fieldA + f_off];
431             /* Jacobian transform - g2, g3 */
432             for (PetscInt fieldA = 0; fieldA < loc_Nf; ++fieldA) {
433               for (d = 0; d < dim; ++d) {
434                 g2[fieldA][d] = 0.0;
435                 for (d2 = 0; d2 < dim; ++d2) {
436                   g2[fieldA][d] += invJj[d * dim + d2] * gg2[fieldA][d2];
437                   g3[fieldA][d][d2] = 0.0;
438                   for (d3 = 0; d3 < dim; ++d3) {
439                     for (dp = 0; dp < dim; ++dp) g3[fieldA][d][d2] += invJj[d * dim + d3] * gg3[fieldA][d3][dp] * invJj[d2 * dim + dp];
440                   }
441                   g3[fieldA][d][d2] *= wj;
442                 }
443                 g2[fieldA][d] *= wj;
444               }
445             }
446           } else { // mass
447             PetscReal wj = ww[jpidx_glb];
448             /* Jacobian transform - g0 */
449             for (PetscInt fieldA = 0; fieldA < loc_Nf; ++fieldA) {
450               if (dim == 2) {
451                 g0[fieldA] = wj * shift * 2. * PETSC_PI; // move this to below and remove g0
452               } else {
453                 g0[fieldA] = wj * shift; // move this to below and remove g0
454               }
455             }
456           }
457           /* FE matrix construction */
458           {
459             PetscInt         fieldA, d, f, d2, g;
460             const PetscReal *BJq = &BB[qj * Nb], *DIq = &DD[qj * Nb * dim];
461             /* assemble - on the diagonal (I,I) */
462             for (fieldA = 0; fieldA < loc_Nf; fieldA++) {
463               for (f = 0; f < Nb; f++) {
464                 const PetscInt i = fieldA * Nb + f; /* Element matrix row */
465                 for (g = 0; g < Nb; ++g) {
466                   const PetscInt j    = fieldA * Nb + g; /* Element matrix column */
467                   const PetscInt fOff = i * totDim + j;
468                   if (shift == 0.0) {
469                     for (d = 0; d < dim; ++d) {
470                       elemMat[fOff] += DIq[f * dim + d] * g2[fieldA][d] * BJq[g];
471                       for (d2 = 0; d2 < dim; ++d2) elemMat[fOff] += DIq[f * dim + d] * g3[fieldA][d][d2] * DIq[g * dim + d2];
472                     }
473                   } else { // mass
474                     elemMat[fOff] += BJq[f] * g0[fieldA] * BJq[g];
475                   }
476                 }
477               }
478             }
479           }
480         }                   /* qj loop */
481         if (shift == 0.0) { // Jacobian
482           PetscCall(PetscLogEventEnd(ctx->events[4], 0, 0, 0, 0));
483         } else {
484           PetscCall(PetscLogEventEnd(ctx->events[16], 0, 0, 0, 0));
485         }
486 #if defined(PETSC_HAVE_THREADSAFETY)
487         endtime = MPI_Wtime();
488         if (ctx->stage) ctx->times[LANDAU_KERNEL] += (endtime - starttime);
489 #endif
490         /* assemble matrix */
491         if (!container) {
492           PetscInt cStart;
493           PetscCall(PetscLogEventBegin(ctx->events[6], 0, 0, 0, 0));
494           PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, NULL));
495           PetscCall(DMPlexMatSetClosure(ctx->plex[grid], section[grid], globsection[grid], subJ[LAND_PACK_IDX(b_id, grid)], loc_elem + cStart, elemMat, ADD_VALUES));
496           PetscCall(PetscLogEventEnd(ctx->events[6], 0, 0, 0, 0));
497         } else { // GPU like assembly for debugging
498           PetscInt    fieldA, q, f, g, d, nr, nc, rows0[LANDAU_MAX_Q_FACE] = {0}, cols0[LANDAU_MAX_Q_FACE] = {0}, rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE];
499           PetscScalar vals[LANDAU_MAX_Q_FACE * LANDAU_MAX_Q_FACE] = {0}, row_scale[LANDAU_MAX_Q_FACE] = {0}, col_scale[LANDAU_MAX_Q_FACE] = {0};
500           LandauIdx *coo_elem_offsets = (LandauIdx *)ctx->SData_d.coo_elem_offsets, *coo_elem_fullNb = (LandauIdx *)ctx->SData_d.coo_elem_fullNb, (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = (LandauIdx(*)[LANDAU_MAX_NQ + 1]) ctx->SData_d.coo_elem_point_offsets;
501           /* assemble - from the diagonal (I,I) in this format for DMPlexMatSetClosure */
502           for (fieldA = 0; fieldA < loc_Nf; fieldA++) {
503             LandauIdx *const Idxs = &maps[grid].gIdx[loc_elem][fieldA][0];
504             for (f = 0; f < Nb; f++) {
505               PetscInt idx = Idxs[f];
506               if (idx >= 0) {
507                 nr           = 1;
508                 rows0[0]     = idx;
509                 row_scale[0] = 1.;
510               } else {
511                 idx = -idx - 1;
512                 for (q = 0, nr = 0; q < maps[grid].num_face; q++, nr++) {
513                   if (maps[grid].c_maps[idx][q].gid < 0) break;
514                   rows0[q]     = maps[grid].c_maps[idx][q].gid;
515                   row_scale[q] = maps[grid].c_maps[idx][q].scale;
516                 }
517               }
518               for (g = 0; g < Nb; ++g) {
519                 idx = Idxs[g];
520                 if (idx >= 0) {
521                   nc           = 1;
522                   cols0[0]     = idx;
523                   col_scale[0] = 1.;
524                 } else {
525                   idx = -idx - 1;
526                   nc  = maps[grid].num_face;
527                   for (q = 0, nc = 0; q < maps[grid].num_face; q++, nc++) {
528                     if (maps[grid].c_maps[idx][q].gid < 0) break;
529                     cols0[q]     = maps[grid].c_maps[idx][q].gid;
530                     col_scale[q] = maps[grid].c_maps[idx][q].scale;
531                   }
532                 }
533                 const PetscInt    i   = fieldA * Nb + f; /* Element matrix row */
534                 const PetscInt    j   = fieldA * Nb + g; /* Element matrix column */
535                 const PetscScalar Aij = elemMat[i * totDim + j];
536                 if (coo_vals) { // mirror (i,j) in CreateStaticGPUData
537                   const int fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb;
538                   const int idx0 = b_id * coo_elem_offsets[elem_offset[num_grids]] + coo_elem_offsets[glb_elem_idx] + fieldA * fullNb2 + fullNb * coo_elem_point_offsets[glb_elem_idx][f] + nr * coo_elem_point_offsets[glb_elem_idx][g];
539                   for (int q = 0, idx2 = idx0; q < nr; q++) {
540                     for (int d = 0; d < nc; d++, idx2++) coo_vals[idx2] = row_scale[q] * col_scale[d] * Aij;
541                   }
542                 } else {
543                   for (q = 0; q < nr; q++) rows[q] = rows0[q] + moffset;
544                   for (d = 0; d < nc; d++) cols[d] = cols0[d] + moffset;
545                   for (q = 0; q < nr; q++) {
546                     for (d = 0; d < nc; d++) vals[q * nc + d] = row_scale[q] * col_scale[d] * Aij;
547                   }
548                   PetscCall(MatSetValues(JacP, nr, rows, nc, cols, vals, ADD_VALUES));
549                 }
550               }
551             }
552           }
553         }
554         if (loc_elem == -1) {
555           PetscCall(PetscPrintf(ctx->comm, "CPU Element matrix\n"));
556           for (int d = 0; d < totDim; ++d) {
557             for (int f = 0; f < totDim; ++f) PetscCall(PetscPrintf(ctx->comm, " %12.5e", (double)PetscRealPart(elemMat[d * totDim + f])));
558             PetscCall(PetscPrintf(ctx->comm, "\n"));
559           }
560           exit(12);
561         }
562         PetscCall(PetscFree(elemMat));
563       }                 /* grid */
564     }                   /* outer element & batch loop */
565     if (shift == 0.0) { // mass
566       PetscCall(PetscFree4(ff, dudx, dudy, dudz));
567     }
568     if (!container) {                                         // 'CPU' assembly move nest matrix to global JacP
569       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // OpenMP
570         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
571           const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); // b_id*b_N + ctx->mat_offset[grid];
572           PetscInt           nloc, nzl, colbuf[1024], row;
573           const PetscInt    *cols;
574           const PetscScalar *vals;
575           Mat                B = subJ[LAND_PACK_IDX(b_id, grid)];
576           PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
577           PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
578           PetscCall(MatGetSize(B, &nloc, NULL));
579           for (int i = 0; i < nloc; i++) {
580             PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
581             PetscCheck(nzl <= 1024, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "Row too big: %" PetscInt_FMT, nzl);
582             for (int j = 0; j < nzl; j++) colbuf[j] = moffset + cols[j];
583             row = moffset + i;
584             PetscCall(MatSetValues(JacP, 1, &row, nzl, colbuf, vals, ADD_VALUES));
585             PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
586           }
587           PetscCall(MatDestroy(&B));
588         }
589       }
590     }
591     if (coo_vals) {
592       PetscCall(MatSetValuesCOO(JacP, coo_vals, ADD_VALUES));
593       PetscCall(PetscFree(coo_vals));
594     }
595   } /* CPU version */
596   PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
597   PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
598   /* clean up */
599   if (cellClosure) PetscCall(PetscFree(cellClosure));
600   if (xdata) PetscCall(VecRestoreArrayReadAndMemType(a_X, &xdata));
601   PetscFunctionReturn(PETSC_SUCCESS);
602 }
603 
604 static PetscErrorCode GeometryDMLandau(DM base, PetscInt point, PetscInt dim, const PetscReal abc[], PetscReal xyz[], void *a_ctx)
605 {
606   PetscReal r = abc[0], z = abc[1];
607 
608   PetscFunctionBegin;
609   xyz[0] = r;
610   xyz[1] = z;
611   if (dim == 3) xyz[2] = abc[2];
612 
613   PetscFunctionReturn(PETSC_SUCCESS);
614 }
615 
616 /* create DMComposite of meshes for each species group */
617 static PetscErrorCode LandauDMCreateVMeshes(MPI_Comm comm_self, const PetscInt dim, const char prefix[], LandauCtx *ctx, DM pack)
618 {
619   PetscFunctionBegin;
620   { /* p4est, quads */
621     /* Create plex mesh of Landau domain */
622     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
623       PetscReal par_radius = ctx->radius_par[grid], perp_radius = ctx->radius_perp[grid];
624       if (!ctx->sphere && !ctx->simplex) { // 2 or 3D (only 3D option)
625         PetscReal      lo[] = {-perp_radius, -par_radius, -par_radius}, hi[] = {perp_radius, par_radius, par_radius};
626         DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, dim == 2 ? DM_BOUNDARY_NONE : DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
627         if (dim == 2) lo[0] = 0;
628         else {
629           lo[1] = -perp_radius;
630           hi[1] = perp_radius; // 3D y is a perp
631         }
632         PetscCall(DMPlexCreateBoxMesh(comm_self, dim, PETSC_FALSE, ctx->cells0, lo, hi, periodicity, PETSC_TRUE, &ctx->plex[grid])); // todo: make composite and create dm[grid] here
633         PetscCall(DMLocalizeCoordinates(ctx->plex[grid]));                                                                           /* needed for periodic */
634         if (dim == 3) PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "cube"));
635         else PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "half-plane"));
636       } else if (dim == 2) {
637         char      filename[PETSC_MAX_PATH_LEN];
638         PetscBool flg;
639         PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_landau_filename", filename, sizeof(filename), &flg));
640         if (flg) {
641           char str[] = "-dm_landau_view_file_0";
642           str[21] += grid;
643           PetscCall(DMPlexCreateFromFile(comm_self, filename, "plexland.c", PETSC_TRUE, &ctx->plex[grid]));
644           PetscCall(PetscInfo(ctx->plex[grid], "%d) Read %s mesh file (%s)", (int)grid, filename, str));
645           PetscCall(DMViewFromOptions(ctx->plex[grid], NULL, str));
646         } else {
647           PetscInt       numCells = ctx->simplex ? 12 : 6, cell_size = ctx->simplex ? 3 : 4, j;
648           const PetscInt numVerts    = 11;
649           PetscInt       cellsT[][4] = {
650             {0,  1, 6, 5 },
651             {1,  2, 7, 6 },
652             {2,  3, 8, 7 },
653             {3,  4, 9, 8 },
654             {5,  6, 7, 10},
655             {10, 7, 8, 9 }
656           };
657           PetscInt cellsS[][3] = {
658             {0,  1, 6 },
659             {1,  2, 6 },
660             {6,  2, 7 },
661             {7,  2, 8 },
662             {8,  2, 3 },
663             {8,  3, 4 },
664             {0,  6, 5 },
665             {5,  6, 7 },
666             {5,  7, 10},
667             {10, 7, 9 },
668             {9,  7, 8 },
669             {9,  8, 4 }
670           };
671           const PetscInt *pcell = (const PetscInt *)(ctx->simplex ? &cellsS[0][0] : &cellsT[0][0]);
672           PetscReal       coords[11][2], *flatCoords = (PetscReal *)&coords[0][0];
673           PetscReal       rad = ctx->radius[grid];
674           for (j = 0; j < 5; j++) { // outside edge
675             PetscReal z, r, theta = -PETSC_PI / 2 + (j % 5) * PETSC_PI / 4;
676             r            = rad * PetscCosReal(theta);
677             coords[j][0] = r;
678             z            = rad * PetscSinReal(theta);
679             coords[j][1] = z;
680           }
681           coords[j][0]   = 0;
682           coords[j++][1] = -rad * ctx->sphere_inner_radius_90degree;
683           coords[j][0]   = rad * ctx->sphere_inner_radius_45degree;
684           coords[j++][1] = -rad * ctx->sphere_inner_radius_45degree;
685           coords[j][0]   = rad * ctx->sphere_inner_radius_90degree;
686           coords[j++][1] = 0;
687           coords[j][0]   = rad * ctx->sphere_inner_radius_45degree;
688           coords[j++][1] = rad * ctx->sphere_inner_radius_45degree;
689           coords[j][0]   = 0;
690           coords[j++][1] = rad * ctx->sphere_inner_radius_90degree;
691           coords[j][0]   = 0;
692           coords[j++][1] = 0;
693           PetscCall(DMPlexCreateFromCellListPetsc(comm_self, 2, numCells, numVerts, cell_size, ctx->interpolate, pcell, 2, flatCoords, &ctx->plex[grid]));
694           PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "semi-circle"));
695           PetscCall(PetscInfo(ctx->plex[grid], "\t%" PetscInt_FMT ") Make circle %s mesh", grid, ctx->simplex ? "simplex" : "tensor"));
696         }
697       } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Velocity space meshes does not support 3V cubed sphere or simplex");
698       PetscCall(DMSetFromOptions(ctx->plex[grid]));
699     } // grid loop
700     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pack, prefix));
701     { /* convert to p4est (or whatever), wait for discretization to create pack */
702       char      convType[256];
703       PetscBool flg;
704 
705       PetscOptionsBegin(ctx->comm, prefix, "Mesh conversion options", "DMPLEX");
706       PetscCall(PetscOptionsFList("-dm_landau_type", "Convert DMPlex to another format (p4est)", "plexland.c", DMList, DMPLEX, convType, 256, &flg));
707       PetscOptionsEnd();
708       if (flg) {
709         ctx->use_p4est = PETSC_TRUE; /* flag for Forest */
710         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
711           DM dmforest;
712           PetscCall(DMConvert(ctx->plex[grid], convType, &dmforest));
713           if (dmforest) {
714             PetscBool isForest;
715             PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dmforest, prefix));
716             PetscCall(DMIsForest(dmforest, &isForest));
717             if (isForest) {
718               if (ctx->sphere) PetscCall(DMForestSetBaseCoordinateMapping(dmforest, GeometryDMLandau, ctx));
719               PetscCall(DMDestroy(&ctx->plex[grid]));
720               ctx->plex[grid] = dmforest; // Forest for adaptivity
721             } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Converted to non Forest?");
722           } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Convert failed?");
723         }
724       } else ctx->use_p4est = PETSC_FALSE; /* flag for Forest */
725     }
726   } /* non-file */
727   PetscCall(DMSetDimension(pack, dim));
728   PetscCall(PetscObjectSetName((PetscObject)pack, "Mesh"));
729   PetscCall(DMSetApplicationContext(pack, ctx));
730 
731   PetscFunctionReturn(PETSC_SUCCESS);
732 }
733 
734 static PetscErrorCode SetupDS(DM pack, PetscInt dim, PetscInt grid, LandauCtx *ctx)
735 {
736   PetscInt     ii, i0;
737   char         buf[256];
738   PetscSection section;
739 
740   PetscFunctionBegin;
741   for (ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
742     if (ii == 0) PetscCall(PetscSNPrintf(buf, sizeof(buf), "e"));
743     else PetscCall(PetscSNPrintf(buf, sizeof(buf), "i%" PetscInt_FMT, ii));
744     /* Setup Discretization - FEM */
745     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, ctx->simplex, NULL, PETSC_DECIDE, &ctx->fe[ii]));
746     PetscCall(PetscObjectSetName((PetscObject)ctx->fe[ii], buf));
747     PetscCall(DMSetField(ctx->plex[grid], i0, NULL, (PetscObject)ctx->fe[ii]));
748   }
749   PetscCall(DMCreateDS(ctx->plex[grid]));
750   PetscCall(DMGetSection(ctx->plex[grid], &section));
751   for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
752     if (ii == 0) PetscCall(PetscSNPrintf(buf, sizeof(buf), "se"));
753     else PetscCall(PetscSNPrintf(buf, sizeof(buf), "si%" PetscInt_FMT, ii));
754     PetscCall(PetscSectionSetComponentName(section, i0, 0, buf));
755   }
756   PetscFunctionReturn(PETSC_SUCCESS);
757 }
758 
759 /* Define a Maxwellian function for testing out the operator. */
760 
761 /* Using cartesian velocity space coordinates, the particle */
762 /* density, [1/m^3], is defined according to */
763 
764 /* $$ n=\int_{R^3} dv^3 \left(\frac{m}{2\pi T}\right)^{3/2}\exp [- mv^2/(2T)] $$ */
765 
766 /* Using some constant, c, we normalize the velocity vector into a */
767 /* dimensionless variable according to v=c*x. Thus the density, $n$, becomes */
768 
769 /* $$ n=\int_{R^3} dx^3 \left(\frac{mc^2}{2\pi T}\right)^{3/2}\exp [- mc^2/(2T)*x^2] $$ */
770 
771 /* Defining $\theta=2T/mc^2$, we thus find that the probability density */
772 /* for finding the particle within the interval in a box dx^3 around x is */
773 
774 /* f(x;\theta)=\left(\frac{1}{\pi\theta}\right)^{3/2} \exp [ -x^2/\theta ] */
775 
776 typedef struct {
777   PetscReal v_0;
778   PetscReal kT_m;
779   PetscReal n;
780   PetscReal shift;
781 } MaxwellianCtx;
782 
783 static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
784 {
785   MaxwellianCtx *mctx = (MaxwellianCtx *)actx;
786   PetscInt       i;
787   PetscReal      v2 = 0, theta = 2 * mctx->kT_m / (mctx->v_0 * mctx->v_0), shift; /* theta = 2kT/mc^2 */
788   PetscFunctionBegin;
789   /* compute the exponents, v^2 */
790   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
791   /* evaluate the Maxwellian */
792   if (mctx->shift < 0) shift = -mctx->shift;
793   else {
794     u[0]  = mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
795     shift = mctx->shift;
796   }
797   if (shift != 0.) {
798     v2 = 0;
799     for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
800     v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
801     /* evaluate the shifted Maxwellian */
802     u[0] += mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
803   }
804   PetscFunctionReturn(PETSC_SUCCESS);
805 }
806 
807 /*@
808  DMPlexLandauAddMaxwellians - Add a Maxwellian distribution to a state
809 
810  Collective
811 
812  Input Parameters:
813  .   dm - The mesh (local)
814  +   time - Current time
815  -   temps - Temperatures of each species (global)
816  .   ns - Number density of each species (global)
817  -   grid - index into current grid - just used for offset into temp and ns
818  .   b_id - batch index
819  -   n_batch - number of batches
820  +   actx - Landau context
821 
822  Output Parameter:
823  .   X  - The state (local to this grid)
824 
825  Level: beginner
826 
827  .keywords: mesh
828 .seealso: `DMPlexLandauCreateVelocitySpace()`
829  @*/
830 PetscErrorCode DMPlexLandauAddMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], PetscInt grid, PetscInt b_id, PetscInt n_batch, void *actx)
831 {
832   LandauCtx *ctx = (LandauCtx *)actx;
833   PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
834   PetscInt       dim;
835   MaxwellianCtx *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES];
836 
837   PetscFunctionBegin;
838   PetscCall(DMGetDimension(dm, &dim));
839   if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx));
840   for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
841     mctxs[i0]      = &data[i0];
842     data[i0].v_0   = ctx->v_0;                             // v_0 same for all grids
843     data[i0].kT_m  = ctx->k * temps[ii] / ctx->masses[ii]; /* kT/m */
844     data[i0].n     = ns[ii];
845     initu[i0]      = maxwellian;
846     data[i0].shift = 0;
847   }
848   data[0].shift = ctx->electronShift;
849   /* need to make ADD_ALL_VALUES work - TODO */
850   PetscCall(DMProjectFunction(dm, time, initu, (void **)mctxs, INSERT_ALL_VALUES, X));
851   PetscFunctionReturn(PETSC_SUCCESS);
852 }
853 
854 /*
855  LandauSetInitialCondition - Adds Maxwellians with context
856 
857  Collective
858 
859  Input Parameters:
860  .   dm - The mesh
861  -   grid - index into current grid - just used for offset into temp and ns
862  .   b_id - batch index
863  -   n_batch - number of batches
864  +   actx - Landau context with T and n
865 
866  Output Parameter:
867  .   X  - The state
868 
869  Level: beginner
870 
871  .keywords: mesh
872 .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauAddMaxwellians()`
873  */
874 static PetscErrorCode LandauSetInitialCondition(DM dm, Vec X, PetscInt grid, PetscInt b_id, PetscInt n_batch, void *actx)
875 {
876   LandauCtx *ctx = (LandauCtx *)actx;
877   PetscFunctionBegin;
878   if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx));
879   PetscCall(VecZeroEntries(X));
880   PetscCall(DMPlexLandauAddMaxwellians(dm, X, 0.0, ctx->thermal_temps, ctx->n, grid, b_id, n_batch, ctx));
881   PetscFunctionReturn(PETSC_SUCCESS);
882 }
883 
884 // adapt a level once. Forest in/out
885 #if defined(PETSC_USE_INFO)
886 static const char *s_refine_names[] = {"RE", "Z1", "Origin", "Z2", "Uniform"};
887 #endif
888 static PetscErrorCode adaptToleranceFEM(PetscFE fem, Vec sol, PetscInt type, PetscInt grid, LandauCtx *ctx, DM *newForest)
889 {
890   DM              forest, plex, adaptedDM = NULL;
891   PetscDS         prob;
892   PetscBool       isForest;
893   PetscQuadrature quad;
894   PetscInt        Nq, *Nb, cStart, cEnd, c, dim, qj, k;
895   DMLabel         adaptLabel = NULL;
896 
897   PetscFunctionBegin;
898   forest = ctx->plex[grid];
899   PetscCall(DMCreateDS(forest));
900   PetscCall(DMGetDS(forest, &prob));
901   PetscCall(DMGetDimension(forest, &dim));
902   PetscCall(DMIsForest(forest, &isForest));
903   PetscCheck(isForest, ctx->comm, PETSC_ERR_ARG_WRONG, "! Forest");
904   PetscCall(DMConvert(forest, DMPLEX, &plex));
905   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
906   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
907   PetscCall(PetscFEGetQuadrature(fem, &quad));
908   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
909   PetscCheck(Nq <= LANDAU_MAX_NQ, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQ (%d)", Nq, LANDAU_MAX_NQ);
910   PetscCall(PetscDSGetDimensions(prob, &Nb));
911   PetscCall(PetscInfo(sol, "%" PetscInt_FMT ") Refine phase: %s\n", grid, s_refine_names[type]));
912   if (type == 4) {
913     for (c = cStart; c < cEnd; c++) PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));
914   } else if (type == 2) {
915     PetscInt  rCellIdx[8], nr = 0, nrmax = (dim == 3) ? 8 : 2;
916     PetscReal minRad = PETSC_INFINITY, r;
917     for (c = cStart; c < cEnd; c++) {
918       PetscReal tt, v0[LANDAU_MAX_NQ * 3], detJ[LANDAU_MAX_NQ];
919       PetscCall(DMPlexComputeCellGeometryFEM(plex, c, quad, v0, NULL, NULL, detJ));
920       for (qj = 0; qj < Nq; ++qj) {
921         tt = PetscSqr(v0[dim * qj + 0]) + PetscSqr(v0[dim * qj + 1]) + PetscSqr(((dim == 3) ? v0[dim * qj + 2] : 0));
922         r  = PetscSqrtReal(tt);
923         if (r < minRad - PETSC_SQRT_MACHINE_EPSILON * 10.) {
924           minRad         = r;
925           nr             = 0;
926           rCellIdx[nr++] = c;
927           PetscCall(PetscInfo(sol, "\t\t%" PetscInt_FMT ") Found first inner r=%e, cell %" PetscInt_FMT ", qp %" PetscInt_FMT "/%" PetscInt_FMT "\n", grid, (double)r, c, qj + 1, Nq));
928         } else if ((r - minRad) < PETSC_SQRT_MACHINE_EPSILON * 100. && nr < nrmax) {
929           for (k = 0; k < nr; k++)
930             if (c == rCellIdx[k]) break;
931           if (k == nr) {
932             rCellIdx[nr++] = c;
933             PetscCall(PetscInfo(sol, "\t\t\t%" PetscInt_FMT ") Found another inner r=%e, cell %" PetscInt_FMT ", qp %" PetscInt_FMT "/%" PetscInt_FMT ", d=%e\n", grid, (double)r, c, qj + 1, Nq, (double)(r - minRad)));
934           }
935         }
936       }
937     }
938     for (k = 0; k < nr; k++) PetscCall(DMLabelSetValue(adaptLabel, rCellIdx[k], DM_ADAPT_REFINE));
939     PetscCall(PetscInfo(sol, "\t\t\t%" PetscInt_FMT ") Refined %" PetscInt_FMT " origin cells %" PetscInt_FMT ",%" PetscInt_FMT " r=%g\n", grid, nr, rCellIdx[0], rCellIdx[1], (double)minRad));
940   } else if (type == 0 || type == 1 || type == 3) { /* refine along r=0 axis */
941     PetscScalar *coef = NULL;
942     Vec          coords;
943     PetscInt     csize, Nv, d, nz, nrefined = 0;
944     DM           cdm;
945     PetscSection cs;
946     PetscCall(DMGetCoordinatesLocal(forest, &coords));
947     PetscCall(DMGetCoordinateDM(forest, &cdm));
948     PetscCall(DMGetLocalSection(cdm, &cs));
949     for (c = cStart; c < cEnd; c++) {
950       PetscInt doit = 0, outside = 0;
951       PetscCall(DMPlexVecGetClosure(cdm, cs, coords, c, &csize, &coef));
952       Nv = csize / dim;
953       for (nz = d = 0; d < Nv; d++) {
954         PetscReal z = PetscRealPart(coef[d * dim + (dim - 1)]), x = PetscSqr(PetscRealPart(coef[d * dim + 0])) + ((dim == 3) ? PetscSqr(PetscRealPart(coef[d * dim + 1])) : 0);
955         x = PetscSqrtReal(x);
956         if (type == 0) {
957           if (ctx->re_radius > PETSC_SQRT_MACHINE_EPSILON && (z < -PETSC_MACHINE_EPSILON * 10. || z > ctx->re_radius + PETSC_MACHINE_EPSILON * 10.)) outside++; /* first pass don't refine bottom */
958         } else if (type == 1 && (z > ctx->vperp0_radius1 || z < -ctx->vperp0_radius1)) {
959           outside++; /* don't refine outside electron refine radius */
960           PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") (debug) found %s cells\n", grid, s_refine_names[type]));
961         } else if (type == 3 && (z > ctx->vperp0_radius2 || z < -ctx->vperp0_radius2)) {
962           outside++; /* refine r=0 cells on refinement front */
963           PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") (debug) found %s cells\n", grid, s_refine_names[type]));
964         }
965         if (x < PETSC_MACHINE_EPSILON * 10. && (type != 0 || ctx->re_radius > PETSC_SQRT_MACHINE_EPSILON)) nz++;
966       }
967       PetscCall(DMPlexVecRestoreClosure(cdm, cs, coords, c, &csize, &coef));
968       if (doit || (outside < Nv && nz)) {
969         PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));
970         nrefined++;
971       }
972     }
973     PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") Refined %" PetscInt_FMT " cells\n", grid, nrefined));
974   }
975   PetscCall(DMDestroy(&plex));
976   PetscCall(DMAdaptLabel(forest, adaptLabel, &adaptedDM));
977   PetscCall(DMLabelDestroy(&adaptLabel));
978   *newForest = adaptedDM;
979   if (adaptedDM) {
980     if (isForest) {
981       PetscCall(DMForestSetAdaptivityForest(adaptedDM, NULL)); // ????
982     }
983     PetscCall(DMConvert(adaptedDM, DMPLEX, &plex));
984     PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
985     PetscCall(PetscInfo(sol, "\t\t\t\t%" PetscInt_FMT ") %" PetscInt_FMT " cells, %" PetscInt_FMT " total quadrature points\n", grid, cEnd - cStart, Nq * (cEnd - cStart)));
986     PetscCall(DMDestroy(&plex));
987   } else *newForest = NULL;
988   PetscFunctionReturn(PETSC_SUCCESS);
989 }
990 
991 // forest goes in (ctx->plex[grid]), plex comes out
992 static PetscErrorCode adapt(PetscInt grid, LandauCtx *ctx, Vec *uu)
993 {
994   PetscInt adaptIter;
995 
996   PetscFunctionBegin;
997   PetscInt type, limits[5] = {(grid == 0) ? ctx->numRERefine : 0, (grid == 0) ? ctx->nZRefine1 : 0, ctx->numAMRRefine[grid], (grid == 0) ? ctx->nZRefine2 : 0, ctx->postAMRRefine[grid]};
998   for (type = 0; type < 5; type++) {
999     for (adaptIter = 0; adaptIter < limits[type]; adaptIter++) {
1000       DM newForest = NULL;
1001       PetscCall(adaptToleranceFEM(ctx->fe[0], *uu, type, grid, ctx, &newForest));
1002       if (newForest) {
1003         PetscCall(DMDestroy(&ctx->plex[grid]));
1004         PetscCall(VecDestroy(uu));
1005         PetscCall(DMCreateGlobalVector(newForest, uu));
1006         PetscCall(PetscObjectSetName((PetscObject)*uu, "uAMR"));
1007         PetscCall(LandauSetInitialCondition(newForest, *uu, grid, 0, 1, ctx));
1008         ctx->plex[grid] = newForest;
1009       } else {
1010         PetscCall(PetscInfo(*uu, "No refinement\n"));
1011       }
1012     }
1013   }
1014   PetscFunctionReturn(PETSC_SUCCESS);
1015 }
1016 
1017 // make log(Lambdas) from NRL Plasma formulary
1018 static PetscErrorCode makeLambdas(LandauCtx *ctx)
1019 {
1020   PetscFunctionBegin;
1021   for (PetscInt gridi = 0; gridi < ctx->num_grids; gridi++) {
1022     int       iii   = ctx->species_offset[gridi];
1023     PetscReal Ti_ev = (ctx->thermal_temps[iii] / 1.1604525e7) * 1000; // convert (back) to eV
1024     PetscReal ni    = ctx->n[iii] * ctx->n_0;
1025     for (PetscInt gridj = gridi; gridj < ctx->num_grids; gridj++) {
1026       PetscInt  jjj = ctx->species_offset[gridj];
1027       PetscReal Zj  = ctx->charges[jjj] / 1.6022e-19;
1028       if (gridi == 0) {
1029         if (gridj == 0) { // lam_ee
1030           ctx->lambdas[gridi][gridj] = 23.5 - PetscLogReal(PetscSqrtReal(ni) * PetscPowReal(Ti_ev, -1.25)) - PetscSqrtReal(1e-5 + PetscSqr(PetscLogReal(Ti_ev) - 2) / 16);
1031         } else { // lam_ei == lam_ie
1032           if (10 * Zj * Zj > Ti_ev) {
1033             ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 23 - PetscLogReal(PetscSqrtReal(ni) * Zj * PetscPowReal(Ti_ev, -1.5));
1034           } else {
1035             ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 24 - PetscLogReal(PetscSqrtReal(ni) / Ti_ev);
1036           }
1037         }
1038       } else { // lam_ii'
1039         PetscReal mui = ctx->masses[iii] / 1.6720e-27, Zi = ctx->charges[iii] / 1.6022e-19;
1040         PetscReal Tj_ev            = (ctx->thermal_temps[jjj] / 1.1604525e7) * 1000; // convert (back) to eV
1041         PetscReal muj              = ctx->masses[jjj] / 1.6720e-27;
1042         PetscReal nj               = ctx->n[jjj] * ctx->n_0;
1043         ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 23 - PetscLogReal(Zi * Zj * (mui + muj) / (mui * Tj_ev + muj * Ti_ev) * PetscSqrtReal(ni * Zi * Zi / Ti_ev + nj * Zj * Zj / Tj_ev));
1044       }
1045     }
1046   }
1047   //PetscReal v0 = PetscSqrtReal(ctx->k * ctx->thermal_temps[iii] / ctx->masses[iii]); /* arbitrary units for non-dimensionalization: plasma formulary def */
1048   PetscFunctionReturn(PETSC_SUCCESS);
1049 }
1050 
1051 static PetscErrorCode ProcessOptions(LandauCtx *ctx, const char prefix[])
1052 {
1053   PetscBool flg;
1054   PetscInt  ii, nt, nm, nc, num_species_grid[LANDAU_MAX_GRIDS], non_dim_grid;
1055   PetscReal v0_grid[LANDAU_MAX_GRIDS], lnLam = 10;
1056   DM        dummy;
1057 
1058   PetscFunctionBegin;
1059   PetscCall(DMCreate(ctx->comm, &dummy));
1060   /* get options - initialize context */
1061   ctx->verbose        = 1; // should be 0 for silent compliance
1062   ctx->batch_sz       = 1;
1063   ctx->batch_view_idx = 0;
1064   ctx->interpolate    = PETSC_TRUE;
1065   ctx->gpu_assembly   = PETSC_TRUE;
1066   ctx->norm_state     = 0;
1067   ctx->electronShift  = 0;
1068   ctx->M              = NULL;
1069   ctx->J              = NULL;
1070   /* geometry and grids */
1071   ctx->sphere    = PETSC_FALSE;
1072   ctx->use_p4est = PETSC_FALSE;
1073   ctx->simplex   = PETSC_FALSE;
1074   for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) {
1075     ctx->radius[grid]             = 5.; /* thermal radius (velocity) */
1076     ctx->radius_perp[grid]        = 5.; /* thermal radius (velocity) */
1077     ctx->radius_par[grid]         = 5.; /* thermal radius (velocity) */
1078     ctx->numAMRRefine[grid]       = 0;
1079     ctx->postAMRRefine[grid]      = 0;
1080     ctx->species_offset[grid + 1] = 1; // one species default
1081     num_species_grid[grid]        = 0;
1082     ctx->plex[grid]               = NULL; /* cache as expensive to Convert */
1083   }
1084   ctx->species_offset[0] = 0;
1085   ctx->re_radius         = 0.;
1086   ctx->vperp0_radius1    = 0;
1087   ctx->vperp0_radius2    = 0;
1088   ctx->nZRefine1         = 0;
1089   ctx->nZRefine2         = 0;
1090   ctx->numRERefine       = 0;
1091   num_species_grid[0]    = 1; // one species default
1092   /* species - [0] electrons, [1] one ion species eg, duetarium, [2] heavy impurity ion, ... */
1093   ctx->charges[0]       = -1;                       /* electron charge (MKS) */
1094   ctx->masses[0]        = 1 / 1835.469965278441013; /* temporary value in proton mass */
1095   ctx->n[0]             = 1;
1096   ctx->v_0              = 1; /* thermal velocity, we could start with a scale != 1 */
1097   ctx->thermal_temps[0] = 1;
1098   /* constants, etc. */
1099   ctx->epsilon0 = 8.8542e-12;     /* permittivity of free space (MKS) F/m */
1100   ctx->k        = 1.38064852e-23; /* Boltzmann constant (MKS) J/K */
1101   ctx->n_0      = 1.e20;          /* typical plasma n, but could set it to 1 */
1102   ctx->Ez       = 0;
1103   for (PetscInt grid = 0; grid < LANDAU_NUM_TIMERS; grid++) ctx->times[grid] = 0;
1104   for (PetscInt ii = 0; ii < LANDAU_DIM; ii++) ctx->cells0[ii] = 2;
1105   if (LANDAU_DIM == 2) ctx->cells0[0] = 1;
1106   ctx->use_matrix_mass                = PETSC_FALSE;
1107   ctx->use_relativistic_corrections   = PETSC_FALSE;
1108   ctx->use_energy_tensor_trick        = PETSC_FALSE; /* Use Eero's trick for energy conservation v --> grad(v^2/2) */
1109   ctx->SData_d.w                      = NULL;
1110   ctx->SData_d.x                      = NULL;
1111   ctx->SData_d.y                      = NULL;
1112   ctx->SData_d.z                      = NULL;
1113   ctx->SData_d.invJ                   = NULL;
1114   ctx->jacobian_field_major_order     = PETSC_FALSE;
1115   ctx->SData_d.coo_elem_offsets       = NULL;
1116   ctx->SData_d.coo_elem_point_offsets = NULL;
1117   ctx->coo_assembly                   = PETSC_FALSE;
1118   ctx->SData_d.coo_elem_fullNb        = NULL;
1119   ctx->SData_d.coo_size               = 0;
1120   PetscOptionsBegin(ctx->comm, prefix, "Options for Fokker-Plank-Landau collision operator", "none");
1121   {
1122     char opstring[256];
1123 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1124     ctx->deviceType = LANDAU_KOKKOS;
1125     PetscCall(PetscStrncpy(opstring, "kokkos", sizeof(opstring)));
1126 #elif defined(PETSC_HAVE_CUDA)
1127     ctx->deviceType = LANDAU_CUDA;
1128     PetscCall(PetscStrncpy(opstring, "cuda", sizeof(opstring)));
1129 #else
1130     ctx->deviceType = LANDAU_CPU;
1131     PetscCall(PetscStrncpy(opstring, "cpu", sizeof(opstring)));
1132 #endif
1133     PetscCall(PetscOptionsString("-dm_landau_device_type", "Use kernels on 'cpu', 'cuda', or 'kokkos'", "plexland.c", opstring, opstring, sizeof(opstring), NULL));
1134     PetscCall(PetscStrcmp("cpu", opstring, &flg));
1135     if (flg) {
1136       ctx->deviceType = LANDAU_CPU;
1137     } else {
1138       PetscCall(PetscStrcmp("cuda", opstring, &flg));
1139       if (flg) {
1140         ctx->deviceType = LANDAU_CUDA;
1141       } else {
1142         PetscCall(PetscStrcmp("kokkos", opstring, &flg));
1143         if (flg) ctx->deviceType = LANDAU_KOKKOS;
1144         else SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_device_type %s", opstring);
1145       }
1146     }
1147   }
1148   PetscCall(PetscOptionsReal("-dm_landau_electron_shift", "Shift in thermal velocity of electrons", "none", ctx->electronShift, &ctx->electronShift, NULL));
1149   PetscCall(PetscOptionsInt("-dm_landau_verbose", "Level of verbosity output", "plexland.c", ctx->verbose, &ctx->verbose, NULL));
1150   PetscCall(PetscOptionsInt("-dm_landau_batch_size", "Number of 'vertices' to batch", "ex2.c", ctx->batch_sz, &ctx->batch_sz, NULL));
1151   PetscCheck(LANDAU_MAX_BATCH_SZ >= ctx->batch_sz, ctx->comm, PETSC_ERR_ARG_WRONG, "LANDAU_MAX_BATCH_SZ %" PetscInt_FMT " < ctx->batch_sz %" PetscInt_FMT, (PetscInt)LANDAU_MAX_BATCH_SZ, ctx->batch_sz);
1152   PetscCall(PetscOptionsInt("-dm_landau_batch_view_idx", "Index of batch for diagnostics like plotting", "ex2.c", ctx->batch_view_idx, &ctx->batch_view_idx, NULL));
1153   PetscCheck(ctx->batch_view_idx < ctx->batch_sz, ctx->comm, PETSC_ERR_ARG_WRONG, "-ctx->batch_view_idx %" PetscInt_FMT " > ctx->batch_sz %" PetscInt_FMT, ctx->batch_view_idx, ctx->batch_sz);
1154   PetscCall(PetscOptionsReal("-dm_landau_Ez", "Initial parallel electric field in unites of Conner-Hastie critical field", "plexland.c", ctx->Ez, &ctx->Ez, NULL));
1155   PetscCall(PetscOptionsReal("-dm_landau_n_0", "Normalization constant for number density", "plexland.c", ctx->n_0, &ctx->n_0, NULL));
1156   PetscCall(PetscOptionsBool("-dm_landau_use_mataxpy_mass", "Use fast but slightly fragile MATAXPY to add mass term", "plexland.c", ctx->use_matrix_mass, &ctx->use_matrix_mass, NULL));
1157   PetscCall(PetscOptionsBool("-dm_landau_use_relativistic_corrections", "Use relativistic corrections", "plexland.c", ctx->use_relativistic_corrections, &ctx->use_relativistic_corrections, NULL));
1158   PetscCall(PetscOptionsBool("-dm_landau_simplex", "Use simplex elements", "plexland.c", ctx->simplex, &ctx->simplex, NULL));
1159   if (LANDAU_DIM == 2 && ctx->use_relativistic_corrections) ctx->use_relativistic_corrections = PETSC_FALSE; // should warn
1160   PetscCall(PetscOptionsBool("-dm_landau_use_energy_tensor_trick", "Use Eero's trick of using grad(v^2/2) instead of v as args to Landau tensor to conserve energy with relativistic corrections and Q1 elements", "plexland.c", ctx->use_energy_tensor_trick,
1161                              &ctx->use_energy_tensor_trick, NULL));
1162 
1163   /* get num species with temperature, set defaults */
1164   for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) {
1165     ctx->thermal_temps[ii] = 1;
1166     ctx->charges[ii]       = 1;
1167     ctx->masses[ii]        = 1;
1168     ctx->n[ii]             = 1;
1169   }
1170   nt = LANDAU_MAX_SPECIES;
1171   PetscCall(PetscOptionsRealArray("-dm_landau_thermal_temps", "Temperature of each species [e,i_0,i_1,...] in keV (must be set to set number of species)", "plexland.c", ctx->thermal_temps, &nt, &flg));
1172   if (flg) {
1173     PetscCall(PetscInfo(dummy, "num_species set to number of thermal temps provided (%" PetscInt_FMT ")\n", nt));
1174     ctx->num_species = nt;
1175   } else SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_thermal_temps ,t1,t2,.. must be provided to set the number of species");
1176   for (ii = 0; ii < ctx->num_species; ii++) ctx->thermal_temps[ii] *= 1.1604525e7; /* convert to Kelvin */
1177   nm = LANDAU_MAX_SPECIES - 1;
1178   PetscCall(PetscOptionsRealArray("-dm_landau_ion_masses", "Mass of each species in units of proton mass [i_0=2,i_1=40...]", "plexland.c", &ctx->masses[1], &nm, &flg));
1179   PetscCheck(!flg || nm == ctx->num_species - 1, ctx->comm, PETSC_ERR_ARG_WRONG, "num ion masses %" PetscInt_FMT " != num species %" PetscInt_FMT, nm, ctx->num_species - 1);
1180   nm = LANDAU_MAX_SPECIES;
1181   PetscCall(PetscOptionsRealArray("-dm_landau_n", "Number density of each species = n_s * n_0", "plexland.c", ctx->n, &nm, &flg));
1182   PetscCheck(!flg || nm == ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "wrong num n: %" PetscInt_FMT " != num species %" PetscInt_FMT, nm, ctx->num_species);
1183   for (ii = 0; ii < LANDAU_MAX_SPECIES; ii++) ctx->masses[ii] *= 1.6720e-27; /* scale by proton mass kg */
1184   ctx->masses[0] = 9.10938356e-31;                                           /* electron mass kg (should be about right already) */
1185   nc             = LANDAU_MAX_SPECIES - 1;
1186   PetscCall(PetscOptionsRealArray("-dm_landau_ion_charges", "Charge of each species in units of proton charge [i_0=2,i_1=18,...]", "plexland.c", &ctx->charges[1], &nc, &flg));
1187   if (flg) PetscCheck(nc == ctx->num_species - 1, ctx->comm, PETSC_ERR_ARG_WRONG, "num charges %" PetscInt_FMT " != num species %" PetscInt_FMT, nc, ctx->num_species - 1);
1188   for (ii = 0; ii < LANDAU_MAX_SPECIES; ii++) ctx->charges[ii] *= 1.6022e-19; /* electron/proton charge (MKS) */
1189   /* geometry and grids */
1190   nt = LANDAU_MAX_GRIDS;
1191   PetscCall(PetscOptionsIntArray("-dm_landau_num_species_grid", "Number of species on each grid: [ 1, ....] or [S, 0 ....] for single grid", "plexland.c", num_species_grid, &nt, &flg));
1192   if (flg) {
1193     ctx->num_grids = nt;
1194     for (ii = nt = 0; ii < ctx->num_grids; ii++) nt += num_species_grid[ii];
1195     PetscCheck(ctx->num_species == nt, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_num_species_grid: sum %" PetscInt_FMT " != num_species = %" PetscInt_FMT ". %" PetscInt_FMT " grids (check that number of grids <= LANDAU_MAX_GRIDS = %d)", nt, ctx->num_species,
1196                ctx->num_grids, LANDAU_MAX_GRIDS);
1197   } else {
1198     ctx->num_grids      = 1; // go back to a single grid run
1199     num_species_grid[0] = ctx->num_species;
1200   }
1201   for (ctx->species_offset[0] = ii = 0; ii < ctx->num_grids; ii++) ctx->species_offset[ii + 1] = ctx->species_offset[ii] + num_species_grid[ii];
1202   PetscCheck(ctx->species_offset[ctx->num_grids] == ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "ctx->species_offset[ctx->num_grids] %" PetscInt_FMT " != ctx->num_species = %" PetscInt_FMT " ???????????", ctx->species_offset[ctx->num_grids],
1203              ctx->num_species);
1204   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1205     int iii       = ctx->species_offset[grid];                                          // normalize with first (arbitrary) species on grid
1206     v0_grid[grid] = PetscSqrtReal(ctx->k * ctx->thermal_temps[iii] / ctx->masses[iii]); /* arbitrary units for non-dimensionalization: plasma formulary def */
1207   }
1208   // get lambdas here because we need them for t_0 etc
1209   PetscCall(PetscOptionsReal("-dm_landau_ln_lambda", "Universal cross section parameter. Default uses NRL formulas", "plexland.c", lnLam, &lnLam, &flg));
1210   if (flg) {
1211     for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) {
1212       for (PetscInt gridj = 0; gridj < LANDAU_MAX_GRIDS; gridj++) ctx->lambdas[gridj][grid] = lnLam; /* cross section ratio large - small angle collisions */
1213     }
1214   } else {
1215     PetscCall(makeLambdas(ctx));
1216   }
1217   non_dim_grid = 0;
1218   PetscCall(PetscOptionsInt("-dm_landau_normalization_grid", "Index of grid to use for setting v_0, m_0, t_0. (Not recommended)", "plexland.c", non_dim_grid, &non_dim_grid, &flg));
1219   if (non_dim_grid != 0) PetscCall(PetscInfo(dummy, "Normalization grid set to %" PetscInt_FMT ", but non-default not well verified\n", non_dim_grid));
1220   PetscCheck(non_dim_grid >= 0 && non_dim_grid < ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "Normalization grid wrong: %" PetscInt_FMT, non_dim_grid);
1221   ctx->v_0 = v0_grid[non_dim_grid];     /* arbitrary units for non dimensionalization: global mean velocity in 1D of electrons */
1222   ctx->m_0 = ctx->masses[non_dim_grid]; /* arbitrary reference mass, electrons */
1223   ctx->t_0 = 8 * PETSC_PI * PetscSqr(ctx->epsilon0 * ctx->m_0 / PetscSqr(ctx->charges[non_dim_grid])) / ctx->lambdas[non_dim_grid][non_dim_grid] / ctx->n_0 * PetscPowReal(ctx->v_0, 3); /* note, this t_0 makes nu[non_dim_grid,non_dim_grid]=1 */
1224   /* domain */
1225   nt = LANDAU_MAX_GRIDS;
1226   PetscCall(PetscOptionsRealArray("-dm_landau_domain_radius", "Phase space size in units of thermal velocity of grid", "plexland.c", ctx->radius, &nt, &flg));
1227   if (flg) {
1228     PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_radius: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1229     while (nt--) ctx->radius_par[nt] = ctx->radius_perp[nt] = ctx->radius[nt];
1230   } else {
1231     nt = LANDAU_MAX_GRIDS;
1232     PetscCall(PetscOptionsRealArray("-dm_landau_domain_max_par", "Parallel velocity domain size in units of thermal velocity of grid", "plexland.c", ctx->radius_par, &nt, &flg));
1233     if (flg) PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_max_par: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1234     PetscCall(PetscOptionsRealArray("-dm_landau_domain_max_perp", "Perpendicular velocity domain size in units of thermal velocity of grid", "plexland.c", ctx->radius_perp, &nt, &flg));
1235     if (flg) PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_max_perp: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1236   }
1237   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1238     if (flg && ctx->radius[grid] <= 0) { /* negative is ratio of c - need to set par and perp with this -- todo */
1239       if (ctx->radius[grid] == 0) ctx->radius[grid] = 0.75;
1240       else ctx->radius[grid] = -ctx->radius[grid];
1241       ctx->radius[grid] = ctx->radius[grid] * SPEED_OF_LIGHT / ctx->v_0; // use any species on grid to normalize (v_0 same for all on grid)
1242       PetscCall(PetscInfo(dummy, "Change domain radius to %g for grid %" PetscInt_FMT "\n", (double)ctx->radius[grid], grid));
1243     }
1244     ctx->radius[grid] *= v0_grid[grid] / ctx->v_0;      // scale domain by thermal radius relative to v_0
1245     ctx->radius_perp[grid] *= v0_grid[grid] / ctx->v_0; // scale domain by thermal radius relative to v_0
1246     ctx->radius_par[grid] *= v0_grid[grid] / ctx->v_0;  // scale domain by thermal radius relative to v_0
1247   }
1248   /* amr parameters */
1249   nt = LANDAU_MAX_GRIDS;
1250   PetscCall(PetscOptionsIntArray("-dm_landau_amr_levels_max", "Number of AMR levels of refinement around origin, after (RE) refinements along z", "plexland.c", ctx->numAMRRefine, &nt, &flg));
1251   PetscCheck(!flg || nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_amr_levels_max: given %" PetscInt_FMT " != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1252   nt = LANDAU_MAX_GRIDS;
1253   PetscCall(PetscOptionsIntArray("-dm_landau_amr_post_refine", "Number of levels to uniformly refine after AMR", "plexland.c", ctx->postAMRRefine, &nt, &flg));
1254   for (ii = 1; ii < ctx->num_grids; ii++) ctx->postAMRRefine[ii] = ctx->postAMRRefine[0]; // all grids the same now
1255   PetscCall(PetscOptionsInt("-dm_landau_amr_re_levels", "Number of levels to refine along v_perp=0, z>0", "plexland.c", ctx->numRERefine, &ctx->numRERefine, &flg));
1256   PetscCall(PetscOptionsInt("-dm_landau_amr_z_refine_pre", "Number of levels to refine along v_perp=0 before origin refine", "plexland.c", ctx->nZRefine1, &ctx->nZRefine1, &flg));
1257   PetscCall(PetscOptionsInt("-dm_landau_amr_z_refine_post", "Number of levels to refine along v_perp=0 after origin refine", "plexland.c", ctx->nZRefine2, &ctx->nZRefine2, &flg));
1258   PetscCall(PetscOptionsReal("-dm_landau_re_radius", "velocity range to refine on positive (z>0) r=0 axis for runaways", "plexland.c", ctx->re_radius, &ctx->re_radius, &flg));
1259   PetscCall(PetscOptionsReal("-dm_landau_z_radius_pre", "velocity range to refine r=0 axis (for electrons)", "plexland.c", ctx->vperp0_radius1, &ctx->vperp0_radius1, &flg));
1260   PetscCall(PetscOptionsReal("-dm_landau_z_radius_post", "velocity range to refine r=0 axis (for electrons) after origin AMR", "plexland.c", ctx->vperp0_radius2, &ctx->vperp0_radius2, &flg));
1261   /* spherical domain (not used) */
1262   PetscCall(PetscOptionsBool("-dm_landau_sphere", "use sphere/semi-circle domain instead of rectangle", "plexland.c", ctx->sphere, &ctx->sphere, NULL));
1263   if (ctx->sphere || ctx->simplex) {
1264     ctx->sphere_inner_radius_90degree = 0.40;
1265     ctx->sphere_inner_radius_45degree = 0.35;
1266     PetscCall(PetscOptionsReal("-dm_landau_sphere_inner_radius_90degree_scale", "Scaling of radius for inner circle on 90 degree grid", "plexland.c", ctx->sphere_inner_radius_90degree, &ctx->sphere_inner_radius_90degree, NULL));
1267     PetscCall(PetscOptionsReal("-dm_landau_sphere_inner_radius_45degree_scale", "Scaling of radius for inner circle on 45 degree grid", "plexland.c", ctx->sphere_inner_radius_45degree, &ctx->sphere_inner_radius_45degree, NULL));
1268   } else {
1269     nt = LANDAU_DIM;
1270     PetscCall(PetscOptionsIntArray("-dm_landau_num_cells", "Number of cells in each dimension of base grid", "plexland.c", ctx->cells0, &nt, &flg));
1271   }
1272   /* processing options */
1273   PetscCall(PetscOptionsBool("-dm_landau_gpu_assembly", "Assemble Jacobian on GPU", "plexland.c", ctx->gpu_assembly, &ctx->gpu_assembly, NULL));
1274   if (ctx->deviceType == LANDAU_CPU || ctx->deviceType == LANDAU_KOKKOS) { // make Kokkos
1275     PetscCall(PetscOptionsBool("-dm_landau_coo_assembly", "Assemble Jacobian with Kokkos on 'device'", "plexland.c", ctx->coo_assembly, &ctx->coo_assembly, NULL));
1276     if (ctx->coo_assembly) PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "COO assembly requires 'gpu assembly' even if Kokkos 'CPU' back-end %d", ctx->coo_assembly);
1277   }
1278   PetscCall(PetscOptionsBool("-dm_landau_jacobian_field_major_order", "Reorder Jacobian for GPU assembly with field major, or block diagonal, ordering (DEPRECATED)", "plexland.c", ctx->jacobian_field_major_order, &ctx->jacobian_field_major_order, NULL));
1279   if (ctx->jacobian_field_major_order) PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order requires -dm_landau_gpu_assembly");
1280   PetscCheck(!ctx->jacobian_field_major_order, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order DEPRECATED");
1281   PetscOptionsEnd();
1282 
1283   for (ii = ctx->num_species; ii < LANDAU_MAX_SPECIES; ii++) ctx->masses[ii] = ctx->thermal_temps[ii] = ctx->charges[ii] = 0;
1284   if (ctx->verbose != 0) {
1285     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "masses:        e=%10.3e; ions in proton mass units:   %10.3e %10.3e ...\n", (double)ctx->masses[0], (double)(ctx->masses[1] / 1.6720e-27), (double)(ctx->num_species > 2 ? ctx->masses[2] / 1.6720e-27 : 0)));
1286     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "charges:       e=%10.3e; charges in elementary units: %10.3e %10.3e\n", (double)ctx->charges[0], (double)(-ctx->charges[1] / ctx->charges[0]), (double)(ctx->num_species > 2 ? -ctx->charges[2] / ctx->charges[0] : 0)));
1287     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "n:             e: %10.3e                           i: %10.3e %10.3e\n", (double)ctx->n[0], (double)ctx->n[1], (double)(ctx->num_species > 2 ? ctx->n[2] : 0)));
1288     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "thermal T (K): e=%10.3e i=%10.3e %10.3e. Normalization grid %d: v_0=%10.3e (%10.3ec) n_0=%10.3e t_0=%10.3e %" PetscInt_FMT " batched, view batch %" PetscInt_FMT "\n", (double)ctx->thermal_temps[0],
1289                           (double)ctx->thermal_temps[1], (double)((ctx->num_species > 2) ? ctx->thermal_temps[2] : 0), (int)non_dim_grid, (double)ctx->v_0, (double)(ctx->v_0 / SPEED_OF_LIGHT), (double)ctx->n_0, (double)ctx->t_0, ctx->batch_sz, ctx->batch_view_idx));
1290     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain radius (AMR levels) grid %d: par=%10.3e perp=%10.3e (%" PetscInt_FMT ") ", 0, (double)ctx->radius_par[0], (double)ctx->radius_perp[0], ctx->numAMRRefine[0]));
1291     for (ii = 1; ii < ctx->num_grids; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %" PetscInt_FMT ": par=%10.3e perp=%10.3e (%" PetscInt_FMT ") ", ii, (double)ctx->radius_par[ii], (double)ctx->radius_perp[ii], ctx->numAMRRefine[ii]));
1292     if (ctx->use_relativistic_corrections) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nUse relativistic corrections\n"));
1293     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1294   }
1295   PetscCall(DMDestroy(&dummy));
1296   {
1297     PetscMPIInt rank;
1298     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1299     ctx->stage = 0;
1300     PetscCall(PetscLogEventRegister("Landau Create", DM_CLASSID, &ctx->events[13]));   /* 13 */
1301     PetscCall(PetscLogEventRegister(" GPU ass. setup", DM_CLASSID, &ctx->events[2]));  /* 2 */
1302     PetscCall(PetscLogEventRegister(" Build matrix", DM_CLASSID, &ctx->events[12]));   /* 12 */
1303     PetscCall(PetscLogEventRegister(" Assembly maps", DM_CLASSID, &ctx->events[15]));  /* 15 */
1304     PetscCall(PetscLogEventRegister("Landau Mass mat", DM_CLASSID, &ctx->events[14])); /* 14 */
1305     PetscCall(PetscLogEventRegister("Landau Operator", DM_CLASSID, &ctx->events[11])); /* 11 */
1306     PetscCall(PetscLogEventRegister("Landau Jacobian", DM_CLASSID, &ctx->events[0]));  /* 0 */
1307     PetscCall(PetscLogEventRegister("Landau Mass", DM_CLASSID, &ctx->events[9]));      /* 9 */
1308     PetscCall(PetscLogEventRegister(" Preamble", DM_CLASSID, &ctx->events[10]));       /* 10 */
1309     PetscCall(PetscLogEventRegister(" static IP Data", DM_CLASSID, &ctx->events[7]));  /* 7 */
1310     PetscCall(PetscLogEventRegister(" dynamic IP-Jac", DM_CLASSID, &ctx->events[1]));  /* 1 */
1311     PetscCall(PetscLogEventRegister(" Kernel-init", DM_CLASSID, &ctx->events[3]));     /* 3 */
1312     PetscCall(PetscLogEventRegister(" Jac-f-df (GPU)", DM_CLASSID, &ctx->events[8]));  /* 8 */
1313     PetscCall(PetscLogEventRegister(" J Kernel (GPU)", DM_CLASSID, &ctx->events[4]));  /* 4 */
1314     PetscCall(PetscLogEventRegister(" M Kernel (GPU)", DM_CLASSID, &ctx->events[16])); /* 16 */
1315     PetscCall(PetscLogEventRegister(" Copy to CPU", DM_CLASSID, &ctx->events[5]));     /* 5 */
1316     PetscCall(PetscLogEventRegister(" CPU assemble", DM_CLASSID, &ctx->events[6]));    /* 6 */
1317 
1318     if (rank) { /* turn off output stuff for duplicate runs - do we need to add the prefix to all this? */
1319       PetscCall(PetscOptionsClearValue(NULL, "-snes_converged_reason"));
1320       PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
1321       PetscCall(PetscOptionsClearValue(NULL, "-snes_monitor"));
1322       PetscCall(PetscOptionsClearValue(NULL, "-ksp_monitor"));
1323       PetscCall(PetscOptionsClearValue(NULL, "-ts_monitor"));
1324       PetscCall(PetscOptionsClearValue(NULL, "-ts_view"));
1325       PetscCall(PetscOptionsClearValue(NULL, "-ts_adapt_monitor"));
1326       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_amr_dm_view"));
1327       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_amr_vec_view"));
1328       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mass_dm_view"));
1329       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mass_view"));
1330       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_jacobian_view"));
1331       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mat_view"));
1332       PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_converged_reason"));
1333       PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_monitor"));
1334       PetscCall(PetscOptionsClearValue(NULL, "-"));
1335       PetscCall(PetscOptionsClearValue(NULL, "-info"));
1336     }
1337   }
1338   PetscFunctionReturn(PETSC_SUCCESS);
1339 }
1340 
1341 static PetscErrorCode CreateStaticData(PetscInt dim, IS grid_batch_is_inv[], LandauCtx *ctx)
1342 {
1343   PetscSection     section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS];
1344   PetscQuadrature  quad;
1345   const PetscReal *quadWeights;
1346   PetscReal        invMass[LANDAU_MAX_SPECIES], nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES];
1347   PetscInt         numCells[LANDAU_MAX_GRIDS], Nq, Nf[LANDAU_MAX_GRIDS], ncellsTot = 0, MAP_BF_SIZE = 64 * LANDAU_DIM * LANDAU_DIM * LANDAU_MAX_Q_FACE * LANDAU_MAX_SPECIES;
1348   PetscTabulation *Tf;
1349   PetscDS          prob;
1350 
1351   PetscFunctionBegin;
1352   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1353     for (PetscInt ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++) {
1354       invMass[ii]  = ctx->m_0 / ctx->masses[ii];
1355       nu_alpha[ii] = PetscSqr(ctx->charges[ii] / ctx->m_0) * ctx->m_0 / ctx->masses[ii];
1356       nu_beta[ii]  = PetscSqr(ctx->charges[ii] / ctx->epsilon0) / (8 * PETSC_PI) * ctx->t_0 * ctx->n_0 / PetscPowReal(ctx->v_0, 3);
1357     }
1358   }
1359   if (ctx->verbose == 4) {
1360     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nu_alpha: "));
1361     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1362       int iii = ctx->species_offset[grid];
1363       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %e", (double)nu_alpha[ii]));
1364     }
1365     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nnu_beta: "));
1366     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1367       int iii = ctx->species_offset[grid];
1368       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %e", (double)nu_beta[ii]));
1369     }
1370     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nnu_alpha[i]*nu_beta[j]*lambda[i][j]:\n"));
1371     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1372       int iii = ctx->species_offset[grid];
1373       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) {
1374         for (PetscInt gridj = 0; gridj < ctx->num_grids; gridj++) {
1375           int jjj = ctx->species_offset[gridj];
1376           for (PetscInt jj = jjj; jj < ctx->species_offset[gridj + 1]; jj++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %14.9e", (double)(nu_alpha[ii] * nu_beta[jj] * ctx->lambdas[grid][gridj])));
1377         }
1378         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1379       }
1380     }
1381     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "lambda[i][j]:\n"));
1382     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1383       int iii = ctx->species_offset[grid];
1384       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) {
1385         for (PetscInt gridj = 0; gridj < ctx->num_grids; gridj++) {
1386           int jjj = ctx->species_offset[gridj];
1387           for (PetscInt jj = jjj; jj < ctx->species_offset[gridj + 1]; jj++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %14.9e", (double)ctx->lambdas[grid][gridj]));
1388         }
1389         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1390       }
1391     }
1392   }
1393   PetscCall(DMGetDS(ctx->plex[0], &prob));    // same DS for all grids
1394   PetscCall(PetscDSGetTabulation(prob, &Tf)); // Bf, &Df same for all grids
1395   /* DS, Tab and quad is same on all grids */
1396   PetscCheck(ctx->plex[0], ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
1397   PetscCall(PetscFEGetQuadrature(ctx->fe[0], &quad));
1398   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, &quadWeights));
1399   PetscCheck(Nq <= LANDAU_MAX_NQ, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQ (%d)", Nq, LANDAU_MAX_NQ);
1400   /* setup each grid */
1401   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1402     PetscInt cStart, cEnd;
1403     PetscCheck(ctx->plex[grid] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
1404     PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
1405     numCells[grid] = cEnd - cStart; // grids can have different topology
1406     PetscCall(DMGetLocalSection(ctx->plex[grid], &section[grid]));
1407     PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid]));
1408     PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid]));
1409     ncellsTot += numCells[grid];
1410   }
1411   /* create GPU assembly data */
1412   if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
1413     PetscContainer container;
1414     PetscScalar   *elemMatrix, *elMat;
1415     pointInterpolationP4est(*pointMaps)[LANDAU_MAX_Q_FACE];
1416     P4estVertexMaps *maps;
1417     const PetscInt  *plex_batch = NULL, Nb = Nq, elMatSz = Nq * Nq * ctx->num_species * ctx->num_species; // tensor elements;
1418     LandauIdx       *coo_elem_offsets = NULL, *coo_elem_fullNb = NULL, (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = NULL;
1419     /* create GPU assembly data */
1420     PetscCall(PetscInfo(ctx->plex[0], "Make GPU maps %d\n", 1));
1421     PetscCall(PetscLogEventBegin(ctx->events[2], 0, 0, 0, 0));
1422     PetscCall(PetscMalloc(sizeof(*maps) * ctx->num_grids, &maps));
1423     PetscCall(PetscMalloc(sizeof(*pointMaps) * MAP_BF_SIZE, &pointMaps));
1424     PetscCall(PetscMalloc(sizeof(*elemMatrix) * elMatSz, &elemMatrix));
1425 
1426     if (ctx->coo_assembly) {                                                                                                      // setup COO assembly -- put COO metadata directly in ctx->SData_d
1427       PetscCall(PetscMalloc3(ncellsTot + 1, &coo_elem_offsets, ncellsTot, &coo_elem_fullNb, ncellsTot, &coo_elem_point_offsets)); // array of integer pointers
1428       coo_elem_offsets[0] = 0;                                                                                                    // finish later
1429       PetscCall(PetscInfo(ctx->plex[0], "COO initialization, %" PetscInt_FMT " cells\n", ncellsTot));
1430       ctx->SData_d.coo_n_cellsTot         = ncellsTot;
1431       ctx->SData_d.coo_elem_offsets       = (void *)coo_elem_offsets;
1432       ctx->SData_d.coo_elem_fullNb        = (void *)coo_elem_fullNb;
1433       ctx->SData_d.coo_elem_point_offsets = (void *)coo_elem_point_offsets;
1434     } else {
1435       ctx->SData_d.coo_elem_offsets = ctx->SData_d.coo_elem_fullNb = NULL;
1436       ctx->SData_d.coo_elem_point_offsets                          = NULL;
1437       ctx->SData_d.coo_n_cellsTot                                  = 0;
1438     }
1439 
1440     ctx->SData_d.coo_max_fullnb = 0;
1441     for (PetscInt grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) {
1442       PetscInt cStart, cEnd, Nfloc = Nf[grid], totDim = Nfloc * Nq;
1443       if (grid_batch_is_inv[grid]) PetscCall(ISGetIndices(grid_batch_is_inv[grid], &plex_batch));
1444       PetscCheck(!plex_batch, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order DEPRECATED");
1445       PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
1446       // make maps
1447       maps[grid].d_self       = NULL;
1448       maps[grid].num_elements = numCells[grid];
1449       maps[grid].num_face     = (PetscInt)(pow(Nq, 1. / ((double)dim)) + .001);                 // Q
1450       maps[grid].num_face     = (PetscInt)(pow(maps[grid].num_face, (double)(dim - 1)) + .001); // Q^2
1451       maps[grid].num_reduced  = 0;
1452       maps[grid].deviceType   = ctx->deviceType;
1453       maps[grid].numgrids     = ctx->num_grids;
1454       // count reduced and get
1455       PetscCall(PetscMalloc(maps[grid].num_elements * sizeof(*maps[grid].gIdx), &maps[grid].gIdx));
1456       for (int ej = cStart, eidx = 0; ej < cEnd; ++ej, ++eidx, glb_elem_idx++) {
1457         if (coo_elem_offsets) coo_elem_offsets[glb_elem_idx + 1] = coo_elem_offsets[glb_elem_idx]; // start with last one, then add
1458         for (int fieldA = 0; fieldA < Nf[grid]; fieldA++) {
1459           int fullNb = 0;
1460           for (int q = 0; q < Nb; ++q) {
1461             PetscInt     numindices, *indices;
1462             PetscScalar *valuesOrig = elMat = elemMatrix;
1463             PetscCall(PetscArrayzero(elMat, totDim * totDim));
1464             elMat[(fieldA * Nb + q) * totDim + fieldA * Nb + q] = 1;
1465             PetscCall(DMPlexGetClosureIndices(ctx->plex[grid], section[grid], globsection[grid], ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **)&elMat));
1466             for (PetscInt f = 0; f < numindices; ++f) { // look for a non-zero on the diagonal (is this too complicated for simplices?)
1467               if (PetscAbs(PetscRealPart(elMat[f * numindices + f])) > PETSC_MACHINE_EPSILON) {
1468                 // found it
1469                 if (PetscAbs(PetscRealPart(elMat[f * numindices + f] - 1.)) < PETSC_MACHINE_EPSILON) { // normal vertex 1.0
1470                   if (plex_batch) {
1471                     maps[grid].gIdx[eidx][fieldA][q] = (LandauIdx)plex_batch[indices[f]];
1472                   } else {
1473                     maps[grid].gIdx[eidx][fieldA][q] = (LandauIdx)indices[f];
1474                   }
1475                   fullNb++;
1476                 } else { //found a constraint
1477                   int            jj                = 0;
1478                   PetscReal      sum               = 0;
1479                   const PetscInt ff                = f;
1480                   maps[grid].gIdx[eidx][fieldA][q] = -maps[grid].num_reduced - 1; // store (-)index: id = -(idx+1): idx = -id - 1
1481                   PetscCheck(!ctx->simplex, ctx->comm, PETSC_ERR_ARG_WRONG, "No constraints with simplex");
1482                   do {                                                                                              // constraints are continuous in Plex - exploit that here
1483                     int ii;                                                                                         // get 'scale'
1484                     for (ii = 0, pointMaps[maps[grid].num_reduced][jj].scale = 0; ii < maps[grid].num_face; ii++) { // sum row of outer product to recover vector value
1485                       if (ff + ii < numindices) {                                                                   // 3D has Q and Q^2 interps so might run off end. We could test that elMat[f*numindices + ff + ii] > 0, and break if not
1486                         pointMaps[maps[grid].num_reduced][jj].scale += PetscRealPart(elMat[f * numindices + ff + ii]);
1487                       }
1488                     }
1489                     sum += pointMaps[maps[grid].num_reduced][jj].scale; // diagnostic
1490                     // get 'gid'
1491                     if (pointMaps[maps[grid].num_reduced][jj].scale == 0) pointMaps[maps[grid].num_reduced][jj].gid = -1; // 3D has Q and Q^2 interps
1492                     else {
1493                       if (plex_batch) {
1494                         pointMaps[maps[grid].num_reduced][jj].gid = plex_batch[indices[f]];
1495                       } else {
1496                         pointMaps[maps[grid].num_reduced][jj].gid = indices[f];
1497                       }
1498                       fullNb++;
1499                     }
1500                   } while (++jj < maps[grid].num_face && ++f < numindices); // jj is incremented if we hit the end
1501                   while (jj < maps[grid].num_face) {
1502                     pointMaps[maps[grid].num_reduced][jj].scale = 0;
1503                     pointMaps[maps[grid].num_reduced][jj].gid   = -1;
1504                     jj++;
1505                   }
1506                   if (PetscAbs(sum - 1.0) > 10 * PETSC_MACHINE_EPSILON) { // debug
1507                     int       d, f;
1508                     PetscReal tmp = 0;
1509                     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\t\t%d.%d.%d) ERROR total I = %22.16e (LANDAU_MAX_Q_FACE=%d, #face=%d)\n", eidx, q, fieldA, (double)sum, LANDAU_MAX_Q_FACE, maps[grid].num_face));
1510                     for (d = 0, tmp = 0; d < numindices; ++d) {
1511                       if (tmp != 0 && PetscAbs(tmp - 1.0) > 10 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3d) %3" PetscInt_FMT ": ", d, indices[d]));
1512                       for (f = 0; f < numindices; ++f) tmp += PetscRealPart(elMat[d * numindices + f]);
1513                       if (tmp != 0) PetscCall(PetscPrintf(ctx->comm, " | %22.16e\n", (double)tmp));
1514                     }
1515                   }
1516                   maps[grid].num_reduced++;
1517                   PetscCheck(maps[grid].num_reduced < MAP_BF_SIZE, PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps[grid].num_reduced %d > %" PetscInt_FMT, maps[grid].num_reduced, MAP_BF_SIZE);
1518                 }
1519                 break;
1520               }
1521             }
1522             // cleanup
1523             PetscCall(DMPlexRestoreClosureIndices(ctx->plex[grid], section[grid], globsection[grid], ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **)&elMat));
1524             if (elMat != valuesOrig) PetscCall(DMRestoreWorkArray(ctx->plex[grid], numindices * numindices, MPIU_SCALAR, &elMat));
1525           }
1526           if (ctx->coo_assembly) {                                 // setup COO assembly
1527             coo_elem_offsets[glb_elem_idx + 1] += fullNb * fullNb; // one species block, adds a block for each species, on this element in this grid
1528             if (fieldA == 0) {                                     // cache full Nb for this element, on this grid per species
1529               coo_elem_fullNb[glb_elem_idx] = fullNb;
1530               if (fullNb > ctx->SData_d.coo_max_fullnb) ctx->SData_d.coo_max_fullnb = fullNb;
1531             } else PetscCheck(coo_elem_fullNb[glb_elem_idx] == fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "full element size change with species %d %d", coo_elem_fullNb[glb_elem_idx], fullNb);
1532           }
1533         } // field
1534       }   // cell
1535       // allocate and copy point data maps[grid].gIdx[eidx][field][q]
1536       PetscCall(PetscMalloc(maps[grid].num_reduced * sizeof(*maps[grid].c_maps), &maps[grid].c_maps));
1537       for (int ej = 0; ej < maps[grid].num_reduced; ++ej) {
1538         for (int q = 0; q < maps[grid].num_face; ++q) {
1539           maps[grid].c_maps[ej][q].scale = pointMaps[ej][q].scale;
1540           maps[grid].c_maps[ej][q].gid   = pointMaps[ej][q].gid;
1541         }
1542       }
1543 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1544       if (ctx->deviceType == LANDAU_KOKKOS) {
1545         PetscCall(LandauKokkosCreateMatMaps(maps, pointMaps, Nf, Nq, grid)); // implies Kokkos does
1546       }                                                                      // else could be CUDA
1547 #endif
1548 #if defined(PETSC_HAVE_CUDA)
1549       if (ctx->deviceType == LANDAU_CUDA) PetscCall(LandauCUDACreateMatMaps(maps, pointMaps, Nf, Nq, grid));
1550 #endif
1551       if (plex_batch) {
1552         PetscCall(ISRestoreIndices(grid_batch_is_inv[grid], &plex_batch));
1553         PetscCall(ISDestroy(&grid_batch_is_inv[grid])); // we are done with this
1554       }
1555     } /* grids */
1556     // finish COO
1557     if (ctx->coo_assembly) { // setup COO assembly
1558       PetscInt *oor, *ooc;
1559       ctx->SData_d.coo_size = coo_elem_offsets[ncellsTot] * ctx->batch_sz;
1560       PetscCall(PetscMalloc2(ctx->SData_d.coo_size, &oor, ctx->SData_d.coo_size, &ooc));
1561       for (int i = 0; i < ctx->SData_d.coo_size; i++) oor[i] = ooc[i] = -1;
1562       // get
1563       for (int grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) {
1564         for (int ej = 0; ej < numCells[grid]; ++ej, glb_elem_idx++) {
1565           const int              fullNb           = coo_elem_fullNb[glb_elem_idx];
1566           const LandauIdx *const Idxs             = &maps[grid].gIdx[ej][0][0]; // just use field-0 maps, They should be the same but this is just for COO storage
1567           coo_elem_point_offsets[glb_elem_idx][0] = 0;
1568           for (int f = 0, cnt2 = 0; f < Nb; f++) {
1569             int idx                                     = Idxs[f];
1570             coo_elem_point_offsets[glb_elem_idx][f + 1] = coo_elem_point_offsets[glb_elem_idx][f]; // start at last
1571             if (idx >= 0) {
1572               cnt2++;
1573               coo_elem_point_offsets[glb_elem_idx][f + 1]++; // inc
1574             } else {
1575               idx = -idx - 1;
1576               for (int q = 0; q < maps[grid].num_face; q++) {
1577                 if (maps[grid].c_maps[idx][q].gid < 0) break;
1578                 cnt2++;
1579                 coo_elem_point_offsets[glb_elem_idx][f + 1]++; // inc
1580               }
1581             }
1582             PetscCheck(cnt2 <= fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "wrong count %d < %d", fullNb, cnt2);
1583           }
1584           PetscCheck(coo_elem_point_offsets[glb_elem_idx][Nb] == fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "coo_elem_point_offsets size %d != fullNb=%d", coo_elem_point_offsets[glb_elem_idx][Nb], fullNb);
1585         }
1586       }
1587       // set
1588       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
1589         for (int grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) {
1590           const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
1591           for (int ej = 0; ej < numCells[grid]; ++ej, glb_elem_idx++) {
1592             const int fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb;
1593             // set (i,j)
1594             for (int fieldA = 0; fieldA < Nf[grid]; fieldA++) {
1595               const LandauIdx *const Idxs = &maps[grid].gIdx[ej][fieldA][0];
1596               int                    rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE];
1597               for (int f = 0; f < Nb; ++f) {
1598                 const int nr = coo_elem_point_offsets[glb_elem_idx][f + 1] - coo_elem_point_offsets[glb_elem_idx][f];
1599                 if (nr == 1) rows[0] = Idxs[f];
1600                 else {
1601                   const int idx = -Idxs[f] - 1;
1602                   for (int q = 0; q < nr; q++) rows[q] = maps[grid].c_maps[idx][q].gid;
1603                 }
1604                 for (int g = 0; g < Nb; ++g) {
1605                   const int nc = coo_elem_point_offsets[glb_elem_idx][g + 1] - coo_elem_point_offsets[glb_elem_idx][g];
1606                   if (nc == 1) cols[0] = Idxs[g];
1607                   else {
1608                     const int idx = -Idxs[g] - 1;
1609                     for (int q = 0; q < nc; q++) cols[q] = maps[grid].c_maps[idx][q].gid;
1610                   }
1611                   const int idx0 = b_id * coo_elem_offsets[ncellsTot] + coo_elem_offsets[glb_elem_idx] + fieldA * fullNb2 + fullNb * coo_elem_point_offsets[glb_elem_idx][f] + nr * coo_elem_point_offsets[glb_elem_idx][g];
1612                   for (int q = 0, idx = idx0; q < nr; q++) {
1613                     for (int d = 0; d < nc; d++, idx++) {
1614                       oor[idx] = rows[q] + moffset;
1615                       ooc[idx] = cols[d] + moffset;
1616                     }
1617                   }
1618                 }
1619               }
1620             }
1621           } // cell
1622         }   // grid
1623       }     // batch
1624       PetscCall(MatSetPreallocationCOO(ctx->J, ctx->SData_d.coo_size, oor, ooc));
1625       PetscCall(PetscFree2(oor, ooc));
1626     }
1627     PetscCall(PetscFree(pointMaps));
1628     PetscCall(PetscFree(elemMatrix));
1629     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
1630     PetscCall(PetscContainerSetPointer(container, (void *)maps));
1631     PetscCall(PetscContainerSetUserDestroy(container, LandauGPUMapsDestroy));
1632     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "assembly_maps", (PetscObject)container));
1633     PetscCall(PetscContainerDestroy(&container));
1634     PetscCall(PetscLogEventEnd(ctx->events[2], 0, 0, 0, 0));
1635   } // end GPU assembly
1636   { /* create static point data, Jacobian called first, only one vertex copy */
1637     PetscReal     *invJe, *ww, *xx, *yy, *zz = NULL, *invJ_a;
1638     PetscInt       outer_ipidx, outer_ej, grid, nip_glb = 0;
1639     PetscFE        fe;
1640     const PetscInt Nb = Nq;
1641     PetscCall(PetscLogEventBegin(ctx->events[7], 0, 0, 0, 0));
1642     PetscCall(PetscInfo(ctx->plex[0], "Initialize static data\n"));
1643     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) nip_glb += Nq * numCells[grid];
1644     /* collect f data, first time is for Jacobian, but make mass now */
1645     if (ctx->verbose != 0) {
1646       PetscInt ncells = 0, N;
1647       PetscCall(MatGetSize(ctx->J, &N, NULL));
1648       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ncells += numCells[grid];
1649       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d) %s %" PetscInt_FMT " IPs, %" PetscInt_FMT " cells total, Nb=%" PetscInt_FMT ", Nq=%" PetscInt_FMT ", dim=%" PetscInt_FMT ", Tab: Nb=%" PetscInt_FMT " Nf=%" PetscInt_FMT " Np=%" PetscInt_FMT " cdim=%" PetscInt_FMT " N=%" PetscInt_FMT "\n", 0, "FormLandau", nip_glb, ncells, Nb, Nq, dim, Nb,
1650                             ctx->num_species, Nb, dim, N));
1651     }
1652     PetscCall(PetscMalloc4(nip_glb, &ww, nip_glb, &xx, nip_glb, &yy, nip_glb * dim * dim, &invJ_a));
1653     if (dim == 3) PetscCall(PetscMalloc1(nip_glb, &zz));
1654     if (ctx->use_energy_tensor_trick) {
1655       PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, ctx->simplex, NULL, PETSC_DECIDE, &fe));
1656       PetscCall(PetscObjectSetName((PetscObject)fe, "energy"));
1657     }
1658     /* init each grids static data - no batch */
1659     for (grid = 0, outer_ipidx = 0, outer_ej = 0; grid < ctx->num_grids; grid++) { // OpenMP (once)
1660       Vec          v2_2 = NULL;                                                    // projected function: v^2/2 for non-relativistic, gamma... for relativistic
1661       PetscSection e_section;
1662       DM           dmEnergy;
1663       PetscInt     cStart, cEnd, ej;
1664 
1665       PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
1666       // prep energy trick, get v^2 / 2 vector
1667       if (ctx->use_energy_tensor_trick) {
1668         PetscErrorCode (*energyf[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {ctx->use_relativistic_corrections ? gamma_m1_f : energy_f};
1669         Vec        glob_v2;
1670         PetscReal *c2_0[1], data[1] = {PetscSqr(C_0(ctx->v_0))};
1671 
1672         PetscCall(DMClone(ctx->plex[grid], &dmEnergy));
1673         PetscCall(PetscObjectSetName((PetscObject)dmEnergy, "energy"));
1674         PetscCall(DMSetField(dmEnergy, 0, NULL, (PetscObject)fe));
1675         PetscCall(DMCreateDS(dmEnergy));
1676         PetscCall(DMGetSection(dmEnergy, &e_section));
1677         PetscCall(DMGetGlobalVector(dmEnergy, &glob_v2));
1678         PetscCall(PetscObjectSetName((PetscObject)glob_v2, "trick"));
1679         c2_0[0] = &data[0];
1680         PetscCall(DMProjectFunction(dmEnergy, 0., energyf, (void **)c2_0, INSERT_ALL_VALUES, glob_v2));
1681         PetscCall(DMGetLocalVector(dmEnergy, &v2_2));
1682         PetscCall(VecZeroEntries(v2_2)); /* zero BCs so don't set */
1683         PetscCall(DMGlobalToLocalBegin(dmEnergy, glob_v2, INSERT_VALUES, v2_2));
1684         PetscCall(DMGlobalToLocalEnd(dmEnergy, glob_v2, INSERT_VALUES, v2_2));
1685         PetscCall(DMViewFromOptions(dmEnergy, NULL, "-energy_dm_view"));
1686         PetscCall(VecViewFromOptions(glob_v2, NULL, "-energy_vec_view"));
1687         PetscCall(DMRestoreGlobalVector(dmEnergy, &glob_v2));
1688       }
1689       /* append part of the IP data for each grid */
1690       for (ej = 0; ej < numCells[grid]; ++ej, ++outer_ej) {
1691         PetscScalar *coefs = NULL;
1692         PetscReal    vj[LANDAU_MAX_NQ * LANDAU_DIM], detJj[LANDAU_MAX_NQ], Jdummy[LANDAU_MAX_NQ * LANDAU_DIM * LANDAU_DIM], c0 = C_0(ctx->v_0), c02 = PetscSqr(c0);
1693         invJe = invJ_a + outer_ej * Nq * dim * dim;
1694         PetscCall(DMPlexComputeCellGeometryFEM(ctx->plex[grid], ej + cStart, quad, vj, Jdummy, invJe, detJj));
1695         if (ctx->use_energy_tensor_trick) PetscCall(DMPlexVecGetClosure(dmEnergy, e_section, v2_2, ej + cStart, NULL, &coefs));
1696         /* create static point data */
1697         for (PetscInt qj = 0; qj < Nq; qj++, outer_ipidx++) {
1698           const PetscInt   gidx = outer_ipidx;
1699           const PetscReal *invJ = &invJe[qj * dim * dim];
1700           ww[gidx]              = detJj[qj] * quadWeights[qj];
1701           if (dim == 2) ww[gidx] *= vj[qj * dim + 0]; /* cylindrical coordinate, w/o 2pi */
1702           // get xx, yy, zz
1703           if (ctx->use_energy_tensor_trick) {
1704             double                 refSpaceDer[3], eGradPhi[3];
1705             const PetscReal *const DD = Tf[0]->T[1];
1706             const PetscReal       *Dq = &DD[qj * Nb * dim];
1707             for (int d = 0; d < 3; ++d) refSpaceDer[d] = eGradPhi[d] = 0.0;
1708             for (int b = 0; b < Nb; ++b) {
1709               for (int d = 0; d < dim; ++d) refSpaceDer[d] += Dq[b * dim + d] * PetscRealPart(coefs[b]);
1710             }
1711             xx[gidx] = 1e10;
1712             if (ctx->use_relativistic_corrections) {
1713               double dg2_c2 = 0;
1714               //for (int d = 0; d < dim; ++d) refSpaceDer[d] *= c02;
1715               for (int d = 0; d < dim; ++d) dg2_c2 += PetscSqr(refSpaceDer[d]);
1716               dg2_c2 *= (double)c02;
1717               if (dg2_c2 >= .999) {
1718                 xx[gidx] = vj[qj * dim + 0]; /* coordinate */
1719                 yy[gidx] = vj[qj * dim + 1];
1720                 if (dim == 3) zz[gidx] = vj[qj * dim + 2];
1721                 PetscCall(PetscPrintf(ctx->comm, "Error: %12.5e %" PetscInt_FMT ".%" PetscInt_FMT ") dg2/c02 = %12.5e x= %12.5e %12.5e %12.5e\n", (double)PetscSqrtReal(xx[gidx] * xx[gidx] + yy[gidx] * yy[gidx] + zz[gidx] * zz[gidx]), ej, qj, dg2_c2, (double)xx[gidx], (double)yy[gidx], (double)zz[gidx]));
1722               } else {
1723                 PetscReal fact = c02 / PetscSqrtReal(1. - dg2_c2);
1724                 for (int d = 0; d < dim; ++d) refSpaceDer[d] *= fact;
1725                 // could test with other point u' that (grad - grad') * U (refSpaceDer, refSpaceDer') == 0
1726               }
1727             }
1728             if (xx[gidx] == 1e10) {
1729               for (int d = 0; d < dim; ++d) {
1730                 for (int e = 0; e < dim; ++e) eGradPhi[d] += invJ[e * dim + d] * refSpaceDer[e];
1731               }
1732               xx[gidx] = eGradPhi[0];
1733               yy[gidx] = eGradPhi[1];
1734               if (dim == 3) zz[gidx] = eGradPhi[2];
1735             }
1736           } else {
1737             xx[gidx] = vj[qj * dim + 0]; /* coordinate */
1738             yy[gidx] = vj[qj * dim + 1];
1739             if (dim == 3) zz[gidx] = vj[qj * dim + 2];
1740           }
1741         } /* q */
1742         if (ctx->use_energy_tensor_trick) PetscCall(DMPlexVecRestoreClosure(dmEnergy, e_section, v2_2, ej + cStart, NULL, &coefs));
1743       } /* ej */
1744       if (ctx->use_energy_tensor_trick) {
1745         PetscCall(DMRestoreLocalVector(dmEnergy, &v2_2));
1746         PetscCall(DMDestroy(&dmEnergy));
1747       }
1748     } /* grid */
1749     if (ctx->use_energy_tensor_trick) PetscCall(PetscFEDestroy(&fe));
1750     /* cache static data */
1751     if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) {
1752 #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_KOKKOS_KERNELS)
1753       if (ctx->deviceType == LANDAU_CUDA) {
1754   #if defined(PETSC_HAVE_CUDA)
1755         PetscCall(LandauCUDAStaticDataSet(ctx->plex[0], Nq, ctx->batch_sz, ctx->num_grids, numCells, ctx->species_offset, ctx->mat_offset, nu_alpha, nu_beta, invMass, (PetscReal *)ctx->lambdas, invJ_a, xx, yy, zz, ww, &ctx->SData_d));
1756   #else
1757         SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type cuda not built");
1758   #endif
1759       } else if (ctx->deviceType == LANDAU_KOKKOS) {
1760   #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1761         PetscCall(LandauKokkosStaticDataSet(ctx->plex[0], Nq, ctx->batch_sz, ctx->num_grids, numCells, ctx->species_offset, ctx->mat_offset, nu_alpha, nu_beta, invMass, (PetscReal *)ctx->lambdas, invJ_a, xx, yy, zz, ww, &ctx->SData_d));
1762   #else
1763         SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type kokkos not built");
1764   #endif
1765       }
1766 #endif
1767       /* free */
1768       PetscCall(PetscFree4(ww, xx, yy, invJ_a));
1769       if (dim == 3) PetscCall(PetscFree(zz));
1770     } else {                                                                                                                                                                   /* CPU version, just copy in, only use part */
1771       PetscReal *nu_alpha_p = (PetscReal *)ctx->SData_d.alpha, *nu_beta_p = (PetscReal *)ctx->SData_d.beta, *invMass_p = (PetscReal *)ctx->SData_d.invMass, *lambdas_p = NULL; // why set these ?
1772       ctx->SData_d.w    = (void *)ww;
1773       ctx->SData_d.x    = (void *)xx;
1774       ctx->SData_d.y    = (void *)yy;
1775       ctx->SData_d.z    = (void *)zz;
1776       ctx->SData_d.invJ = (void *)invJ_a;
1777       PetscCall(PetscMalloc4(ctx->num_species, &nu_alpha_p, ctx->num_species, &nu_beta_p, ctx->num_species, &invMass_p, LANDAU_MAX_GRIDS * LANDAU_MAX_GRIDS, &lambdas_p));
1778       for (PetscInt ii = 0; ii < ctx->num_species; ii++) {
1779         nu_alpha_p[ii] = nu_alpha[ii];
1780         nu_beta_p[ii]  = nu_beta[ii];
1781         invMass_p[ii]  = invMass[ii];
1782       }
1783       ctx->SData_d.alpha   = (void *)nu_alpha_p;
1784       ctx->SData_d.beta    = (void *)nu_beta_p;
1785       ctx->SData_d.invMass = (void *)invMass_p;
1786       ctx->SData_d.lambdas = (void *)lambdas_p;
1787       for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) {
1788         PetscReal(*lambdas)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS] = (PetscReal(*)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS])ctx->SData_d.lambdas;
1789         for (PetscInt gridj = 0; gridj < LANDAU_MAX_GRIDS; gridj++) { (*lambdas)[grid][gridj] = ctx->lambdas[grid][gridj]; }
1790       }
1791     }
1792     PetscCall(PetscLogEventEnd(ctx->events[7], 0, 0, 0, 0));
1793   } // initialize
1794   PetscFunctionReturn(PETSC_SUCCESS);
1795 }
1796 
1797 /* < v, u > */
1798 static void g0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1799 {
1800   g0[0] = 1.;
1801 }
1802 
1803 /* < v, u > */
1804 static void g0_fake(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1805 {
1806   static double ttt = 1e-12;
1807   g0[0]             = ttt++;
1808 }
1809 
1810 /* < v, u > */
1811 static void g0_r(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1812 {
1813   g0[0] = 2. * PETSC_PI * x[0];
1814 }
1815 
1816 static PetscErrorCode MatrixNfDestroy(void *ptr)
1817 {
1818   PetscInt *nf = (PetscInt *)ptr;
1819   PetscFunctionBegin;
1820   PetscCall(PetscFree(nf));
1821   PetscFunctionReturn(PETSC_SUCCESS);
1822 }
1823 
1824 /*
1825  LandauCreateJacobianMatrix - creates ctx->J with without real data. Hard to keep sparse.
1826   - Like DMPlexLandauCreateMassMatrix. Should remove one and combine
1827   - has old support for field major ordering
1828  */
1829 static PetscErrorCode LandauCreateJacobianMatrix(MPI_Comm comm, Vec X, IS grid_batch_is_inv[LANDAU_MAX_GRIDS], LandauCtx *ctx)
1830 {
1831   PetscInt *idxs = NULL;
1832   Mat       subM[LANDAU_MAX_GRIDS];
1833 
1834   PetscFunctionBegin;
1835   if (!ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
1836     PetscFunctionReturn(PETSC_SUCCESS);
1837   }
1838   // get the RCM for this grid to separate out species into blocks -- create 'idxs' & 'ctx->batch_is' -- not used
1839   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) PetscCall(PetscMalloc1(ctx->mat_offset[ctx->num_grids] * ctx->batch_sz, &idxs));
1840   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1841     const PetscInt *values, n = ctx->mat_offset[grid + 1] - ctx->mat_offset[grid];
1842     Mat             gMat;
1843     DM              massDM;
1844     PetscDS         prob;
1845     Vec             tvec;
1846     // get "mass" matrix for reordering
1847     PetscCall(DMClone(ctx->plex[grid], &massDM));
1848     PetscCall(DMCopyFields(ctx->plex[grid], massDM));
1849     PetscCall(DMCreateDS(massDM));
1850     PetscCall(DMGetDS(massDM, &prob));
1851     for (int ix = 0, ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++, ix++) PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_fake, NULL, NULL, NULL));
1852     PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only")); // this trick is need to both sparsify the matrix and avoid runtime error
1853     PetscCall(DMCreateMatrix(massDM, &gMat));
1854     PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false"));
1855     PetscCall(MatSetOption(gMat, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE));
1856     PetscCall(MatSetOption(gMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
1857     PetscCall(DMCreateLocalVector(ctx->plex[grid], &tvec));
1858     PetscCall(DMPlexSNESComputeJacobianFEM(massDM, tvec, gMat, gMat, ctx));
1859     PetscCall(MatViewFromOptions(gMat, NULL, "-dm_landau_reorder_mat_view"));
1860     PetscCall(DMDestroy(&massDM));
1861     PetscCall(VecDestroy(&tvec));
1862     subM[grid] = gMat;
1863     if (ctx->gpu_assembly && ctx->jacobian_field_major_order) {
1864       MatOrderingType rtype = MATORDERINGRCM;
1865       IS              isrow, isicol;
1866       PetscCall(MatGetOrdering(gMat, rtype, &isrow, &isicol));
1867       PetscCall(ISInvertPermutation(isrow, PETSC_DECIDE, &grid_batch_is_inv[grid]));
1868       PetscCall(ISGetIndices(isrow, &values));
1869       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid
1870 #if !defined(LANDAU_SPECIES_MAJOR)
1871         PetscInt N = ctx->mat_offset[ctx->num_grids], n0 = ctx->mat_offset[grid] + b_id * N;
1872         for (int ii = 0; ii < n; ++ii) idxs[n0 + ii] = values[ii] + n0;
1873 #else
1874         PetscInt n0 = ctx->mat_offset[grid] * ctx->batch_sz + b_id * n;
1875         for (int ii = 0; ii < n; ++ii) idxs[n0 + ii] = values[ii] + n0;
1876 #endif
1877       }
1878       PetscCall(ISRestoreIndices(isrow, &values));
1879       PetscCall(ISDestroy(&isrow));
1880       PetscCall(ISDestroy(&isicol));
1881     }
1882   }
1883   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) PetscCall(ISCreateGeneral(comm, ctx->mat_offset[ctx->num_grids] * ctx->batch_sz, idxs, PETSC_OWN_POINTER, &ctx->batch_is));
1884   // get a block matrix
1885   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1886     Mat      B = subM[grid];
1887     PetscInt nloc, nzl, *colbuf, row, COL_BF_SIZE = 1024;
1888     PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
1889     PetscCall(MatGetSize(B, &nloc, NULL));
1890     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
1891       const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
1892       const PetscInt    *cols;
1893       const PetscScalar *vals;
1894       for (int i = 0; i < nloc; i++) {
1895         PetscCall(MatGetRow(B, i, &nzl, NULL, NULL));
1896         if (nzl > COL_BF_SIZE) {
1897           PetscCall(PetscFree(colbuf));
1898           PetscCall(PetscInfo(ctx->plex[grid], "Realloc buffer %" PetscInt_FMT " to %" PetscInt_FMT " (row size %" PetscInt_FMT ") \n", COL_BF_SIZE, 2 * COL_BF_SIZE, nzl));
1899           COL_BF_SIZE = nzl;
1900           PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
1901         }
1902         PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
1903         for (int j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset;
1904         row = i + moffset;
1905         PetscCall(MatSetValues(ctx->J, 1, &row, nzl, colbuf, vals, INSERT_VALUES));
1906         PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
1907       }
1908     }
1909     PetscCall(PetscFree(colbuf));
1910   }
1911   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(MatDestroy(&subM[grid]));
1912   PetscCall(MatAssemblyBegin(ctx->J, MAT_FINAL_ASSEMBLY));
1913   PetscCall(MatAssemblyEnd(ctx->J, MAT_FINAL_ASSEMBLY));
1914 
1915   // debug
1916   PetscCall(MatViewFromOptions(ctx->J, NULL, "-dm_landau_mat_view"));
1917   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) {
1918     Mat mat_block_order;
1919     PetscCall(MatCreateSubMatrix(ctx->J, ctx->batch_is, ctx->batch_is, MAT_INITIAL_MATRIX, &mat_block_order)); // use MatPermute
1920     PetscCall(MatViewFromOptions(mat_block_order, NULL, "-dm_landau_mat_view"));
1921     PetscCall(MatDestroy(&mat_block_order));
1922     PetscCall(VecScatterCreate(X, ctx->batch_is, X, NULL, &ctx->plex_batch));
1923     PetscCall(VecDuplicate(X, &ctx->work_vec));
1924   }
1925 
1926   PetscFunctionReturn(PETSC_SUCCESS);
1927 }
1928 
1929 PetscErrorCode DMPlexLandauCreateMassMatrix(DM pack, Mat *Amat);
1930 /*@C
1931  DMPlexLandauCreateVelocitySpace - Create a DMPlex velocity space mesh
1932 
1933  Collective
1934 
1935  Input Parameters:
1936  +   comm  - The MPI communicator
1937  .   dim - velocity space dimension (2 for axisymmetric, 3 for full 3X + 3V solver)
1938  -   prefix - prefix for options (not tested)
1939 
1940  Output Parameter:
1941  .   pack  - The DM object representing the mesh
1942  +   X - A vector (user destroys)
1943  -   J - Optional matrix (object destroys)
1944 
1945  Level: beginner
1946 
1947  .keywords: mesh
1948 .seealso: `DMPlexCreate()`, `DMPlexLandauDestroyVelocitySpace()`
1949  @*/
1950 PetscErrorCode DMPlexLandauCreateVelocitySpace(MPI_Comm comm, PetscInt dim, const char prefix[], Vec *X, Mat *J, DM *pack)
1951 {
1952   LandauCtx *ctx;
1953   Vec        Xsub[LANDAU_MAX_GRIDS];
1954   IS         grid_batch_is_inv[LANDAU_MAX_GRIDS];
1955 
1956   PetscFunctionBegin;
1957   PetscCheck(dim == 2 || dim == 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Only 2D and 3D supported");
1958   PetscCheck(LANDAU_DIM == dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM);
1959   PetscCall(PetscNew(&ctx));
1960   ctx->comm = comm; /* used for diagnostics and global errors */
1961   /* process options */
1962   PetscCall(ProcessOptions(ctx, prefix));
1963   if (dim == 2) ctx->use_relativistic_corrections = PETSC_FALSE;
1964   /* Create Mesh */
1965   PetscCall(DMCompositeCreate(PETSC_COMM_SELF, pack));
1966   PetscCall(PetscLogEventBegin(ctx->events[13], 0, 0, 0, 0));
1967   PetscCall(PetscLogEventBegin(ctx->events[15], 0, 0, 0, 0));
1968   PetscCall(LandauDMCreateVMeshes(PETSC_COMM_SELF, dim, prefix, ctx, *pack)); // creates grids (Forest of AMR)
1969   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1970     /* create FEM */
1971     PetscCall(SetupDS(ctx->plex[grid], dim, grid, ctx));
1972     /* set initial state */
1973     PetscCall(DMCreateGlobalVector(ctx->plex[grid], &Xsub[grid]));
1974     PetscCall(PetscObjectSetName((PetscObject)Xsub[grid], "u_orig"));
1975     /* initial static refinement, no solve */
1976     PetscCall(LandauSetInitialCondition(ctx->plex[grid], Xsub[grid], grid, 0, 1, ctx));
1977     /* forest refinement - forest goes in (if forest), plex comes out */
1978     if (ctx->use_p4est) {
1979       DM plex;
1980       PetscCall(adapt(grid, ctx, &Xsub[grid]));                                      // forest goes in, plex comes out
1981       PetscCall(DMViewFromOptions(ctx->plex[grid], NULL, "-dm_landau_amr_dm_view")); // need to differentiate - todo
1982       PetscCall(VecViewFromOptions(Xsub[grid], NULL, "-dm_landau_amr_vec_view"));
1983       // convert to plex, all done with this level
1984       PetscCall(DMConvert(ctx->plex[grid], DMPLEX, &plex));
1985       PetscCall(DMDestroy(&ctx->plex[grid]));
1986       ctx->plex[grid] = plex;
1987     }
1988 #if !defined(LANDAU_SPECIES_MAJOR)
1989     PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid]));
1990 #else
1991     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid
1992       PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid]));
1993     }
1994 #endif
1995     PetscCall(DMSetApplicationContext(ctx->plex[grid], ctx));
1996   }
1997 #if !defined(LANDAU_SPECIES_MAJOR)
1998   // stack the batched DMs, could do it all here!!! b_id=0
1999   for (PetscInt b_id = 1; b_id < ctx->batch_sz; b_id++) {
2000     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid]));
2001   }
2002 #endif
2003   // create ctx->mat_offset
2004   ctx->mat_offset[0] = 0;
2005   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2006     PetscInt n;
2007     PetscCall(VecGetLocalSize(Xsub[grid], &n));
2008     ctx->mat_offset[grid + 1] = ctx->mat_offset[grid] + n;
2009   }
2010   // creat DM & Jac
2011   PetscCall(DMSetApplicationContext(*pack, ctx));
2012   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only"));
2013   PetscCall(DMCreateMatrix(*pack, &ctx->J));
2014   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false"));
2015   PetscCall(MatSetOption(ctx->J, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE));
2016   PetscCall(MatSetOption(ctx->J, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
2017   PetscCall(PetscObjectSetName((PetscObject)ctx->J, "Jac"));
2018   // construct initial conditions in X
2019   PetscCall(DMCreateGlobalVector(*pack, X));
2020   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2021     PetscInt n;
2022     PetscCall(VecGetLocalSize(Xsub[grid], &n));
2023     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
2024       PetscScalar const *values;
2025       const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
2026       PetscCall(LandauSetInitialCondition(ctx->plex[grid], Xsub[grid], grid, b_id, ctx->batch_sz, ctx));
2027       PetscCall(VecGetArrayRead(Xsub[grid], &values)); // Drop whole grid in Plex ordering
2028       for (int i = 0, idx = moffset; i < n; i++, idx++) PetscCall(VecSetValue(*X, idx, values[i], INSERT_VALUES));
2029       PetscCall(VecRestoreArrayRead(Xsub[grid], &values));
2030     }
2031   }
2032   // cleanup
2033   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(VecDestroy(&Xsub[grid]));
2034   /* check for correct matrix type */
2035   if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
2036     PetscBool flg;
2037     if (ctx->deviceType == LANDAU_CUDA) {
2038       PetscCall(PetscObjectTypeCompareAny((PetscObject)ctx->J, &flg, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, MATAIJCUSPARSE, ""));
2039       PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must use '-dm_mat_type aijcusparse -dm_vec_type cuda' for GPU assembly and Cuda or use '-dm_landau_device_type cpu'");
2040     } else if (ctx->deviceType == LANDAU_KOKKOS) {
2041       PetscCall(PetscObjectTypeCompareAny((PetscObject)ctx->J, &flg, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, MATAIJKOKKOS, ""));
2042 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
2043       PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for GPU assembly and Kokkos or use '-dm_landau_device_type cpu'");
2044 #else
2045       PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must configure with '--download-kokkos-kernels' for GPU assembly and Kokkos or use '-dm_landau_device_type cpu'");
2046 #endif
2047     }
2048   }
2049   PetscCall(PetscLogEventEnd(ctx->events[15], 0, 0, 0, 0));
2050 
2051   // create field major ordering
2052   ctx->work_vec   = NULL;
2053   ctx->plex_batch = NULL;
2054   ctx->batch_is   = NULL;
2055   for (int i = 0; i < LANDAU_MAX_GRIDS; i++) grid_batch_is_inv[i] = NULL;
2056   PetscCall(PetscLogEventBegin(ctx->events[12], 0, 0, 0, 0));
2057   PetscCall(LandauCreateJacobianMatrix(comm, *X, grid_batch_is_inv, ctx));
2058   PetscCall(PetscLogEventEnd(ctx->events[12], 0, 0, 0, 0));
2059 
2060   // create AMR GPU assembly maps and static GPU data
2061   PetscCall(CreateStaticData(dim, grid_batch_is_inv, ctx));
2062 
2063   PetscCall(PetscLogEventEnd(ctx->events[13], 0, 0, 0, 0));
2064 
2065   // create mass matrix
2066   PetscCall(DMPlexLandauCreateMassMatrix(*pack, NULL));
2067 
2068   if (J) *J = ctx->J;
2069 
2070   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) {
2071     PetscContainer container;
2072     // cache ctx for KSP with batch/field major Jacobian ordering -ksp_type gmres/etc -dm_landau_jacobian_field_major_order
2073     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
2074     PetscCall(PetscContainerSetPointer(container, (void *)ctx));
2075     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "LandauCtx", (PetscObject)container));
2076     PetscCall(PetscContainerDestroy(&container));
2077     // batch solvers need to map -- can batch solvers work
2078     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
2079     PetscCall(PetscContainerSetPointer(container, (void *)ctx->plex_batch));
2080     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "plex_batch_is", (PetscObject)container));
2081     PetscCall(PetscContainerDestroy(&container));
2082   }
2083   // for batch solvers
2084   {
2085     PetscContainer container;
2086     PetscInt      *pNf;
2087     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
2088     PetscCall(PetscMalloc1(sizeof(*pNf), &pNf));
2089     *pNf = ctx->batch_sz;
2090     PetscCall(PetscContainerSetPointer(container, (void *)pNf));
2091     PetscCall(PetscContainerSetUserDestroy(container, MatrixNfDestroy));
2092     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "batch size", (PetscObject)container));
2093     PetscCall(PetscContainerDestroy(&container));
2094   }
2095 
2096   PetscFunctionReturn(PETSC_SUCCESS);
2097 }
2098 
2099 /*@
2100  DMPlexLandauAccess - Access to the distribution function with user callback
2101 
2102  Collective
2103 
2104  Input Parameters:
2105  .   pack - the DMComposite
2106  +   func - call back function
2107  .   user_ctx - user context
2108 
2109  Input/Output Parameter:
2110  .   X - Vector to data to
2111 
2112  Level: advanced
2113 
2114  .keywords: mesh
2115 .seealso: `DMPlexLandauCreateVelocitySpace()`
2116  @*/
2117 PetscErrorCode DMPlexLandauAccess(DM pack, Vec X, PetscErrorCode (*func)(DM, Vec, PetscInt, PetscInt, PetscInt, void *), void *user_ctx)
2118 {
2119   LandauCtx *ctx;
2120   PetscFunctionBegin;
2121   PetscCall(DMGetApplicationContext(pack, &ctx)); // uses ctx->num_grids; ctx->plex[grid]; ctx->batch_sz; ctx->mat_offset
2122   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2123     PetscInt dim, n;
2124     PetscCall(DMGetDimension(pack, &dim));
2125     for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) {
2126       Vec      vec;
2127       PetscInt vf[1] = {i0};
2128       IS       vis;
2129       DM       vdm;
2130       PetscCall(DMCreateSubDM(ctx->plex[grid], 1, vf, &vis, &vdm));
2131       PetscCall(DMSetApplicationContext(vdm, ctx)); // the user might want this
2132       PetscCall(DMCreateGlobalVector(vdm, &vec));
2133       PetscCall(VecGetSize(vec, &n));
2134       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
2135         const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
2136         PetscCall(VecZeroEntries(vec));
2137         /* Add your data with 'dm' for species 'sp' to 'vec' */
2138         PetscCall(func(vdm, vec, i0, grid, b_id, user_ctx));
2139         /* add to global */
2140         PetscScalar const *values;
2141         const PetscInt    *offsets;
2142         PetscCall(VecGetArrayRead(vec, &values));
2143         PetscCall(ISGetIndices(vis, &offsets));
2144         for (int i = 0; i < n; i++) PetscCall(VecSetValue(X, moffset + offsets[i], values[i], ADD_VALUES));
2145         PetscCall(VecRestoreArrayRead(vec, &values));
2146         PetscCall(ISRestoreIndices(vis, &offsets));
2147       } // batch
2148       PetscCall(VecDestroy(&vec));
2149       PetscCall(ISDestroy(&vis));
2150       PetscCall(DMDestroy(&vdm));
2151     }
2152   } // grid
2153   PetscFunctionReturn(PETSC_SUCCESS);
2154 }
2155 
2156 /*@
2157  DMPlexLandauDestroyVelocitySpace - Destroy a DMPlex velocity space mesh
2158 
2159  Collective
2160 
2161  Input/Output Parameters:
2162  .   dm - the dm to destroy
2163 
2164  Level: beginner
2165 
2166  .keywords: mesh
2167 .seealso: `DMPlexLandauCreateVelocitySpace()`
2168  @*/
2169 PetscErrorCode DMPlexLandauDestroyVelocitySpace(DM *dm)
2170 {
2171   LandauCtx *ctx;
2172   PetscFunctionBegin;
2173   PetscCall(DMGetApplicationContext(*dm, &ctx));
2174   PetscCall(MatDestroy(&ctx->M));
2175   PetscCall(MatDestroy(&ctx->J));
2176   for (PetscInt ii = 0; ii < ctx->num_species; ii++) PetscCall(PetscFEDestroy(&ctx->fe[ii]));
2177   PetscCall(ISDestroy(&ctx->batch_is));
2178   PetscCall(VecDestroy(&ctx->work_vec));
2179   PetscCall(VecScatterDestroy(&ctx->plex_batch));
2180   if (ctx->deviceType == LANDAU_CUDA) {
2181 #if defined(PETSC_HAVE_CUDA)
2182     PetscCall(LandauCUDAStaticDataClear(&ctx->SData_d));
2183 #else
2184     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "cuda");
2185 #endif
2186   } else if (ctx->deviceType == LANDAU_KOKKOS) {
2187 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
2188     PetscCall(LandauKokkosStaticDataClear(&ctx->SData_d));
2189 #else
2190     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "kokkos");
2191 #endif
2192   } else {
2193     if (ctx->SData_d.x) { /* in a CPU run */
2194       PetscReal *invJ = (PetscReal *)ctx->SData_d.invJ, *xx = (PetscReal *)ctx->SData_d.x, *yy = (PetscReal *)ctx->SData_d.y, *zz = (PetscReal *)ctx->SData_d.z, *ww = (PetscReal *)ctx->SData_d.w;
2195       LandauIdx *coo_elem_offsets = (LandauIdx *)ctx->SData_d.coo_elem_offsets, *coo_elem_fullNb = (LandauIdx *)ctx->SData_d.coo_elem_fullNb, (*coo_elem_point_offsets)[LANDAU_MAX_NQ + 1] = (LandauIdx(*)[LANDAU_MAX_NQ + 1]) ctx->SData_d.coo_elem_point_offsets;
2196       PetscCall(PetscFree4(ww, xx, yy, invJ));
2197       if (zz) PetscCall(PetscFree(zz));
2198       if (coo_elem_offsets) {
2199         PetscCall(PetscFree3(coo_elem_offsets, coo_elem_fullNb, coo_elem_point_offsets)); // could be NULL
2200       }
2201       PetscCall(PetscFree4(ctx->SData_d.alpha, ctx->SData_d.beta, ctx->SData_d.invMass, ctx->SData_d.lambdas));
2202     }
2203   }
2204 
2205   if (ctx->times[LANDAU_MATRIX_TOTAL] > 0) { // OMP timings
2206     PetscCall(PetscPrintf(ctx->comm, "TSStep               N  1.0 %10.3e\n", ctx->times[LANDAU_EX2_TSSOLVE]));
2207     PetscCall(PetscPrintf(ctx->comm, "2:           Solve:  %10.3e with %" PetscInt_FMT " threads\n", ctx->times[LANDAU_EX2_TSSOLVE] - ctx->times[LANDAU_MATRIX_TOTAL], ctx->batch_sz));
2208     PetscCall(PetscPrintf(ctx->comm, "3:          Landau:  %10.3e\n", ctx->times[LANDAU_MATRIX_TOTAL]));
2209     PetscCall(PetscPrintf(ctx->comm, "Landau Jacobian       %" PetscInt_FMT " 1.0 %10.3e\n", (PetscInt)ctx->times[LANDAU_JACOBIAN_COUNT], ctx->times[LANDAU_JACOBIAN]));
2210     PetscCall(PetscPrintf(ctx->comm, "Landau Operator       N 1.0  %10.3e\n", ctx->times[LANDAU_OPERATOR]));
2211     PetscCall(PetscPrintf(ctx->comm, "Landau Mass           N 1.0  %10.3e\n", ctx->times[LANDAU_MASS]));
2212     PetscCall(PetscPrintf(ctx->comm, " Jac-f-df (GPU)       N 1.0  %10.3e\n", ctx->times[LANDAU_F_DF]));
2213     PetscCall(PetscPrintf(ctx->comm, " Kernel (GPU)         N 1.0  %10.3e\n", ctx->times[LANDAU_KERNEL]));
2214     PetscCall(PetscPrintf(ctx->comm, "MatLUFactorNum        X 1.0 %10.3e\n", ctx->times[KSP_FACTOR]));
2215     PetscCall(PetscPrintf(ctx->comm, "MatSolve              X 1.0 %10.3e\n", ctx->times[KSP_SOLVE]));
2216   }
2217   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMDestroy(&ctx->plex[grid]));
2218   PetscCall(PetscFree(ctx));
2219   PetscCall(DMDestroy(dm));
2220   PetscFunctionReturn(PETSC_SUCCESS);
2221 }
2222 
2223 /* < v, ru > */
2224 static void f0_s_den(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2225 {
2226   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2227   f0[0]       = u[ii];
2228 }
2229 
2230 /* < v, ru > */
2231 static void f0_s_mom(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2232 {
2233   PetscInt ii = (PetscInt)PetscRealPart(constants[0]), jj = (PetscInt)PetscRealPart(constants[1]);
2234   f0[0] = x[jj] * u[ii]; /* x momentum */
2235 }
2236 
2237 static void f0_s_v2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2238 {
2239   PetscInt i, ii = (PetscInt)PetscRealPart(constants[0]);
2240   double   tmp1 = 0.;
2241   for (i = 0; i < dim; ++i) tmp1 += x[i] * x[i];
2242   f0[0] = tmp1 * u[ii];
2243 }
2244 
2245 static PetscErrorCode gamma_n_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *actx)
2246 {
2247   const PetscReal *c2_0_arr = ((PetscReal *)actx);
2248   const PetscReal  c02      = c2_0_arr[0];
2249 
2250   PetscFunctionBegin;
2251   for (int s = 0; s < Nf; s++) {
2252     PetscReal tmp1 = 0.;
2253     for (int i = 0; i < dim; ++i) tmp1 += x[i] * x[i];
2254 #if defined(PETSC_USE_DEBUG)
2255     u[s] = PetscSqrtReal(1. + tmp1 / c02); //  u[0] = PetscSqrtReal(1. + xx);
2256 #else
2257     {
2258       PetscReal xx = tmp1 / c02;
2259       u[s]         = xx / (PetscSqrtReal(1. + xx) + 1.); // better conditioned = xx/(PetscSqrtReal(1. + xx) + 1.)
2260     }
2261 #endif
2262   }
2263   PetscFunctionReturn(PETSC_SUCCESS);
2264 }
2265 
2266 /* < v, ru > */
2267 static void f0_s_rden(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2268 {
2269   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2270   f0[0]       = 2. * PETSC_PI * x[0] * u[ii];
2271 }
2272 
2273 /* < v, ru > */
2274 static void f0_s_rmom(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2275 {
2276   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2277   f0[0]       = 2. * PETSC_PI * x[0] * x[1] * u[ii];
2278 }
2279 
2280 static void f0_s_rv2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2281 {
2282   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2283   f0[0]       = 2. * PETSC_PI * x[0] * (x[0] * x[0] + x[1] * x[1]) * u[ii];
2284 }
2285 
2286 /*@
2287  DMPlexLandauPrintNorms - collects moments and prints them
2288 
2289  Collective
2290 
2291  Input Parameters:
2292  +   X  - the state
2293  -   stepi - current step to print
2294 
2295  Level: beginner
2296 
2297  .keywords: mesh
2298 .seealso: `DMPlexLandauCreateVelocitySpace()`
2299  @*/
2300 PetscErrorCode DMPlexLandauPrintNorms(Vec X, PetscInt stepi)
2301 {
2302   LandauCtx  *ctx;
2303   PetscDS     prob;
2304   DM          pack;
2305   PetscInt    cStart, cEnd, dim, ii, i0, nDMs;
2306   PetscScalar xmomentumtot = 0, ymomentumtot = 0, zmomentumtot = 0, energytot = 0, densitytot = 0, tt[LANDAU_MAX_SPECIES];
2307   PetscScalar xmomentum[LANDAU_MAX_SPECIES], ymomentum[LANDAU_MAX_SPECIES], zmomentum[LANDAU_MAX_SPECIES], energy[LANDAU_MAX_SPECIES], density[LANDAU_MAX_SPECIES];
2308   Vec        *globXArray;
2309 
2310   PetscFunctionBegin;
2311   PetscCall(VecGetDM(X, &pack));
2312   PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Vector has no DM");
2313   PetscCall(DMGetDimension(pack, &dim));
2314   PetscCheck(dim == 2 || dim == 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " not in [2,3]", dim);
2315   PetscCall(DMGetApplicationContext(pack, &ctx));
2316   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2317   /* print momentum and energy */
2318   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
2319   PetscCheck(nDMs == ctx->num_grids * ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "#DM wrong %" PetscInt_FMT " %" PetscInt_FMT, nDMs, ctx->num_grids * ctx->batch_sz);
2320   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
2321   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
2322   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2323     Vec Xloc = globXArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)];
2324     PetscCall(DMGetDS(ctx->plex[grid], &prob));
2325     for (ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
2326       PetscScalar user[2] = {(PetscScalar)i0, (PetscScalar)ctx->charges[ii]};
2327       PetscCall(PetscDSSetConstants(prob, 2, user));
2328       if (dim == 2) { /* 2/3X + 3V (cylindrical coordinates) */
2329         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rden));
2330         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2331         density[ii] = tt[0] * ctx->n_0 * ctx->charges[ii];
2332         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rmom));
2333         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2334         zmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2335         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rv2));
2336         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2337         energy[ii] = tt[0] * 0.5 * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ii];
2338         zmomentumtot += zmomentum[ii];
2339         energytot += energy[ii];
2340         densitytot += density[ii];
2341         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") species-%" PetscInt_FMT ": charge density= %20.13e z-momentum= %20.13e energy= %20.13e", stepi, ii, (double)PetscRealPart(density[ii]), (double)PetscRealPart(zmomentum[ii]), (double)PetscRealPart(energy[ii])));
2342       } else { /* 2/3Xloc + 3V */
2343         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_den));
2344         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2345         density[ii] = tt[0] * ctx->n_0 * ctx->charges[ii];
2346         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_mom));
2347         user[1] = 0;
2348         PetscCall(PetscDSSetConstants(prob, 2, user));
2349         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2350         xmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2351         user[1]       = 1;
2352         PetscCall(PetscDSSetConstants(prob, 2, user));
2353         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2354         ymomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2355         user[1]       = 2;
2356         PetscCall(PetscDSSetConstants(prob, 2, user));
2357         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2358         zmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2359         if (ctx->use_relativistic_corrections) {
2360           /* gamma * M * f */
2361           if (ii == 0 && grid == 0) { // do all at once
2362             Vec Mf, globGamma, *globMfArray, *globGammaArray;
2363             PetscErrorCode (*gammaf[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {gamma_n_f};
2364             PetscReal *c2_0[1], data[1];
2365 
2366             PetscCall(VecDuplicate(X, &globGamma));
2367             PetscCall(VecDuplicate(X, &Mf));
2368             PetscCall(PetscMalloc(sizeof(*globMfArray) * nDMs, &globMfArray));
2369             PetscCall(PetscMalloc(sizeof(*globMfArray) * nDMs, &globGammaArray));
2370             /* M * f */
2371             PetscCall(MatMult(ctx->M, X, Mf));
2372             /* gamma */
2373             PetscCall(DMCompositeGetAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2374             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // yes a grid loop in a grid loop to print nice, need to fix for batching
2375               Vec v1  = globGammaArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)];
2376               data[0] = PetscSqr(C_0(ctx->v_0));
2377               c2_0[0] = &data[0];
2378               PetscCall(DMProjectFunction(ctx->plex[grid], 0., gammaf, (void **)c2_0, INSERT_ALL_VALUES, v1));
2379             }
2380             PetscCall(DMCompositeRestoreAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2381             /* gamma * Mf */
2382             PetscCall(DMCompositeGetAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2383             PetscCall(DMCompositeGetAccessArray(pack, Mf, nDMs, NULL, globMfArray));
2384             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // yes a grid loop in a grid loop to print nice
2385               PetscInt Nf    = ctx->species_offset[grid + 1] - ctx->species_offset[grid], N, bs;
2386               Vec      Mfsub = globMfArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)], Gsub = globGammaArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)], v1, v2;
2387               // get each component
2388               PetscCall(VecGetSize(Mfsub, &N));
2389               PetscCall(VecCreate(ctx->comm, &v1));
2390               PetscCall(VecSetSizes(v1, PETSC_DECIDE, N / Nf));
2391               PetscCall(VecCreate(ctx->comm, &v2));
2392               PetscCall(VecSetSizes(v2, PETSC_DECIDE, N / Nf));
2393               PetscCall(VecSetFromOptions(v1)); // ???
2394               PetscCall(VecSetFromOptions(v2));
2395               // get each component
2396               PetscCall(VecGetBlockSize(Gsub, &bs));
2397               PetscCheck(bs == Nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " != num_species %" PetscInt_FMT " in Gsub", bs, Nf);
2398               PetscCall(VecGetBlockSize(Mfsub, &bs));
2399               PetscCheck(bs == Nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " != num_species %" PetscInt_FMT, bs, Nf);
2400               for (int i = 0, ix = ctx->species_offset[grid]; i < Nf; i++, ix++) {
2401                 PetscScalar val;
2402                 PetscCall(VecStrideGather(Gsub, i, v1, INSERT_VALUES)); // this is not right -- TODO
2403                 PetscCall(VecStrideGather(Mfsub, i, v2, INSERT_VALUES));
2404                 PetscCall(VecDot(v1, v2, &val));
2405                 energy[ix] = PetscRealPart(val) * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ix];
2406               }
2407               PetscCall(VecDestroy(&v1));
2408               PetscCall(VecDestroy(&v2));
2409             } /* grids */
2410             PetscCall(DMCompositeRestoreAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2411             PetscCall(DMCompositeRestoreAccessArray(pack, Mf, nDMs, NULL, globMfArray));
2412             PetscCall(PetscFree(globGammaArray));
2413             PetscCall(PetscFree(globMfArray));
2414             PetscCall(VecDestroy(&globGamma));
2415             PetscCall(VecDestroy(&Mf));
2416           }
2417         } else {
2418           PetscCall(PetscDSSetObjective(prob, 0, &f0_s_v2));
2419           PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2420           energy[ii] = 0.5 * tt[0] * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ii];
2421         }
2422         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") species %" PetscInt_FMT ": density=%20.13e, x-momentum=%20.13e, y-momentum=%20.13e, z-momentum=%20.13e, energy=%21.13e", stepi, ii, (double)PetscRealPart(density[ii]), (double)PetscRealPart(xmomentum[ii]), (double)PetscRealPart(ymomentum[ii]), (double)PetscRealPart(zmomentum[ii]), (double)PetscRealPart(energy[ii])));
2423         xmomentumtot += xmomentum[ii];
2424         ymomentumtot += ymomentum[ii];
2425         zmomentumtot += zmomentum[ii];
2426         energytot += energy[ii];
2427         densitytot += density[ii];
2428       }
2429       if (ctx->num_species > 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
2430     }
2431   }
2432   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
2433   PetscCall(PetscFree(globXArray));
2434   /* totals */
2435   PetscCall(DMPlexGetHeightStratum(ctx->plex[0], 0, &cStart, &cEnd));
2436   if (ctx->num_species > 1) {
2437     if (dim == 2) {
2438       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t%3" PetscInt_FMT ") Total: charge density=%21.13e, momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %" PetscInt_FMT " cells on electron grid)", stepi, (double)PetscRealPart(densitytot), (double)PetscRealPart(zmomentumtot), (double)PetscRealPart(energytot),
2439                             (double)(ctx->masses[1] / ctx->masses[0]), cEnd - cStart));
2440     } else {
2441       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t%3" PetscInt_FMT ") Total: charge density=%21.13e, x-momentum=%21.13e, y-momentum=%21.13e, z-momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %" PetscInt_FMT " cells)", stepi, (double)PetscRealPart(densitytot), (double)PetscRealPart(xmomentumtot), (double)PetscRealPart(ymomentumtot), (double)PetscRealPart(zmomentumtot), (double)PetscRealPart(energytot),
2442                             (double)(ctx->masses[1] / ctx->masses[0]), cEnd - cStart));
2443     }
2444   } else PetscCall(PetscPrintf(PETSC_COMM_WORLD, " -- %" PetscInt_FMT " cells", cEnd - cStart));
2445   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
2446   PetscFunctionReturn(PETSC_SUCCESS);
2447 }
2448 
2449 /*@
2450  DMPlexLandauCreateMassMatrix - Create mass matrix for Landau in Plex space (not field major order of Jacobian)
2451   - puts mass matrix into ctx->M
2452 
2453  Collective
2454 
2455  Input Parameter:
2456 . pack     - the DM object. Puts matrix in Landau context M field
2457 
2458  Output Parameter:
2459 . Amat - The mass matrix (optional), mass matrix is added to the DM context
2460 
2461  Level: beginner
2462 
2463  .keywords: mesh
2464 .seealso: `DMPlexLandauCreateVelocitySpace()`
2465  @*/
2466 PetscErrorCode DMPlexLandauCreateMassMatrix(DM pack, Mat *Amat)
2467 {
2468   DM         mass_pack, massDM[LANDAU_MAX_GRIDS];
2469   PetscDS    prob;
2470   PetscInt   ii, dim, N1 = 1, N2;
2471   LandauCtx *ctx;
2472   Mat        packM, subM[LANDAU_MAX_GRIDS];
2473 
2474   PetscFunctionBegin;
2475   PetscValidHeaderSpecific(pack, DM_CLASSID, 1);
2476   if (Amat) PetscValidPointer(Amat, 2);
2477   PetscCall(DMGetApplicationContext(pack, &ctx));
2478   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2479   PetscCall(PetscLogEventBegin(ctx->events[14], 0, 0, 0, 0));
2480   PetscCall(DMGetDimension(pack, &dim));
2481   PetscCall(DMCompositeCreate(PetscObjectComm((PetscObject)pack), &mass_pack));
2482   /* create pack mass matrix */
2483   for (PetscInt grid = 0, ix = 0; grid < ctx->num_grids; grid++) {
2484     PetscCall(DMClone(ctx->plex[grid], &massDM[grid]));
2485     PetscCall(DMCopyFields(ctx->plex[grid], massDM[grid]));
2486     PetscCall(DMCreateDS(massDM[grid]));
2487     PetscCall(DMGetDS(massDM[grid], &prob));
2488     for (ix = 0, ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++, ix++) {
2489       if (dim == 3) PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_1, NULL, NULL, NULL));
2490       else PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_r, NULL, NULL, NULL));
2491     }
2492 #if !defined(LANDAU_SPECIES_MAJOR)
2493     PetscCall(DMCompositeAddDM(mass_pack, massDM[grid]));
2494 #else
2495     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid
2496       PetscCall(DMCompositeAddDM(mass_pack, massDM[grid]));
2497     }
2498 #endif
2499     PetscCall(DMCreateMatrix(massDM[grid], &subM[grid]));
2500   }
2501 #if !defined(LANDAU_SPECIES_MAJOR)
2502   // stack the batched DMs
2503   for (PetscInt b_id = 1; b_id < ctx->batch_sz; b_id++) {
2504     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCompositeAddDM(mass_pack, massDM[grid]));
2505   }
2506 #endif
2507   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only"));
2508   PetscCall(DMCreateMatrix(mass_pack, &packM));
2509   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false"));
2510   PetscCall(MatSetOption(packM, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE));
2511   PetscCall(MatSetOption(packM, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
2512   PetscCall(DMDestroy(&mass_pack));
2513   /* make mass matrix for each block */
2514   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2515     Vec locX;
2516     DM  plex = massDM[grid];
2517     PetscCall(DMGetLocalVector(plex, &locX));
2518     /* Mass matrix is independent of the input, so no need to fill locX */
2519     PetscCall(DMPlexSNESComputeJacobianFEM(plex, locX, subM[grid], subM[grid], ctx));
2520     PetscCall(DMRestoreLocalVector(plex, &locX));
2521     PetscCall(DMDestroy(&massDM[grid]));
2522   }
2523   PetscCall(MatGetSize(ctx->J, &N1, NULL));
2524   PetscCall(MatGetSize(packM, &N2, NULL));
2525   PetscCheck(N1 == N2, PetscObjectComm((PetscObject)pack), PETSC_ERR_PLIB, "Incorrect matrix sizes: |Jacobian| = %" PetscInt_FMT ", |Mass|=%" PetscInt_FMT, N1, N2);
2526   /* assemble block diagonals */
2527   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2528     Mat      B = subM[grid];
2529     PetscInt nloc, nzl, *colbuf, COL_BF_SIZE = 1024, row;
2530     PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
2531     PetscCall(MatGetSize(B, &nloc, NULL));
2532     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
2533       const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
2534       const PetscInt    *cols;
2535       const PetscScalar *vals;
2536       for (int i = 0; i < nloc; i++) {
2537         PetscCall(MatGetRow(B, i, &nzl, NULL, NULL));
2538         if (nzl > COL_BF_SIZE) {
2539           PetscCall(PetscFree(colbuf));
2540           PetscCall(PetscInfo(pack, "Realloc buffer %" PetscInt_FMT " to %" PetscInt_FMT " (row size %" PetscInt_FMT ") \n", COL_BF_SIZE, 2 * COL_BF_SIZE, nzl));
2541           COL_BF_SIZE = nzl;
2542           PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
2543         }
2544         PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
2545         for (int j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset;
2546         row = i + moffset;
2547         PetscCall(MatSetValues(packM, 1, &row, nzl, colbuf, vals, INSERT_VALUES));
2548         PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
2549       }
2550     }
2551     PetscCall(PetscFree(colbuf));
2552   }
2553   // cleanup
2554   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(MatDestroy(&subM[grid]));
2555   PetscCall(MatAssemblyBegin(packM, MAT_FINAL_ASSEMBLY));
2556   PetscCall(MatAssemblyEnd(packM, MAT_FINAL_ASSEMBLY));
2557   PetscCall(PetscObjectSetName((PetscObject)packM, "mass"));
2558   PetscCall(MatViewFromOptions(packM, NULL, "-dm_landau_mass_view"));
2559   ctx->M = packM;
2560   if (Amat) *Amat = packM;
2561   PetscCall(PetscLogEventEnd(ctx->events[14], 0, 0, 0, 0));
2562   PetscFunctionReturn(PETSC_SUCCESS);
2563 }
2564 
2565 /*@
2566  DMPlexLandauIFunction - TS residual calculation, confusingly this computes the Jacobian w/o mass
2567 
2568  Collective
2569 
2570  Input Parameters:
2571 +   TS  - The time stepping context
2572 .   time_dummy - current time (not used)
2573 .   X - Current state
2574 .   X_t - Time derivative of current state
2575 -   actx - Landau context
2576 
2577  Output Parameter:
2578 .   F  - The residual
2579 
2580  Level: beginner
2581 
2582  .keywords: mesh
2583 .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauIJacobian()`
2584  @*/
2585 PetscErrorCode DMPlexLandauIFunction(TS ts, PetscReal time_dummy, Vec X, Vec X_t, Vec F, void *actx)
2586 {
2587   LandauCtx *ctx = (LandauCtx *)actx;
2588   PetscInt   dim;
2589   DM         pack;
2590 #if defined(PETSC_HAVE_THREADSAFETY)
2591   double starttime, endtime;
2592 #endif
2593   PetscObjectState state;
2594 
2595   PetscFunctionBegin;
2596   PetscCall(TSGetDM(ts, &pack));
2597   PetscCall(DMGetApplicationContext(pack, &ctx));
2598   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2599   if (ctx->stage) PetscCall(PetscLogStagePush(ctx->stage));
2600   PetscCall(PetscLogEventBegin(ctx->events[11], 0, 0, 0, 0));
2601   PetscCall(PetscLogEventBegin(ctx->events[0], 0, 0, 0, 0));
2602 #if defined(PETSC_HAVE_THREADSAFETY)
2603   starttime = MPI_Wtime();
2604 #endif
2605   PetscCall(DMGetDimension(pack, &dim));
2606   PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state));
2607   if (state != ctx->norm_state) {
2608     PetscCall(PetscInfo(ts, "Create Landau Jacobian t=%g J.state %" PetscInt64_FMT " --> %" PetscInt64_FMT "\n", (double)time_dummy, ctx->norm_state, state));
2609     PetscCall(MatZeroEntries(ctx->J));
2610     PetscCall(LandauFormJacobian_Internal(X, ctx->J, dim, 0.0, (void *)ctx));
2611     PetscCall(MatViewFromOptions(ctx->J, NULL, "-dm_landau_jacobian_view"));
2612     PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state));
2613     ctx->norm_state = state;
2614   } else {
2615     PetscCall(PetscInfo(ts, "WARNING Skip forming Jacobian, has not changed %" PetscInt64_FMT "\n", state));
2616   }
2617   /* mat vec for op */
2618   PetscCall(MatMult(ctx->J, X, F)); /* C*f */
2619   /* add time term */
2620   if (X_t) PetscCall(MatMultAdd(ctx->M, X_t, F, F));
2621 #if defined(PETSC_HAVE_THREADSAFETY)
2622   if (ctx->stage) {
2623     endtime = MPI_Wtime();
2624     ctx->times[LANDAU_OPERATOR] += (endtime - starttime);
2625     ctx->times[LANDAU_JACOBIAN] += (endtime - starttime);
2626     ctx->times[LANDAU_MATRIX_TOTAL] += (endtime - starttime);
2627     ctx->times[LANDAU_JACOBIAN_COUNT] += 1;
2628   }
2629 #endif
2630   PetscCall(PetscLogEventEnd(ctx->events[0], 0, 0, 0, 0));
2631   PetscCall(PetscLogEventEnd(ctx->events[11], 0, 0, 0, 0));
2632   if (ctx->stage) PetscCall(PetscLogStagePop());
2633   PetscFunctionReturn(PETSC_SUCCESS);
2634 }
2635 
2636 /*@
2637  DMPlexLandauIJacobian - TS Jacobian construction, confusingly this adds mass
2638 
2639  Collective
2640 
2641  Input Parameters:
2642 +   TS  - The time stepping context
2643 .   time_dummy - current time (not used)
2644 .   X - Current state
2645 .   U_tdummy - Time derivative of current state (not used)
2646 .   shift - shift for du/dt term
2647 -   actx - Landau context
2648 
2649  Output Parameters:
2650 +   Amat  - Jacobian
2651 -   Pmat  - same as Amat
2652 
2653  Level: beginner
2654 
2655  .keywords: mesh
2656 .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauIFunction()`
2657  @*/
2658 PetscErrorCode DMPlexLandauIJacobian(TS ts, PetscReal time_dummy, Vec X, Vec U_tdummy, PetscReal shift, Mat Amat, Mat Pmat, void *actx)
2659 {
2660   LandauCtx *ctx = NULL;
2661   PetscInt   dim;
2662   DM         pack;
2663 #if defined(PETSC_HAVE_THREADSAFETY)
2664   double starttime, endtime;
2665 #endif
2666   PetscObjectState state;
2667 
2668   PetscFunctionBegin;
2669   PetscCall(TSGetDM(ts, &pack));
2670   PetscCall(DMGetApplicationContext(pack, &ctx));
2671   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2672   PetscCheck(Amat == Pmat && Amat == ctx->J, ctx->comm, PETSC_ERR_PLIB, "Amat!=Pmat || Amat!=ctx->J");
2673   PetscCall(DMGetDimension(pack, &dim));
2674   /* get collision Jacobian into A */
2675   if (ctx->stage) PetscCall(PetscLogStagePush(ctx->stage));
2676   PetscCall(PetscLogEventBegin(ctx->events[11], 0, 0, 0, 0));
2677   PetscCall(PetscLogEventBegin(ctx->events[9], 0, 0, 0, 0));
2678 #if defined(PETSC_HAVE_THREADSAFETY)
2679   starttime = MPI_Wtime();
2680 #endif
2681   PetscCall(PetscInfo(ts, "Adding mass to Jacobian t=%g, shift=%g\n", (double)time_dummy, (double)shift));
2682   PetscCheck(shift != 0.0, ctx->comm, PETSC_ERR_PLIB, "zero shift");
2683   PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state));
2684   PetscCheck(state == ctx->norm_state, ctx->comm, PETSC_ERR_PLIB, "wrong state, %" PetscInt64_FMT " %" PetscInt64_FMT "", ctx->norm_state, state);
2685   if (!ctx->use_matrix_mass) {
2686     PetscCall(LandauFormJacobian_Internal(X, ctx->J, dim, shift, (void *)ctx));
2687   } else { /* add mass */
2688     PetscCall(MatAXPY(Pmat, shift, ctx->M, SAME_NONZERO_PATTERN));
2689   }
2690 #if defined(PETSC_HAVE_THREADSAFETY)
2691   if (ctx->stage) {
2692     endtime = MPI_Wtime();
2693     ctx->times[LANDAU_OPERATOR] += (endtime - starttime);
2694     ctx->times[LANDAU_MASS] += (endtime - starttime);
2695     ctx->times[LANDAU_MATRIX_TOTAL] += (endtime - starttime);
2696   }
2697 #endif
2698   PetscCall(PetscLogEventEnd(ctx->events[9], 0, 0, 0, 0));
2699   PetscCall(PetscLogEventEnd(ctx->events[11], 0, 0, 0, 0));
2700   if (ctx->stage) PetscCall(PetscLogStagePop());
2701   PetscFunctionReturn(PETSC_SUCCESS);
2702 }
2703