xref: /petsc/src/ts/tests/ex30.c (revision bb4b53ef092968f72b740b90dbab8a2b6700db0d)
1 static char help[] = "Grid based Landau collision operator with PIC interface with OpenMP setup. (one species per grid)\n";
2 
3 /*
4    Support 2.5V with axisymmetric coordinates
5      - r,z coordinates
6      - Domain and species data input by Landau operator
7      - "radius" for each grid, normalized with electron thermal velocity
8      - Domain: (0,radius) x (-radius,radius), thus first coordinate x[0] is perpendicular velocity and 2pi*x[0] term is added for axisymmetric
9    Supports full 3V
10 
11  */
12 
13 #include <petscdmplex.h>
14 #include <petscds.h>
15 #include <petscdmswarm.h>
16 #include <petscksp.h>
17 #include <petsc/private/petscimpl.h>
18 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
19   #include <omp.h>
20 #endif
21 #include <petsclandau.h>
22 #include <petscdmcomposite.h>
23 
24 typedef struct {
25   Mat MpTrans;
26   Mat Mp;
27   Vec ff;
28   Vec uu;
29 } MatShellCtx;
30 
31 typedef struct {
32   PetscInt   v_target;
33   PetscInt   g_target;
34   PetscInt   global_vertex_id_0;
35   DM        *globSwarmArray;
36   LandauCtx *ctx;
37   DM        *grid_dm;
38   Mat       *g_Mass;
39   Mat       *globMpArray;
40   Vec       *globXArray;
41   PetscBool  print;
42   PetscBool  print_entropy;
43 } PrintCtx;
44 
45 PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy)
46 {
47   MatShellCtx *matshellctx;
48 
49   PetscFunctionBeginUser;
50   PetscCall(MatShellGetContext(MtM, &matshellctx));
51   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
52   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
53   PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
54   PetscFunctionReturn(PETSC_SUCCESS);
55 }
56 
57 PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz)
58 {
59   MatShellCtx *matshellctx;
60 
61   PetscFunctionBeginUser;
62   PetscCall(MatShellGetContext(MtM, &matshellctx));
63   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
64   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
65   PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
66   PetscFunctionReturn(PETSC_SUCCESS);
67 }
68 
69 PetscErrorCode createSwarm(const DM dm, PetscInt dim, DM *sw)
70 {
71   PetscInt Nc = 1;
72 
73   PetscFunctionBeginUser;
74   PetscCall(DMCreate(PETSC_COMM_SELF, sw));
75   PetscCall(DMSetType(*sw, DMSWARM));
76   PetscCall(DMSetDimension(*sw, dim));
77   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
78   PetscCall(DMSwarmSetCellDM(*sw, dm));
79   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_REAL));
80   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
81   PetscCall(DMSetFromOptions(*sw));
82   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particle Grid"));
83   PetscFunctionReturn(PETSC_SUCCESS);
84 }
85 
86 static PetscErrorCode makeSwarm(DM sw, const PetscInt dim, const PetscInt Np, const PetscReal xx[], const PetscReal yy[], const PetscReal zz[])
87 {
88   PetscReal    *coords;
89   PetscDataType dtype;
90   PetscInt      bs, p, zero = 0;
91 
92   PetscFunctionBeginUser;
93   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
94   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
95   for (p = 0; p < Np; p++) {
96     coords[p * dim + 0] = xx[p];
97     coords[p * dim + 1] = yy[p];
98     if (dim == 3) coords[p * dim + 2] = zz[p];
99   }
100   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
101   PetscFunctionReturn(PETSC_SUCCESS);
102 }
103 
104 static PetscErrorCode createMp(const DM dm, DM sw, Mat *Mp_out)
105 {
106   PetscBool removePoints = PETSC_TRUE;
107   Mat       M_p;
108 
109   PetscFunctionBeginUser;
110   // migrate after coords are set
111   PetscCall(DMSwarmMigrate(sw, removePoints));
112   //
113   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
114 
115   /* PetscInt  N,*count,nmin=10000,nmax=0,ntot=0; */
116   /* // count */
117   /* PetscCall(DMSwarmCreatePointPerCellCount(sw, &N, &count)); */
118   /* for (int i=0, n; i< N ; i++) { */
119   /*   if ((n=count[i]) > nmax) nmax = n; */
120   /*   if (n < nmin) nmin = n; */
121   /*   PetscCall(PetscInfo(dm, " %d) %d particles\n", i, n)); */
122   /*   ntot += n; */
123   /* } */
124   /* PetscCall(PetscFree(count)); */
125   /* PetscCall(PetscInfo(dm, " %" PetscInt_FMT " max particle / cell, and %" PetscInt_FMT " min, ratio = %g,  %" PetscInt_FMT " total\n", nmax, nmin, (double)nmax/(double)nmin,ntot)); */
126 
127   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
128   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
129   PetscCall(DMViewFromOptions(sw, NULL, "-ex30_sw_view"));
130   // output
131   *Mp_out = M_p;
132   PetscFunctionReturn(PETSC_SUCCESS);
133 }
134 
135 static PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt a_tid, const PetscInt dim, const PetscReal a_wp[], Vec rho, Mat M_p)
136 {
137   PetscReal    *wq;
138   PetscDataType dtype;
139   Vec           ff;
140   PetscInt      bs, p, Np;
141 
142   PetscFunctionBeginUser;
143   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
144   PetscCall(DMSwarmGetLocalSize(sw, &Np));
145   for (p = 0; p < Np; p++) wq[p] = a_wp[p];
146   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
147   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
148   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff));
149   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
150   PetscCall(MatMultTranspose(M_p, ff, rho));
151   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
152   PetscFunctionReturn(PETSC_SUCCESS);
153 }
154 
155 //
156 // add grid to arg 'sw.w_q'
157 //
158 PetscErrorCode gridToParticles(const DM dm, DM sw, const Vec rhs, Vec work_ferhs, Mat M_p, Mat Mass)
159 {
160   PetscBool    is_lsqr;
161   KSP          ksp;
162   Mat          PM_p = NULL, MtM, D = NULL;
163   Vec          ff;
164   PetscInt     N, M, nzl;
165   MatShellCtx *matshellctx = NULL;
166   PC           pc;
167 
168   PetscFunctionBeginUser;
169   // 1) apply M in, for Moore-Penrose with mass: Mp (Mp' Mp)^-1 M
170   PetscCall(MatMult(Mass, rhs, work_ferhs));
171   // 2) pseudo-inverse, first part: (Mp' Mp)^-1
172   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
173   PetscCall(KSPSetType(ksp, KSPCG));
174   PetscCall(KSPGetPC(ksp, &pc));
175   PetscCall(PCSetType(pc, PCJACOBI));
176   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
177   PetscCall(KSPSetFromOptions(ksp));
178   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
179   if (!is_lsqr) {
180     PetscCall(MatGetLocalSize(M_p, &M, &N));
181     if (N > M) {
182       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") more vertices than particles: revert to lsqr\n", M, N));
183       is_lsqr = PETSC_TRUE;
184       PetscCall(KSPSetType(ksp, KSPLSQR));
185       PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp^T Mp), move projection Mp before solve
186     } else {
187       PetscCall(PetscNew(&matshellctx));
188       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
189       if (0) {
190         PetscCall(MatTransposeMatMult(M_p, M_p, MAT_INITIAL_MATRIX, 4, &MtM));
191         PetscCall(KSPSetOperators(ksp, MtM, MtM));
192         PetscCall(PetscInfo(M_p, "createMtM KSP with explicit Mp'Mp\n"));
193         PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
194       } else {
195         PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
196         PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
197         matshellctx->Mp = M_p;
198         PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
199         PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
200         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
201         PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpT_mat_view"));
202         for (int i = 0; i < N; i++) {
203           const PetscScalar *vals;
204           const PetscInt    *cols;
205           PetscScalar        dot = 0;
206           PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
207           for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
208           if (dot < PETSC_MACHINE_EPSILON) {
209             PetscCall(PetscInfo(ksp, "empty row in pseudo-inverse %d\n", i));
210             is_lsqr = PETSC_TRUE; // empty rows
211             PetscCall(KSPSetType(ksp, KSPLSQR));
212             PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve
213             // clean up
214             PetscCall(MatDestroy(&matshellctx->MpTrans));
215             PetscCall(VecDestroy(&matshellctx->ff));
216             PetscCall(VecDestroy(&matshellctx->uu));
217             PetscCall(MatDestroy(&D));
218             PetscCall(MatDestroy(&MtM));
219             PetscCall(PetscFree(matshellctx));
220             D = NULL;
221             break;
222           }
223           PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
224         }
225         if (D) {
226           PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
227           PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
228           PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
229           PetscCall(KSPSetOperators(ksp, MtM, D));
230           PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
231           PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
232           PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
233           PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
234         }
235       }
236     }
237   }
238   if (is_lsqr) {
239     PC        pc2;
240     PetscBool is_bjac;
241     PetscCall(KSPGetPC(ksp, &pc2));
242     PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac));
243     if (is_bjac) {
244       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
245       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
246     } else {
247       PetscCall(KSPSetOperators(ksp, M_p, M_p));
248     }
249     PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
250   }
251   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access
252   if (!is_lsqr) {
253     PetscErrorCode ierr;
254     ierr = KSPSolve(ksp, work_ferhs, matshellctx->uu);
255     if (!ierr) {
256       // 3) with Moore-Penrose apply Mp: M_p (Mp' Mp)^-1 M
257       PetscCall(MatMult(M_p, matshellctx->uu, ff));
258     } else { // failed
259       PC        pc2;
260       PetscBool is_bjac;
261       PetscCall(PetscInfo(ksp, "Solver failed, probably singular, try lsqr\n"));
262       PetscCall(KSPReset(ksp));
263       PetscCall(KSPSetType(ksp, KSPLSQR));
264       PetscCall(KSPGetPC(ksp, &pc2));
265       PetscCall(PCSetType(pc2, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve
266       PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
267       PetscCall(KSPSetFromOptions(ksp));
268       PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac));
269       if (is_bjac) {
270         PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
271         PetscCall(KSPSetOperators(ksp, M_p, PM_p));
272       } else {
273         PetscCall(KSPSetOperators(ksp, M_p, M_p));
274       }
275       ierr = KSPSolveTranspose(ksp, work_ferhs, ff);
276       if (ierr) { PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "backup LSQR solver failed - need to add N_v > N_p Moore-Penrose pseudo-inverse"); }
277     }
278     if (D) PetscCall(MatDestroy(&D));
279     PetscCall(MatDestroy(&MtM));
280     if (matshellctx->MpTrans) PetscCall(MatDestroy(&matshellctx->MpTrans));
281     PetscCall(VecDestroy(&matshellctx->ff));
282     PetscCall(VecDestroy(&matshellctx->uu));
283     PetscCall(PetscFree(matshellctx));
284   } else {
285     PetscErrorCode ierr;
286     // finally with LSQR apply M_p^\dagger
287     ierr = KSPSolveTranspose(ksp, work_ferhs, ff);
288     if (ierr) { PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "backup LSQR solver failed - need to add N_v > N_p Moore-Penrose pseudo-inverse"); }
289   }
290   PetscCall(KSPDestroy(&ksp));
291   PetscCall(MatDestroy(&PM_p));
292   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
293   PetscFunctionReturn(PETSC_SUCCESS);
294 }
295 
296 #define EX30_MAX_NUM_THRDS 12
297 #define EX30_MAX_BATCH_SZ  1024
298 //
299 // add grid to arg 'globSwarmArray[].w_q'
300 //
301 PetscErrorCode gridToParticles_private(DM grid_dm[], DM globSwarmArray[], const PetscInt dim, const PetscInt v_target, const PetscInt numthreads, const PetscInt num_vertices, const PetscInt global_vertex_id, Mat globMpArray[], Mat g_Mass[], Vec t_fhat[][EX30_MAX_NUM_THRDS], PetscReal moments[], Vec globXArray[], LandauCtx *ctx)
302 {
303   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
304 
305   PetscFunctionBeginUser;
306   // map back to particles
307   for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
308     PetscCall(PetscInfo(grid_dm[0], "g2p: global batch %" PetscInt_FMT " of %" PetscInt_FMT ", Landau batch %" PetscInt_FMT " of %" PetscInt_FMT ": map back to particles\n", global_vertex_id + 1, num_vertices, v_id_0 + 1, ctx->batch_sz));
309     //PetscPragmaOMP(parallel for)
310     for (int tid = 0; tid < numthreads; tid++) {
311       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
312       if (glb_v_id < num_vertices) {
313         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
314           PetscErrorCode ierr_t;
315           ierr_t = PetscInfo(grid_dm[0], "gridToParticles: global batch %" PetscInt_FMT ", local batch b=%" PetscInt_FMT ", grid g=%" PetscInt_FMT ", index(b,g) %" PetscInt_FMT "\n", global_vertex_id, v_id, grid, LAND_PACK_IDX(v_id, grid));
316           ierr_t = gridToParticles(grid_dm[grid], globSwarmArray[LAND_PACK_IDX(v_id, grid)], globXArray[LAND_PACK_IDX(v_id, grid)], t_fhat[grid][tid], globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]);
317           if (ierr_t) ierr = ierr_t;
318         }
319       }
320     }
321     PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
322     /* Get moments */
323     PetscCall(PetscInfo(grid_dm[0], "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", v_id_0, v_id_0 + numthreads));
324     for (int tid = 0; tid < numthreads; tid++) {
325       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
326       if (glb_v_id == v_target) {
327         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
328           PetscDataType dtype;
329           PetscReal    *wp, *coords;
330           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
331           PetscInt      npoints, bs = 1;
332           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
333           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
334           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
335           for (int p = 0; p < npoints; p++) {
336             PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1, w = fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
337             for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]);
338             moments[0] += w;
339             moments[1] += w * ctx->v_0 * coords[p * dim + 1]; // z-momentum
340             moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
341           }
342           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
343           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
344         }
345         const PetscReal N_inv = 1 / moments[0];
346         PetscCall(PetscInfo(grid_dm[0], "gridToParticles_private [%" PetscInt_FMT "], n = %g\n", v_id, (double)moments[0]));
347         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
348           PetscDataType dtype;
349           PetscReal    *wp, *coords;
350           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
351           PetscInt      npoints, bs = 1;
352           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
353           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
354           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
355           for (int p = 0; p < npoints; p++) {
356             const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1, w = fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
357             if (w > PETSC_REAL_MIN) {
358               moments[3] -= ww * PetscLogReal(ww);
359               PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ww (%g) > 1", (double)ww);
360             } else moments[4] -= w; // keep track of density that is lost
361           }
362           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
363           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
364         }
365       }
366     } // thread batch
367   } // batch
368   PetscFunctionReturn(PETSC_SUCCESS);
369 }
370 
371 static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscReal shift, PetscScalar *u)
372 {
373   PetscInt  i;
374   PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */
375 
376   if (shift != 0.) {
377     v2 = 0;
378     for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
379     v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
380     /* evaluate the shifted Maxwellian */
381     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
382   } else {
383     /* compute the exponents, v^2 */
384     for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
385     /* evaluate the Maxwellian */
386     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
387   }
388 }
389 
390 static PetscErrorCode PostStep(TS ts)
391 {
392   PetscInt   n, dim, nDMs, v_id;
393   PetscReal  t;
394   LandauCtx *ctx;
395   Vec        X;
396   PrintCtx  *printCtx;
397   DM         pack;
398   PetscReal  moments[5], e_grid[LANDAU_MAX_GRIDS];
399 
400   PetscFunctionBeginUser;
401   PetscCall(TSGetApplicationContext(ts, &printCtx));
402   if (!printCtx->print && !printCtx->print_entropy) PetscFunctionReturn(PETSC_SUCCESS);
403   ctx = printCtx->ctx;
404   if (printCtx->v_target < printCtx->global_vertex_id_0 || printCtx->v_target >= printCtx->global_vertex_id_0 + ctx->batch_sz) PetscFunctionReturn(PETSC_SUCCESS);
405   for (int i = 0; i < 5; i++) moments[i] = 0;
406   for (int i = 0; i < LANDAU_MAX_GRIDS; i++) e_grid[i] = 0;
407   v_id = printCtx->v_target % ctx->batch_sz;
408   PetscCall(TSGetDM(ts, &pack));
409   PetscCall(DMGetDimension(pack, &dim));
410   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
411   PetscCall(TSGetSolution(ts, &X));
412   PetscCall(TSGetStepNumber(ts, &n));
413   PetscCall(TSGetTime(ts, &t));
414   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
415   if (printCtx->print_entropy && printCtx->v_target >= 0 && 0) {
416     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
417       PetscDataType dtype;
418       PetscReal    *wp, *coords;
419       DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
420       Vec           work, subX = printCtx->globXArray[LAND_PACK_IDX(v_id, grid)];
421       PetscInt      bs, NN;
422       // C-G moments
423       PetscCall(VecDuplicate(subX, &work));
424       PetscCall(gridToParticles(printCtx->grid_dm[grid], sw, subX, work, printCtx->globMpArray[LAND_PACK_IDX(v_id, grid)], printCtx->g_Mass[grid]));
425       PetscCall(VecDestroy(&work));
426       // moments
427       PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
428       PetscCall(DMSwarmGetLocalSize(sw, &NN));
429       PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
430       for (int pp = 0; pp < NN; pp++) {
431         PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
432         for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
433         moments[0] += w;
434         moments[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
435         moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
436         e_grid[grid] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
437       }
438       PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
439       PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
440     }
441     // entropy
442     const PetscReal N_inv = 1 / moments[0];
443     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
444       PetscDataType dtype;
445       PetscReal    *wp, *coords;
446       DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
447       PetscInt      bs, NN;
448       PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
449       PetscCall(DMSwarmGetLocalSize(sw, &NN));
450       PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
451       for (int pp = 0; pp < NN; pp++) {
452         PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
453         if (w > PETSC_REAL_MIN) {
454           moments[3] -= ww * PetscLogReal(ww);
455         } else moments[4] -= w;
456       }
457       PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
458       PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
459     }
460     PetscCall(PetscInfo(X, "%4d) time %e, Landau particle moments: 0: %18.12e 1: %19.12e 2: %18.12e entropy: %e loss %e. energy = %e + %e + %e\n", (int)n, (double)t, (double)moments[0], (double)moments[1], (double)moments[2], (double)moments[3], (double)(moments[4] / moments[0]), (double)e_grid[0], (double)e_grid[1], (double)e_grid[2]));
461   }
462   if (printCtx->print && printCtx->g_target >= 0) {
463     PetscInt         grid   = printCtx->g_target, id;
464     static PetscReal last_t = -100000, period = .5;
465     if (last_t == -100000) last_t = -period + t;
466     if (t >= last_t + period) {
467       last_t = t;
468       PetscCall(DMGetOutputSequenceNumber(ctx->plex[grid], &id, NULL));
469       PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid], id + 1, t));
470       PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid)], NULL, "-ex30_vec_view"));
471       if (ctx->num_grids > grid + 1) {
472         PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid + 1], id + 1, t));
473         PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid + 1)], NULL, "-ex30_vec_view2"));
474       }
475       PetscCall(PetscInfo(X, "%4d) time %e View\n", (int)n, (double)t));
476     }
477   }
478   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
479   PetscFunctionReturn(PETSC_SUCCESS);
480 }
481 
482 PetscErrorCode go(TS ts, Vec X, const PetscInt num_vertices, const PetscInt a_Np, const PetscInt dim, const PetscInt v_target, const PetscInt g_target, PetscReal shift, PetscBool use_uniform_particle_grid)
483 {
484   DM             pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS];
485   Mat           *globMpArray, g_Mass[LANDAU_MAX_GRIDS];
486   KSP            t_ksp[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
487   Vec            t_fhat[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
488   PetscInt       nDMs;
489   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
490 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
491   PetscInt numthreads = PetscNumOMPThreads;
492 #else
493   PetscInt numthreads = 1;
494 #endif
495   LandauCtx *ctx;
496   Vec       *globXArray;
497   PetscReal  moments_0[5], moments_1a[5], moments_1b[5], dt_init;
498   PrintCtx  *printCtx;
499 
500   PetscFunctionBeginUser;
501   PetscCheck(numthreads <= EX30_MAX_NUM_THRDS, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
502   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
503   PetscCall(TSGetDM(ts, &pack));
504   PetscCall(DMGetApplicationContext(pack, &ctx));
505   PetscCheck(ctx->batch_sz % numthreads == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "batch size (-dm_landau_batch_size) %" PetscInt_FMT "  mod #threads %" PetscInt_FMT " must equal zero", ctx->batch_sz, numthreads);
506   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
507   PetscCall(PetscInfo(pack, "Have %" PetscInt_FMT " total grids, with %" PetscInt_FMT " Landau local batched and %" PetscInt_FMT " global items (vertices) %d DMs\n", ctx->num_grids, ctx->batch_sz, num_vertices, (int)nDMs));
508   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
509   PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray));
510   PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray));
511   // print ctx
512   PetscCall(PetscNew(&printCtx));
513   PetscCall(TSSetApplicationContext(ts, printCtx));
514   printCtx->v_target       = v_target;
515   printCtx->g_target       = g_target;
516   printCtx->ctx            = ctx;
517   printCtx->globSwarmArray = globSwarmArray;
518   printCtx->grid_dm        = grid_dm;
519   printCtx->globMpArray    = globMpArray;
520   printCtx->g_Mass         = g_Mass;
521   printCtx->globXArray     = globXArray;
522   printCtx->print_entropy  = PETSC_FALSE;
523   PetscOptionsBegin(PETSC_COMM_SELF, "", "Print Options", "DMPLEX");
524   PetscCall(PetscOptionsBool("-print_entropy", "Print entropy and moments at each time step", "ex30.c", printCtx->print_entropy, &printCtx->print_entropy, NULL));
525   PetscOptionsEnd();
526   // view
527   PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view"));
528   if (ctx->num_grids > g_target + 1) { PetscCall(DMViewFromOptions(ctx->plex[g_target + 1], NULL, "-ex30_dm_view2")); }
529   // create mesh mass matrices
530   PetscCall(VecZeroEntries(X));
531   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate
532   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {               // add same particels for all grids
533     Vec          subX = globXArray[LAND_PACK_IDX(0, grid)];
534     DM           dm   = ctx->plex[grid];
535     PetscSection s;
536     grid_dm[grid] = dm;
537     PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid]));
538     //
539     PetscCall(DMGetLocalSection(dm, &s));
540     PetscCall(DMPlexCreateClosureIndex(dm, s));
541     for (int tid = 0; tid < numthreads; tid++) {
542       PC pc;
543       PetscCall(VecDuplicate(subX, &t_fhat[grid][tid]));
544       PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid]));
545       PetscCall(KSPSetType(t_ksp[grid][tid], KSPCG));
546       PetscCall(KSPGetPC(t_ksp[grid][tid], &pc));
547       PetscCall(PCSetType(pc, PCJACOBI));
548       PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_"));
549       PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid]));
550       PetscCall(KSPSetFromOptions(t_ksp[grid][tid]));
551     }
552   }
553   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
554   PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper
555   // loop over all vertices in chucks that are batched for TSSolve
556   for (int i = 0; i < 5; i++) moments_0[i] = moments_1a[i] = moments_1b[i] = 0;
557   for (PetscInt global_vertex_id_0 = 0; global_vertex_id_0 < num_vertices; global_vertex_id_0 += ctx->batch_sz, shift /= 2) { // outer vertex loop
558     PetscCall(TSSetTime(ts, 0));
559     PetscCall(TSSetStepNumber(ts, 0));
560     PetscCall(TSSetTimeStep(ts, dt_init));
561     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
562     printCtx->global_vertex_id_0 = global_vertex_id_0;
563     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
564       PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "rho"));
565       printCtx->print = PETSC_TRUE;
566     } else printCtx->print = PETSC_FALSE;
567     // create fake particles in batches with threads
568     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
569       PetscReal *xx_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS] /* , radiuses[80000] */;
570       PetscInt   Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
571       // make particles
572       for (int tid = 0; tid < numthreads; tid++) {
573         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
574         if (glb_v_id < num_vertices) {                                                     // the ragged edge (in last batch)
575           PetscInt Npp0 = a_Np + (glb_v_id % (a_Np / 10 + 1)), nTargetP[LANDAU_MAX_GRIDS]; // n of particels in each dim with load imbalance
576           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {                         // add same particels for all grids
577             // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) {
578             const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ctx->species_offset[grid]] / ctx->masses[ctx->species_offset[grid]] / (ctx->v_0 * ctx->v_0);                      /* theta = 2kT/mc^2 per species */
579             PetscReal       lo[3] = {-ctx->radius[grid], -ctx->radius[grid], -ctx->radius[grid]}, hi[3] = {ctx->radius[grid], ctx->radius[grid], ctx->radius[grid]}, hp[3], vole; // would be nice to get box from DM
580             PetscInt        Npi = Npp0, Npj = 2 * Npp0, Npk = 1;
581             PetscRandom     rand;
582             PetscReal       sigma = ctx->thermal_speed[grid] / ctx->thermal_speed[0], p2_shift = grid == 0 ? shift : -shift; // symmetric shift of e vs ions
583             PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
584             PetscCall(PetscRandomSetInterval(rand, 0., 1.));
585             PetscCall(PetscRandomSetFromOptions(rand));
586             if (dim == 2) lo[0] = 0; // Landau coordinate (r,z)
587             else Npi = Npj = Npk = Npp0;
588             // User: use glb_v_id to index into your data
589             const PetscInt NNreal = Npi * Npj * Npk, NN = NNreal + (dim == 2 ? 3 : 6); // make room for bounding box
590             Np_t[grid][tid] = NN;
591             if (glb_v_id == v_target) nTargetP[grid] = NN;
592             PetscCall(PetscMalloc4(NN, &xx_t[grid][tid], NN, &yy_t[grid][tid], NN, &wp_t[grid][tid], dim == 2 ? 1 : NN, &zz_t[grid][tid]));
593             hp[0] = (hi[0] - lo[0]) / Npi;
594             hp[1] = (hi[1] - lo[1]) / Npj;
595             hp[2] = (hi[2] - lo[2]) / Npk;
596             if (dim == 2) hp[2] = 1;
597             PetscCall(PetscInfo(pack, " lo = %14.7e, hi = %14.7e; hp = %14.7e, %14.7e; kT_m = %g; \n", (double)lo[1], (double)hi[1], (double)hp[0], (double)hp[1], (double)kT_m)); // temp
598             vole = hp[0] * hp[1] * hp[2] * ctx->n[grid];                                                                                                                           // fix for multi-species
599             PetscCall(PetscInfo(pack, "Vertex %" PetscInt_FMT ", grid %" PetscInt_FMT " with %" PetscInt_FMT " particles (diagnostic target = %" PetscInt_FMT ")\n", glb_v_id, grid, NN, v_target));
600             for (int pj = 0, pp = 0; pj < Npj; pj++) {
601               for (int pk = 0; pk < Npk; pk++) {
602                 for (int pi = 0; pi < Npi; pi++, pp++) {
603                   PetscReal p_shift   = p2_shift;
604                   wp_t[grid][tid][pp] = 0;
605                   if (use_uniform_particle_grid) {
606                     xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
607                     yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
608                     if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
609                     PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
610                     p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
611                     if (ctx->sphere && PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) {
612                       wp_t[grid][tid][pp] = 0;
613                     } else {
614                       maxwellian(dim, x, kT_m, vole, p_shift, &wp_t[grid][tid][pp]);
615                       if (ctx->num_grids == 1 && shift != 0) {                          // bi-maxwellian, electron plasma
616                         maxwellian(dim, x, kT_m, vole, -p_shift, &wp_t[grid][tid][pp]); // symmetric shift of electron plasma
617                       }
618                     }
619                   } else {
620                     PetscReal u1, u2;
621                     do {
622                       do {
623                         PetscCall(PetscRandomGetValueReal(rand, &u1));
624                       } while (u1 == 0);
625                       PetscCall(PetscRandomGetValueReal(rand, &u2));
626                       //compute z0 and z1
627                       PetscReal mag       = sigma * PetscSqrtReal(-2.0 * PetscLogReal(u1)); // is this the same scale grid Maxwellian? t_therm = sigma
628                       xx_t[grid][tid][pp] = mag * PetscCosReal(2.0 * PETSC_PI * u2);
629                       yy_t[grid][tid][pp] = mag * PetscSinReal(2.0 * PETSC_PI * u2);
630                       if (dim == 2 && xx_t[grid][tid][pp] < lo[0]) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp];
631                       if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
632                       if (!ctx->sphere) {
633                         if (dim == 2 && xx_t[grid][tid][pp] < 0) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp]; // ???
634                         else if (dim == 3) {
635                           while (zz_t[grid][tid][pp] >= hi[2] || zz_t[grid][tid][pp] <= lo[2]) zz_t[grid][tid][pp] *= .9;
636                         }
637                         while (xx_t[grid][tid][pp] >= hi[0] || xx_t[grid][tid][pp] <= lo[0]) xx_t[grid][tid][pp] *= .9;
638                         while (yy_t[grid][tid][pp] >= hi[1] || yy_t[grid][tid][pp] <= lo[1]) yy_t[grid][tid][pp] *= .9;
639                       } else { // 2D
640                         //if (glb_v_id == v_target && pp < 80000) radiuses[pp] = PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp]));
641                         while (PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) { // safety factor for facets of sphere
642                           xx_t[grid][tid][pp] *= .9;
643                           yy_t[grid][tid][pp] *= .9;
644                         }
645                       }
646                       if (ctx->num_grids == 1 && pp % 2 == 0) p_shift = 0; // one species, split bi-max
647                       p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
648                       if (dim == 3) zz_t[grid][tid][pp] += p_shift;
649                       else yy_t[grid][tid][pp] += p_shift;
650                       wp_t[grid][tid][pp] += ctx->n[grid] / NNreal * PetscSqrtReal(ctx->masses[ctx->species_offset[grid]] / ctx->masses[0]);
651                       if (p_shift <= 0) break; // add bi-max for electron plasma only
652                       p_shift = -p_shift;
653                     } while (ctx->num_grids == 1); // add bi-max for electron plasma only
654                   }
655                   {
656                     if (glb_v_id == v_target) {
657                       PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
658                       PetscReal v2 = 0, fact = dim == 2 ? 2.0 * PETSC_PI * x[0] : 1, w = fact * wp_t[grid][tid][pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
659                       for (int i = 0; i < dim; ++i) v2 += PetscSqr(x[i]);
660                       moments_0[0] += w;                   // not thread safe
661                       moments_0[1] += w * ctx->v_0 * x[1]; // z-momentum
662                       moments_0[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
663                     }
664                   }
665                 }
666               }
667             }
668             if (dim == 2) { // fix bounding box
669               int pp                = NNreal;
670               wp_t[grid][tid][pp]   = 0;
671               xx_t[grid][tid][pp]   = 1.e-7;
672               yy_t[grid][tid][pp++] = hi[1] - 5.e-7;
673               wp_t[grid][tid][pp]   = 0;
674               xx_t[grid][tid][pp]   = hi[0] - 5.e-7;
675               yy_t[grid][tid][pp++] = 0;
676               wp_t[grid][tid][pp]   = 0;
677               xx_t[grid][tid][pp]   = 1.e-7;
678               yy_t[grid][tid][pp++] = lo[1] + 5.e-7;
679             } else {
680               const int p0 = NNreal;
681               for (int pj = 0; pj < 6; pj++) { xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0; }
682               xx_t[grid][tid][p0 + 0] = lo[0];
683               xx_t[grid][tid][p0 + 1] = hi[0];
684               yy_t[grid][tid][p0 + 2] = lo[1];
685               yy_t[grid][tid][p0 + 3] = hi[1];
686               zz_t[grid][tid][p0 + 4] = lo[2];
687               zz_t[grid][tid][p0 + 5] = hi[2];
688             }
689             PetscCall(PetscRandomDestroy(&rand));
690           }
691           // entropy init, need global n
692           if (glb_v_id == v_target) {
693             const PetscReal N_inv = 1 / moments_0[0];
694             PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_v_id, nTargetP[0]));
695             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
696               const PetscInt NN = nTargetP[grid];
697               for (int pp = 0; pp < NN; pp++) {
698                 const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * xx_t[grid][tid][pp] : 1, w = fact * ctx->n_0 * ctx->masses[ctx->species_offset[grid]] * wp_t[grid][tid][pp], ww = w * N_inv;
699                 if (w > PETSC_REAL_MIN) {
700                   moments_0[3] -= ww * PetscLogReal(ww);
701                   PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
702                 } else moments_0[4] -= w;
703               }
704             } // grid
705           } // target
706         } // active
707       } // threads
708       /* Create particle swarm */
709       for (int tid = 0; tid < numthreads; tid++) {
710         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
711         if (glb_v_id < num_vertices) {                             // the ragged edge of the last batch
712           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
713             PetscSection section;
714             PetscInt     Nf;
715             DM           dm = grid_dm[grid];
716             PetscCall(DMGetLocalSection(dm, &section));
717             PetscCall(PetscSectionGetNumFields(section, &Nf));
718             PetscCheck(Nf == 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only on species per grid supported -- todo");
719             PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
720             PetscCall(PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)));
721             PetscCall(createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
722           }
723         } // active
724       } // threads
725       PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Only support one species per grid");
726       // make globMpArray
727       PetscPragmaOMP(parallel for)
728       for (int tid = 0; tid < numthreads; tid++) {
729         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
730         if (glb_v_id < num_vertices) {
731           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
732             // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
733             PetscErrorCode ierr_t;
734             DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
735             ierr_t            = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
736             ierr_t            = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
737             if (ierr_t) ierr = ierr_t;
738           }
739         }
740       }
741       for (int tid = 0; tid < numthreads; tid++) {
742         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
743         if (glb_v_id < num_vertices) {
744           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
745             DM dm = grid_dm[grid];
746             DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
747             PetscCall(PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)));
748             PetscCall(createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]));
749             PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_id, grid)], NULL, "-mp_mat_view"));
750           }
751         }
752       }
753       // p --> g: set X
754       // PetscPragmaOMP(parallel for)
755       for (int tid = 0; tid < numthreads; tid++) {
756         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
757         if (glb_v_id < num_vertices) {
758           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
759             PetscErrorCode ierr_t;
760             DM             dm   = grid_dm[grid];
761             DM             sw   = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
762             Vec            subX = globXArray[LAND_PACK_IDX(v_id, grid)], work = t_fhat[grid][tid];
763             ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
764             ierr_t = particlesToGrid(dm, sw, tid, dim, wp_t[grid][tid], subX, globMpArray[LAND_PACK_IDX(v_id, grid)]);
765             if (ierr_t) ierr = ierr_t;
766             // u = M^_1 f_w
767             ierr_t = VecCopy(subX, work);
768             ierr_t = KSPSolve(t_ksp[grid][tid], work, subX);
769             if (ierr_t) ierr = ierr_t;
770           }
771         }
772       }
773       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
774       /* Cleanup */
775       for (int tid = 0; tid < numthreads; tid++) {
776         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
777         if (glb_v_id < num_vertices) {
778           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
779             PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
780           }
781         } // active
782       } // threads
783     } // (fake) particle loop
784     // standard view of initial conditions
785     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
786       PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 0, 0.0));
787       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
788       if (ctx->num_grids > g_target + 1) {
789         PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target + 1], 0, 0.0));
790         PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target + 1)], NULL, "-ex30_vec_view2"));
791       }
792       PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_mass_mat_view"));
793       PetscCall(DMViewFromOptions(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_sw_view"));
794       PetscCall(DMSwarmViewXDMF(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "initial_swarm.xmf")); // writes a file by default!!!
795     }
796     // coarse graining moments_1a, bring f back from grid before advance
797     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
798       const PetscInt v_id = v_target % ctx->batch_sz;
799       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
800         PetscDataType dtype;
801         PetscReal    *wp, *coords;
802         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
803         Vec           work, subX = globXArray[LAND_PACK_IDX(v_id, grid)];
804         PetscInt      bs, NN;
805         // C-G moments
806         PetscCall(VecDuplicate(subX, &work));
807         PetscCall(gridToParticles(grid_dm[grid], sw, subX, work, globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]));
808         PetscCall(VecDestroy(&work));
809         // moments
810         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
811         PetscCall(DMSwarmGetLocalSize(sw, &NN));
812         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
813         for (int pp = 0; pp < NN; pp++) {
814           PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
815           for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
816           moments_1a[0] += w;
817           moments_1a[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
818           moments_1a[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
819         }
820         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
821         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
822       }
823       // entropy
824       const PetscReal N_inv = 1 / moments_1a[0];
825       PetscCall(PetscInfo(pack, "Entropy batch %" PetscInt_FMT " of %" PetscInt_FMT ", n = %g\n", v_target, num_vertices, (double)(1 / N_inv)));
826       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
827         PetscDataType dtype;
828         PetscReal    *wp, *coords;
829         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
830         PetscInt      bs, NN;
831         PetscCall(DMSwarmGetLocalSize(sw, &NN));
832         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
833         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
834         for (int pp = 0; pp < NN; pp++) {
835           PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
836           if (w > PETSC_REAL_MIN) {
837             moments_1a[3] -= ww * PetscLogReal(ww);
838             PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
839           } else moments_1a[4] -= w;
840         }
841         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
842         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
843       }
844     }
845     // restore vector
846     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
847     // view initial grid
848     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { PetscCall(DMPlexLandauPrintNorms(X, 0)); }
849     // advance
850     PetscCall(TSSetSolution(ts, X));
851     PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT "\n", global_vertex_id_0, global_vertex_id_0 + ctx->batch_sz));
852     PetscCall(TSSetPostStep(ts, PostStep));
853     PetscCall(PostStep(ts));
854     PetscCall(TSSolve(ts, X));
855     // view
856     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
857     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
858       /* Visualize original particle field */
859       DM  sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
860       Vec f;
861       PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
862       PetscCall(DMViewFromOptions(grid_dm[g_target], NULL, "-weights_dm_view"));
863       PetscCall(DMViewFromOptions(sw, NULL, "-weights_sw_view"));
864       PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
865       PetscCall(PetscObjectSetName((PetscObject)f, "weights"));
866       PetscCall(VecViewFromOptions(f, NULL, "-weights_vec_view"));
867       PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
868       //
869       PetscCall(DMPlexLandauPrintNorms(X, 1));
870     }
871     if (!use_uniform_particle_grid) { // resample to uniform grid
872       for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
873         PetscReal *xx_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
874         PetscInt   Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
875         for (int tid = 0; tid < numthreads; tid++) {
876           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
877           if (glb_v_id < num_vertices) {
878             // create uniform grid w/o weights & smaller
879             PetscInt Npp0 = (a_Np + (glb_v_id % (a_Np / 10 + 1))) / 2, Nv; // 1/2 of uniform particle grid size
880             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
881               // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++)
882               PetscReal lo[3] = {-ctx->radius[grid], -ctx->radius[grid], -ctx->radius[grid]}, hi[3] = {ctx->radius[grid], ctx->radius[grid], ctx->radius[grid]}, hp[3];
883               PetscInt  Npi = Npp0, Npj = 2 * Npp0, Npk = 1, NN;
884               // delete old particles and particle mass matrix
885               PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
886               PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
887               // create fake particles in batches with threads
888               PetscCall(MatGetLocalSize(g_Mass[grid], &Nv, NULL));
889               if (dim == 2) lo[0] = 0;
890               else Npi = Npj = Npk = Npp0;
891               NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6); // make a regular grid of particles Npp x Npp
892               while (Npi * Npj * Npk < Nv) {             // make stable - no LS
893                 Npi++;
894                 Npj++;
895                 Npk++;
896                 NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6);
897               }
898               Np_t[grid][tid] = NN;
899               PetscCall(PetscMalloc4(NN, &xx_t[grid][tid], NN, &yy_t[grid][tid], NN, &wp_t[grid][tid], dim == 2 ? 1 : NN, &zz_t[grid][tid]));
900               hp[0] = (hi[0] - lo[0]) / Npi;
901               hp[1] = (hi[1] - lo[1]) / Npj;
902               hp[2] = (hi[2] - lo[2]) / Npk;
903               if (dim == 2) hp[2] = 1;
904               PetscCall(PetscInfo(pack, "Resampling %d particles, %d vertices\n", (int)NN, (int)Nv)); // temp
905               for (int pj = 0, pp = 0; pj < Npj; pj++) {
906                 for (int pk = 0; pk < Npk; pk++) {
907                   for (int pi = 0; pi < Npi; pi++, pp++) {
908                     wp_t[grid][tid][pp] = 0;
909                     xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
910                     yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
911                     if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
912                   }
913                 }
914               }
915               if (dim == 2) { // fix bounding box
916                 int pp                = NN - 3;
917                 wp_t[grid][tid][pp]   = 0;
918                 xx_t[grid][tid][pp]   = 1.e-7;
919                 yy_t[grid][tid][pp++] = hi[1] - 5.e-7;
920                 wp_t[grid][tid][pp]   = 0;
921                 xx_t[grid][tid][pp]   = hi[0] - 5.e-7;
922                 yy_t[grid][tid][pp++] = 0;
923                 wp_t[grid][tid][pp]   = 0;
924                 xx_t[grid][tid][pp]   = 1.e-7;
925                 yy_t[grid][tid][pp++] = lo[1] + 5.e-7;
926               } else {
927                 const int p0 = NN - 6;
928                 for (int pj = 0; pj < 6; pj++) { xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0; }
929                 xx_t[grid][tid][p0 + 0] = lo[0];
930                 xx_t[grid][tid][p0 + 1] = hi[0];
931                 yy_t[grid][tid][p0 + 2] = lo[1];
932                 yy_t[grid][tid][p0 + 3] = hi[1];
933                 zz_t[grid][tid][p0 + 4] = lo[2];
934                 zz_t[grid][tid][p0 + 5] = hi[2];
935               }
936             }
937           } // active
938         } // threads
939         /* Create particle swarm */
940         for (int tid = 0; tid < numthreads; tid++) {
941           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
942           if (glb_v_id < num_vertices) {                             // the ragged edge of the last batch
943             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
944               // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
945               PetscErrorCode ierr_t;
946               PetscSection   section;
947               PetscInt       Nf;
948               DM             dm = grid_dm[grid];
949               ierr_t            = DMGetLocalSection(dm, &section);
950               ierr_t            = PetscSectionGetNumFields(section, &Nf);
951               if (Nf != 1) ierr_t = (PetscErrorCode)9999;
952               else {
953                 ierr_t = DMViewFromOptions(dm, NULL, "-dm_view");
954                 ierr_t = PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
955                 ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]);
956               }
957               if (ierr_t) ierr = ierr_t;
958             }
959           } // active
960         } // threads
961         PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid");
962         PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
963         // make globMpArray
964         PetscPragmaOMP(parallel for)
965         for (int tid = 0; tid < numthreads; tid++) {
966           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
967           if (glb_v_id < num_vertices) {
968             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
969               // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
970               PetscErrorCode ierr_t;
971               DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
972               ierr_t            = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
973               ierr_t            = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
974               if (ierr_t) ierr = ierr_t;
975             }
976           } // active
977         } // threads
978         // create particle mass matrices
979         //PetscPragmaOMP(parallel for)
980         for (int tid = 0; tid < numthreads; tid++) {
981           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
982           if (glb_v_id < num_vertices) {
983             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
984               PetscErrorCode ierr_t;
985               DM             dm = grid_dm[grid];
986               DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
987               ierr_t            = PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
988               ierr_t            = createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]);
989               if (ierr_t) ierr = ierr_t;
990             }
991           } // active
992         } // threads
993         PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
994         /* Cleanup */
995         for (int tid = 0; tid < numthreads; tid++) {
996           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
997           if (glb_v_id < num_vertices) {
998             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
999               PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
1000             }
1001           } // active
1002         } // threads
1003       } // batch
1004       // view
1005       if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
1006         /* Visualize particle field */
1007         DM  sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
1008         Vec f;
1009         PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
1010         PetscCall(DMViewFromOptions(sw, NULL, "-resampled_weights_sw_view"));
1011         PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
1012         PetscCall(PetscObjectSetName((PetscObject)f, "resampled_weights"));
1013         PetscCall(VecViewFromOptions(f, NULL, "-resampled_weights_vec_view"));
1014         PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
1015         PetscCall(DMSwarmViewXDMF(sw, "resampled.xmf"));
1016       }
1017     } // !uniform
1018     // particles to grid, compute moments and entropy, for target vertex only
1019     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
1020       PetscReal energy_error_rel;
1021       PetscCall(gridToParticles_private(grid_dm, globSwarmArray, dim, v_target, numthreads, num_vertices, global_vertex_id_0, globMpArray, g_Mass, t_fhat, moments_1b, globXArray, ctx));
1022       energy_error_rel = PetscAbsReal((moments_1b[2] - moments_0[2])) / moments_0[2];
1023       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Particle Moments:\t number density      momentum (par)     energy             entropy            negative weights  : # OMP threads %g\n", (double)numthreads));
1024       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tInitial:         %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_0[0], (double)moments_0[1], (double)moments_0[2], (double)moments_0[3], 100 * (double)(moments_0[4] / moments_0[0])));
1025       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tCoarse-graining: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1a[0], (double)moments_1a[1], (double)moments_1a[2], (double)moments_1a[3], 100 * (double)(moments_1a[4] / moments_0[0])));
1026       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tLandau:          %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1b[0], (double)moments_1b[1], (double)moments_1b[2], (double)moments_1b[3], 100 * (double)(moments_1b[4] / moments_0[0])));
1027       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coarse-graining entropy generation = %e ; Landau entropy generation = %e\n", (double)(moments_1a[3] - moments_0[3]), (double)(moments_1b[3] - moments_0[3])));
1028       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(relative) energy conservation: Coarse-graining = %e, Landau = %e (%g %d)\n", (double)(moments_1a[2] - moments_0[2]) / (double)moments_0[2], (double)energy_error_rel, (double)PetscLog10Real(energy_error_rel), (int)(PetscLog10Real(energy_error_rel) + .5)));
1029     }
1030     // restore vector
1031     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
1032     // cleanup
1033     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
1034       for (int tid = 0; tid < numthreads; tid++) {
1035         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
1036         if (glb_v_id < num_vertices) {
1037           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1038             PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
1039             PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
1040           }
1041         }
1042       }
1043     }
1044   } // user batch, not used
1045   /* Cleanup */
1046   PetscCall(PetscFree(globXArray));
1047   PetscCall(PetscFree(globSwarmArray));
1048   PetscCall(PetscFree(globMpArray));
1049   PetscCall(PetscFree(printCtx));
1050   // clean up mass matrices
1051   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
1052     PetscCall(MatDestroy(&g_Mass[grid]));
1053     for (int tid = 0; tid < numthreads; tid++) {
1054       PetscCall(VecDestroy(&t_fhat[grid][tid]));
1055       PetscCall(KSPDestroy(&t_ksp[grid][tid]));
1056     }
1057   }
1058   PetscFunctionReturn(PETSC_SUCCESS);
1059 }
1060 
1061 int main(int argc, char **argv)
1062 {
1063   DM         pack;
1064   Vec        X;
1065   PetscInt   dim = 2, num_vertices = 1, Np = 10, v_target = 0, g_target = 0;
1066   TS         ts;
1067   Mat        J;
1068   LandauCtx *ctx;
1069   PetscReal  shift                     = 0;
1070   PetscBool  use_uniform_particle_grid = PETSC_TRUE;
1071 
1072   PetscFunctionBeginUser;
1073   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1074   // process args
1075   PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX");
1076   PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL));
1077   PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", num_vertices, &num_vertices, NULL));
1078   PetscCall(PetscOptionsInt("-number_particles_per_dimension", "Number of particles per grid, with slight modification per spatial vertex, in each dimension of base Cartesian grid", "ex30.c", Np, &Np, NULL));
1079   PetscCall(PetscOptionsBool("-use_uniform_particle_grid", "Use uniform particle grid", "ex30.c", use_uniform_particle_grid, &use_uniform_particle_grid, NULL));
1080   PetscCall(PetscOptionsInt("-vertex_view_target", "Global vertex for diagnostics", "ex30.c", v_target, &v_target, NULL));
1081   PetscCall(PetscOptionsReal("-e_shift", "Bi-Maxwellian shift", "ex30.c", shift, &shift, NULL));
1082   PetscCheck(v_target < num_vertices, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Batch to view %" PetscInt_FMT " should be < number of vertices %" PetscInt_FMT, v_target, num_vertices);
1083   PetscCall(PetscOptionsInt("-grid_view_target", "Grid to view with diagnostics", "ex30.c", g_target, &g_target, NULL));
1084   PetscOptionsEnd();
1085   /* Create a mesh */
1086   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
1087   PetscCall(DMGetApplicationContext(pack, &ctx));
1088   PetscCall(DMSetUp(pack));
1089   PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
1090   PetscCheck(g_target < ctx->num_grids, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Grid to view %" PetscInt_FMT " should be < number of grids %" PetscInt_FMT, g_target, ctx->num_grids);
1091   PetscCheck(ctx->batch_view_idx == v_target % ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Global view index %" PetscInt_FMT " mode batch size %" PetscInt_FMT " != ctx->batch_view_idx %" PetscInt_FMT, v_target, ctx->batch_sz, ctx->batch_view_idx);
1092   // PetscCheck(!use_uniform_particle_grid || !ctx->sphere, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Can not use -use_uniform_particle_grid and -dm_landau_sphere");
1093   /* Create timestepping solver context */
1094   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1095   PetscCall(TSSetDM(ts, pack));
1096   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
1097   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
1098   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1099   PetscCall(TSSetFromOptions(ts));
1100   PetscCall(PetscObjectSetName((PetscObject)X, "X"));
1101   // do particle advance
1102   PetscCall(go(ts, X, num_vertices, Np, dim, v_target, g_target, shift, use_uniform_particle_grid));
1103   PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic
1104   /* clean up */
1105   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
1106   PetscCall(TSDestroy(&ts));
1107   PetscCall(VecDestroy(&X));
1108   PetscCall(PetscFinalize());
1109   return 0;
1110 }
1111 
1112 /*TEST
1113 
1114   build:
1115     requires: !complex
1116 
1117   testset:
1118     requires: double defined(PETSC_USE_DMLANDAU_2D)
1119     output_file: output/ex30_0.out
1120     args: -dim 2 -petscspace_degree 3 -dm_landau_num_species_grid 1,1,1 -dm_refine 1 -number_particles_per_dimension 20 \
1121           -dm_landau_batch_size 4 -number_spatial_vertices 6 -vertex_view_target 5 -grid_view_target 1 -dm_landau_batch_view_idx 1 \
1122           -dm_landau_n 1.000018,1,1e-6 -dm_landau_thermal_temps 2,1,1 -dm_landau_ion_masses 2,180 -dm_landau_ion_charges 1,18 \
1123           -ftop_ksp_rtol 1e-10 -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_factor_shift_type nonzero -ftop_sub_pc_type lu -ftop_ksp_error_if_not_converged \
1124           -ksp_type gmres -ksp_error_if_not_converged -dm_landau_verbose 4 -print_entropy \
1125           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged\
1126           -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \
1127           -ts_dt 0.01 -ts_rtol 1e-1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler
1128     test:
1129       suffix: cpu
1130       args: -dm_landau_device_type cpu -pc_type jacobi
1131     test:
1132       suffix: kokkos
1133       # failed on Sunspot@ALCF with sycl
1134       requires: kokkos_kernels !openmp !sycl
1135       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi
1136 
1137   testset:
1138     requires: double !defined(PETSC_USE_DMLANDAU_2D)
1139     output_file: output/ex30_3d.out
1140     args: -dim 3 -petscspace_degree 2 -dm_landau_num_species_grid 1,1 -dm_refine 0 -number_particles_per_dimension 10 -dm_plex_hash_location \
1141           -dm_landau_batch_size 1 -number_spatial_vertices 1 -vertex_view_target 0 -grid_view_target 0 -dm_landau_batch_view_idx 0 \
1142           -dm_landau_n 1.000018,1 -dm_landau_thermal_temps 2,1 -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 \
1143           -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-12 -ftop_ksp_error_if_not_converged -ksp_type preonly -pc_type lu -ksp_error_if_not_converged \
1144           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged \
1145           -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \
1146           -ts_dt 0.1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler -print_entropy
1147     test:
1148       suffix: cpu_3d
1149       args: -dm_landau_device_type cpu
1150     test:
1151       suffix: kokkos_3d
1152       requires: kokkos_kernels !openmp
1153       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi
1154 
1155   test:
1156     suffix: conserve
1157     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
1158     args: -dm_landau_batch_size 4 -dm_refine 0 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_dt .1 \
1159           -ts_max_steps 1 -ksp_type preonly -ksp_error_if_not_converged -snes_rtol 1e-14 -snes_stol 1e-14 -dm_landau_device_type cpu -number_particles_per_dimension 20 \
1160           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-14 -ptof_ksp_error_if_not_converged -pc_type lu -dm_landau_simplex 1 -use_uniform_particle_grid false -dm_landau_sphere -print_entropy -number_particles_per_dimension 50 -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-14
1161 
1162 TEST*/
1163