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