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