xref: /petsc/src/ts/tests/ex28.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "Example application of the Bhatnagar-Gross-Krook (BGK) collision operator.\n\
2 This example is a 0D-1V setting for the kinetic equation\n\
3 https://en.wikipedia.org/wiki/Bhatnagar%E2%80%93Gross%E2%80%93Krook_operator\n";
4 
5 #include <petscdmplex.h>
6 #include <petscdmswarm.h>
7 #include <petscts.h>
8 #include <petscdraw.h>
9 #include <petscviewer.h>
10 
11 typedef struct {
12   PetscInt    particlesPerCell; /* The number of partices per cell */
13   PetscReal   momentTol;        /* Tolerance for checking moment conservation */
14   PetscBool   monitorhg;        /* Flag for using the TS histogram monitor */
15   PetscBool   monitorsp;        /* Flag for using the TS scatter monitor */
16   PetscBool   monitorks;        /* Monitor to perform KS test to the maxwellian */
17   PetscBool   error;            /* Flag for printing the error */
18   PetscInt    ostep;            /* print the energy at each ostep time steps */
19   PetscDraw   draw;             /* The draw object for histogram monitoring */
20   PetscDrawHG drawhg;           /* The histogram draw context for monitoring */
21   PetscDrawSP drawsp;           /* The scatter plot draw context for the monitor */
22   PetscDrawSP drawks;           /* Scatterplot draw context for KS test */
23 } AppCtx;
24 
25 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
26   PetscFunctionBeginUser;
27   options->monitorhg        = PETSC_FALSE;
28   options->monitorsp        = PETSC_FALSE;
29   options->monitorks        = PETSC_FALSE;
30   options->particlesPerCell = 1;
31   options->momentTol        = 100.0 * PETSC_MACHINE_EPSILON;
32   options->ostep            = 100;
33 
34   PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");
35   PetscCall(PetscOptionsBool("-monitorhg", "Flag to use the TS histogram monitor", "ex28.c", options->monitorhg, &options->monitorhg, NULL));
36   PetscCall(PetscOptionsBool("-monitorsp", "Flag to use the TS scatter plot monitor", "ex28.c", options->monitorsp, &options->monitorsp, NULL));
37   PetscCall(PetscOptionsBool("-monitorks", "Flag to plot KS test results", "ex28.c", options->monitorks, &options->monitorks, NULL));
38   PetscCall(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex28.c", options->particlesPerCell, &options->particlesPerCell, NULL));
39   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex28.c", options->ostep, &options->ostep, PETSC_NULL));
40   PetscOptionsEnd();
41 
42   PetscFunctionReturn(0);
43 }
44 
45 /* Create the mesh for velocity space */
46 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) {
47   PetscFunctionBeginUser;
48   PetscCall(DMCreate(comm, dm));
49   PetscCall(DMSetType(*dm, DMPLEX));
50   PetscCall(DMSetFromOptions(*dm));
51   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
52   PetscFunctionReturn(0);
53 }
54 
55 /* Since we are putting the same number of particles in each cell, this amounts to a uniform distribution of v */
56 static PetscErrorCode SetInitialCoordinates(DM sw) {
57   AppCtx        *user;
58   PetscRandom    rnd;
59   DM             dm;
60   DMPolytopeType ct;
61   PetscBool      simplex;
62   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *vals;
63   PetscInt       dim, d, cStart, cEnd, c, Np, p;
64 
65   PetscFunctionBeginUser;
66   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
67   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
68   PetscCall(PetscRandomSetFromOptions(rnd));
69 
70   PetscCall(DMGetApplicationContext(sw, &user));
71   Np = user->particlesPerCell;
72   PetscCall(DMGetDimension(sw, &dim));
73   PetscCall(DMSwarmGetCellDM(sw, &dm));
74   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
75   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
76   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
77   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
78   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
79   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
80   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals));
81   for (c = cStart; c < cEnd; ++c) {
82     if (Np == 1) {
83       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
84       for (d = 0; d < dim; ++d) {
85         coords[c * dim + d] = centroid[d];
86         if ((coords[c * dim + d] >= -1) && (coords[c * dim + d] <= 1)) {
87           vals[c] = 1.0;
88         } else {
89           vals[c] = 0.;
90         }
91       }
92     } else {
93       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
94       for (p = 0; p < Np; ++p) {
95         const PetscInt n   = c * Np + p;
96         PetscReal      sum = 0.0, refcoords[3];
97 
98         for (d = 0; d < dim; ++d) {
99           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
100           sum += refcoords[d];
101         }
102         if (simplex && sum > 0.0)
103           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
104         vals[n] = 1.0;
105         PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim]));
106       }
107     }
108   }
109   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
110   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals));
111   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
112   PetscCall(PetscRandomDestroy(&rnd));
113   PetscFunctionReturn(0);
114 }
115 
116 /* The intiial conditions are just the initial particle weights */
117 static PetscErrorCode SetInitialConditions(DM dmSw, Vec u) {
118   DM           dm;
119   AppCtx      *user;
120   PetscReal   *vals;
121   PetscScalar *initialConditions;
122   PetscInt     dim, d, cStart, cEnd, c, Np, p, n;
123 
124   PetscFunctionBeginUser;
125   PetscCall(VecGetLocalSize(u, &n));
126   PetscCall(DMGetApplicationContext(dmSw, &user));
127   Np = user->particlesPerCell;
128   PetscCall(DMSwarmGetCellDM(dmSw, &dm));
129   PetscCall(DMGetDimension(dm, &dim));
130   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
131   PetscCheck(n == (cEnd - cStart) * Np, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "TS solution local size %" PetscInt_FMT " != %" PetscInt_FMT " nm particles", n, (cEnd - cStart) * Np);
132   PetscCall(DMSwarmGetField(dmSw, "w_q", NULL, NULL, (void **)&vals));
133   PetscCall(VecGetArray(u, &initialConditions));
134   for (c = cStart; c < cEnd; ++c) {
135     for (p = 0; p < Np; ++p) {
136       const PetscInt n = c * Np + p;
137       for (d = 0; d < dim; d++) initialConditions[n] = vals[n];
138     }
139   }
140   PetscCall(VecRestoreArray(u, &initialConditions));
141   PetscCall(DMSwarmRestoreField(dmSw, "w_q", NULL, NULL, (void **)&vals));
142   PetscFunctionReturn(0);
143 }
144 
145 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) {
146   PetscInt *cellid;
147   PetscInt  dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
148 
149   PetscFunctionBeginUser;
150   PetscCall(DMGetDimension(dm, &dim));
151   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
152   PetscCall(DMSetType(*sw, DMSWARM));
153   PetscCall(DMSetDimension(*sw, dim));
154   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
155   PetscCall(DMSwarmSetCellDM(*sw, dm));
156   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL));
157   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
158   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
159   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
160   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
161   PetscCall(DMSetFromOptions(*sw));
162   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
163   for (c = cStart; c < cEnd; ++c) {
164     for (p = 0; p < Np; ++p) {
165       const PetscInt n = c * Np + p;
166       cellid[n]        = c;
167     }
168   }
169   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
170   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
171   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
172   PetscFunctionReturn(0);
173 }
174 
175 /*
176   f_t = 1/\tau \left( f_eq - f \right)
177   n_t = 1/\tau \left( \int f_eq - \int f \right) = 1/\tau (n - n) = 0
178   v_t = 1/\tau \left( \int v f_eq - \int v f \right) = 1/\tau (v - v) = 0
179   E_t = 1/\tau \left( \int v^2 f_eq - \int v^2 f \right) = 1/\tau (T - T) = 0
180 
181   Let's look at a single cell:
182 
183     \int_C f_t             = 1/\tau \left( \int_C f_eq - \int_C f \right)
184     \sum_{x_i \in C} w^i_t = 1/\tau \left( F_eq(C) - \sum_{x_i \in C} w_i \right)
185 */
186 
187 /* This computes the 1D Maxwellian distribution for given mass n, velocity v, and temperature T */
188 static PetscReal ComputePDF(PetscReal m, PetscReal n, PetscReal T, PetscReal v[]) {
189   return (n / PetscSqrtReal(2.0 * PETSC_PI * T / m)) * PetscExpReal(-0.5 * m * PetscSqr(v[0]) / T);
190 }
191 
192 /*
193   erf z = \frac{2}{\sqrt\pi} \int^z_0 dt e^{-t^2} and erf \infty = 1
194 
195   We begin with our distribution
196 
197     \sqrt{\frac{m}{2 \pi T}} e^{-m v^2/2T}
198 
199   Let t = \sqrt{m/2T} v, z = \sqrt{m/2T} w, so that we now have
200 
201       \sqrt{\frac{m}{2 \pi T}} \int^w_0 dv e^{-m v^2/2T}
202     = \sqrt{\frac{m}{2 \pi T}} \int^{\sqrt{m/2T} w}_0 \sqrt{2T/m} dt e^{-t^2}
203     = 1/\sqrt{\pi} \int^{\sqrt{m/2T} w}_0 dt e^{-t^2}
204     = 1/2 erf(\sqrt{m/2T} w)
205 */
206 static PetscReal ComputeCDF(PetscReal m, PetscReal n, PetscReal T, PetscReal va, PetscReal vb) {
207   PetscReal alpha = PetscSqrtReal(0.5 * m / T);
208   PetscReal za    = alpha * va;
209   PetscReal zb    = alpha * vb;
210   PetscReal sum   = 0.0;
211 
212   sum += zb >= 0 ? erf(zb) : -erf(-zb);
213   sum -= za >= 0 ? erf(za) : -erf(-za);
214   return 0.5 * n * sum;
215 }
216 
217 static PetscErrorCode CheckDistribution(DM dm, PetscReal m, PetscReal n, PetscReal T, PetscReal v[]) {
218   PetscSection coordSection;
219   Vec          coordsLocal;
220   PetscReal   *xq, *wq;
221   PetscReal    vmin, vmax, neq, veq, Teq;
222   PetscInt     Nq = 100, q, cStart, cEnd, c;
223 
224   PetscFunctionBeginUser;
225   PetscCall(DMGetBoundingBox(dm, &vmin, &vmax));
226   /* Check analytic over entire line */
227   neq = ComputeCDF(m, n, T, vmin, vmax);
228   PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", (double)neq, (double)n, (double)(neq - n));
229   /* Check analytic over cells */
230   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
231   PetscCall(DMGetCoordinateSection(dm, &coordSection));
232   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
233   neq = 0.0;
234   for (c = cStart; c < cEnd; ++c) {
235     PetscScalar *vcoords = NULL;
236 
237     PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords));
238     neq += ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
239     PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords));
240   }
241   PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell Int f %g != %g mass (%g)", (double)neq, (double)n, (double)(neq - n));
242   /* Check quadrature over entire line */
243   PetscCall(PetscMalloc2(Nq, &xq, Nq, &wq));
244   PetscCall(PetscDTGaussQuadrature(100, vmin, vmax, xq, wq));
245   neq = 0.0;
246   for (q = 0; q < Nq; ++q) neq += ComputePDF(m, n, T, &xq[q]) * wq[q];
247   PetscCheck(PetscAbsReal(neq - n) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int f %g != %g mass (%g)", (double)neq, (double)n, (double)(neq - n));
248   /* Check omemnts with quadrature */
249   veq = 0.0;
250   for (q = 0; q < Nq; ++q) veq += xq[q] * ComputePDF(m, n, T, &xq[q]) * wq[q];
251   veq /= neq;
252   PetscCheck(PetscAbsReal(veq - v[0]) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v f %g != %g velocity (%g)", (double)veq, (double)v[0], (double)(veq - v[0]));
253   Teq = 0.0;
254   for (q = 0; q < Nq; ++q) Teq += PetscSqr(xq[q]) * ComputePDF(m, n, T, &xq[q]) * wq[q];
255   Teq = Teq * m / neq - PetscSqr(veq);
256   PetscCheck(PetscAbsReal(Teq - T) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Int v^2 f %g != %g temperature (%g)", (double)Teq, (double)T, (double)(Teq - T));
257   PetscCall(PetscFree2(xq, wq));
258   PetscFunctionReturn(0);
259 }
260 
261 static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx) {
262   const PetscScalar *u;
263   PetscSection       coordSection;
264   Vec                coordsLocal;
265   PetscScalar       *r, *coords;
266   PetscReal          n = 0.0, v = 0.0, E = 0.0, T = 0.0, m = 1.0, cn = 0.0, cv = 0.0, cE = 0.0, pE = 0.0, eqE = 0.0;
267   PetscInt           dim, d, Np, Ncp, p, cStart, cEnd, c;
268   DM                 dmSw, plex;
269 
270   PetscFunctionBeginUser;
271   PetscCall(VecGetLocalSize(U, &Np));
272   PetscCall(VecGetArrayRead(U, &u));
273   PetscCall(VecGetArray(R, &r));
274   PetscCall(TSGetDM(ts, &dmSw));
275   PetscCall(DMSwarmGetCellDM(dmSw, &plex));
276   PetscCall(DMGetDimension(dmSw, &dim));
277   PetscCall(DMGetCoordinatesLocal(plex, &coordsLocal));
278   PetscCall(DMGetCoordinateSection(plex, &coordSection));
279   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
280   Np /= dim;
281   Ncp = Np / (cEnd - cStart);
282   /* Calculate moments of particle distribution, note that velocity is in the coordinate */
283   PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
284   for (p = 0; p < Np; ++p) {
285     PetscReal m1 = 0.0, m2 = 0.0;
286 
287     for (d = 0; d < dim; ++d) {
288       m1 += PetscRealPart(coords[p * dim + d]);
289       m2 += PetscSqr(PetscRealPart(coords[p * dim + d]));
290     }
291     n += u[p];
292     v += u[p] * m1;
293     E += u[p] * m2;
294   }
295   v /= n;
296   T = E * m / n - v * v;
297   PetscCall(PetscInfo(ts, "Time %.2f: mass %.4f velocity: %+.4f temperature: %.4f\n", (double)t, (double)n, (double)v, (double)T));
298   PetscCall(CheckDistribution(plex, m, n, T, &v));
299   /*
300      Begin cellwise evaluation of the collision operator. Essentially, penalize the weights of the particles
301      in that cell against the maxwellian for the number of particles expected to be in that cell
302   */
303   for (c = cStart; c < cEnd; ++c) {
304     PetscScalar *vcoords    = NULL;
305     PetscReal    relaxation = 1.0, neq;
306     PetscInt     sp         = c * Ncp, q;
307 
308     /* Calculate equilibrium occupation for this velocity cell */
309     PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords));
310     neq = ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
311     PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords));
312     for (q = 0; q < Ncp; ++q) r[sp + q] = (1.0 / relaxation) * (neq - u[sp + q]);
313   }
314   /* Check update */
315   for (p = 0; p < Np; ++p) {
316     PetscReal    m1 = 0.0, m2 = 0.0;
317     PetscScalar *vcoords = NULL;
318 
319     for (d = 0; d < dim; ++d) {
320       m1 += PetscRealPart(coords[p * dim + d]);
321       m2 += PetscSqr(PetscRealPart(coords[p * dim + d]));
322     }
323     cn += r[p];
324     cv += r[p] * m1;
325     cE += r[p] * m2;
326     pE += u[p] * m2;
327     PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords));
328     eqE += ComputeCDF(m, n, T, vcoords[0], vcoords[1]) * m2;
329     PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords));
330   }
331   PetscCall(PetscInfo(ts, "Time %.2f: mass update %.8f velocity update: %+.8f energy update: %.8f (%.8f, %.8f)\n", (double)t, (double)cn, (double)cv, (double)cE, (double)pE, (double)eqE));
332   PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
333   PetscCall(VecRestoreArrayRead(U, &u));
334   PetscCall(VecRestoreArray(R, &r));
335   PetscFunctionReturn(0);
336 }
337 
338 static PetscErrorCode HGMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) {
339   AppCtx            *user = (AppCtx *)ctx;
340   const PetscScalar *u;
341   DM                 sw, dm;
342   PetscInt           dim, Np, p;
343 
344   PetscFunctionBeginUser;
345   if (step < 0) PetscFunctionReturn(0);
346   if (((user->ostep > 0) && (!(step % user->ostep)))) {
347     PetscDrawAxis axis;
348 
349     PetscCall(TSGetDM(ts, &sw));
350     PetscCall(DMSwarmGetCellDM(sw, &dm));
351     PetscCall(DMGetDimension(dm, &dim));
352     PetscCall(PetscDrawHGReset(user->drawhg));
353     PetscCall(PetscDrawHGGetAxis(user->drawhg, &axis));
354     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "f(V)"));
355     PetscCall(PetscDrawAxisSetLimits(axis, -3, 3, 0, 100));
356     PetscCall(PetscDrawHGSetLimits(user->drawhg, -3.0, 3.0, 0, 10));
357     PetscCall(VecGetLocalSize(U, &Np));
358     Np /= dim;
359     PetscCall(VecGetArrayRead(U, &u));
360     /* get points from solution vector */
361     for (p = 0; p < Np; ++p) PetscCall(PetscDrawHGAddValue(user->drawhg, u[p]));
362     PetscCall(PetscDrawHGDraw(user->drawhg));
363     PetscCall(VecRestoreArrayRead(U, &u));
364   }
365   PetscFunctionReturn(0);
366 }
367 
368 static PetscErrorCode SPMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) {
369   AppCtx            *user = (AppCtx *)ctx;
370   const PetscScalar *u;
371   PetscReal         *v, *coords;
372   PetscInt           Np, p;
373   DM                 dmSw;
374 
375   PetscFunctionBeginUser;
376 
377   if (step < 0) PetscFunctionReturn(0);
378   if (((user->ostep > 0) && (!(step % user->ostep)))) {
379     PetscDrawAxis axis;
380 
381     PetscCall(TSGetDM(ts, &dmSw));
382     PetscCall(PetscDrawSPReset(user->drawsp));
383     PetscCall(PetscDrawSPGetAxis(user->drawsp, &axis));
384     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "w"));
385     PetscCall(VecGetLocalSize(U, &Np));
386     PetscCall(VecGetArrayRead(U, &u));
387     /* get points from solution vector */
388     PetscCall(PetscMalloc1(Np, &v));
389     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
390     PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
391     for (p = 0; p < Np - 1; ++p) PetscCall(PetscDrawSPAddPoint(user->drawsp, &coords[p], &v[p]));
392     PetscCall(PetscDrawSPDraw(user->drawsp, PETSC_TRUE));
393     PetscCall(VecRestoreArrayRead(U, &u));
394     PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
395     PetscCall(PetscFree(v));
396   }
397   PetscFunctionReturn(0);
398 }
399 
400 static PetscErrorCode KSConv(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) {
401   AppCtx            *user = (AppCtx *)ctx;
402   const PetscScalar *u;
403   PetscScalar        sup;
404   PetscReal         *v, *coords, T = 0., vel = 0., step_cast, w_sum;
405   PetscInt           dim, Np, p, cStart, cEnd;
406   DM                 sw, plex;
407 
408   PetscFunctionBeginUser;
409   if (step < 0) PetscFunctionReturn(0);
410   if (((user->ostep > 0) && (!(step % user->ostep)))) {
411     PetscDrawAxis axis;
412     PetscCall(PetscDrawSPGetAxis(user->drawks, &axis));
413     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "t", "D_n"));
414     PetscCall(PetscDrawSPSetLimits(user->drawks, 0., 100, 0., 3.5));
415     PetscCall(TSGetDM(ts, &sw));
416     PetscCall(VecGetLocalSize(U, &Np));
417     PetscCall(VecGetArrayRead(U, &u));
418     /* get points from solution vector */
419     PetscCall(PetscMalloc1(Np, &v));
420     PetscCall(DMSwarmGetCellDM(sw, &plex));
421     PetscCall(DMGetDimension(sw, &dim));
422     PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
423     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
424     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
425     w_sum = 0.;
426     for (p = 0; p < Np; ++p) {
427       w_sum += u[p];
428       T += u[p] * coords[p] * coords[p];
429       vel += u[p] * coords[p];
430     }
431     vel /= w_sum;
432     T   = T / w_sum - vel * vel;
433     sup = 0.0;
434     for (p = 0; p < Np; ++p) {
435       PetscReal tmp = 0.;
436       tmp           = PetscAbs(u[p] - ComputePDF(1.0, w_sum, T, &coords[p * dim]));
437       if (tmp > sup) sup = tmp;
438     }
439     step_cast = (PetscReal)step;
440     PetscCall(PetscDrawSPAddPoint(user->drawks, &step_cast, &sup));
441     PetscCall(PetscDrawSPDraw(user->drawks, PETSC_TRUE));
442     PetscCall(VecRestoreArrayRead(U, &u));
443     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
444     PetscCall(PetscFree(v));
445   }
446   PetscFunctionReturn(0);
447 }
448 
449 static PetscErrorCode InitializeSolve(TS ts, Vec u) {
450   DM      dm;
451   AppCtx *user;
452 
453   PetscFunctionBeginUser;
454   PetscCall(TSGetDM(ts, &dm));
455   PetscCall(DMGetApplicationContext(dm, &user));
456   PetscCall(SetInitialCoordinates(dm));
457   PetscCall(SetInitialConditions(dm, u));
458   PetscFunctionReturn(0);
459 }
460 /*
461      A single particle is generated for each velocity space cell using the dmswarmpicfield_coor and is used to evaluate collisions in that cell.
462      0 weight ghost particles are initialized outside of a small velocity domain to ensure the tails of the amxwellian are resolved.
463    */
464 int main(int argc, char **argv) {
465   TS       ts;     /* nonlinear solver */
466   DM       dm, sw; /* Velocity space mesh and Particle Swarm */
467   Vec      u, w;   /* swarm vector */
468   MPI_Comm comm;
469   AppCtx   user;
470 
471   PetscFunctionBeginUser;
472   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
473   comm = PETSC_COMM_WORLD;
474   PetscCall(ProcessOptions(comm, &user));
475   PetscCall(CreateMesh(comm, &dm, &user));
476   PetscCall(CreateParticles(dm, &sw, &user));
477   PetscCall(DMSetApplicationContext(sw, &user));
478   PetscCall(TSCreate(comm, &ts));
479   PetscCall(TSSetDM(ts, sw));
480   PetscCall(TSSetMaxTime(ts, 10.0));
481   PetscCall(TSSetTimeStep(ts, 0.01));
482   PetscCall(TSSetMaxSteps(ts, 100000));
483   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
484   if (user.monitorhg) {
485     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
486     PetscCall(PetscDrawSetFromOptions(user.draw));
487     PetscCall(PetscDrawHGCreate(user.draw, 20, &user.drawhg));
488     PetscCall(PetscDrawHGSetColor(user.drawhg, 3));
489     PetscCall(TSMonitorSet(ts, HGMonitor, &user, NULL));
490   } else if (user.monitorsp) {
491     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
492     PetscCall(PetscDrawSetFromOptions(user.draw));
493     PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawsp));
494     PetscCall(TSMonitorSet(ts, SPMonitor, &user, NULL));
495   } else if (user.monitorks) {
496     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
497     PetscCall(PetscDrawSetFromOptions(user.draw));
498     PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawks));
499     PetscCall(TSMonitorSet(ts, KSConv, &user, NULL));
500   }
501   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user));
502   PetscCall(TSSetFromOptions(ts));
503   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
504   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w));
505   PetscCall(VecDuplicate(w, &u));
506   PetscCall(VecCopy(w, u));
507   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w));
508   PetscCall(TSComputeInitialCondition(ts, u));
509   PetscCall(TSSolve(ts, u));
510   if (user.monitorhg) {
511     PetscCall(PetscDrawSave(user.draw));
512     PetscCall(PetscDrawHGDestroy(&user.drawhg));
513     PetscCall(PetscDrawDestroy(&user.draw));
514   }
515   if (user.monitorsp) {
516     PetscCall(PetscDrawSave(user.draw));
517     PetscCall(PetscDrawSPDestroy(&user.drawsp));
518     PetscCall(PetscDrawDestroy(&user.draw));
519   }
520   if (user.monitorks) {
521     PetscCall(PetscDrawSave(user.draw));
522     PetscCall(PetscDrawSPDestroy(&user.drawks));
523     PetscCall(PetscDrawDestroy(&user.draw));
524   }
525   PetscCall(VecDestroy(&u));
526   PetscCall(TSDestroy(&ts));
527   PetscCall(DMDestroy(&sw));
528   PetscCall(DMDestroy(&dm));
529   PetscCall(PetscFinalize());
530   return 0;
531 }
532 
533 /*TEST
534    build:
535      requires: double !complex
536    test:
537      suffix: 1
538      args: -particles_per_cell 1 -output_step 10 -ts_type euler -dm_plex_dim 1 -dm_plex_box_faces 200 -dm_plex_box_lower -10 -dm_plex_box_upper 10 -dm_view -monitorsp
539    test:
540      suffix: 2
541      args: -particles_per_cell 1 -output_step 50 -ts_type euler -dm_plex_dim 1 -dm_plex_box_faces 200 -dm_plex_box_lower -10 -dm_plex_box_upper 10 -dm_view -monitorks
542 TEST*/
543