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