xref: /petsc/src/ts/tests/ex28.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
1d0c080abSJoseph Pusztay static char help[] = "Example application of the Bhatnagar-Gross-Krook (BGK) collision operator.\n\
2d0c080abSJoseph Pusztay This example is a 0D-1V setting for the kinetic equation\n\
3d0c080abSJoseph Pusztay https://en.wikipedia.org/wiki/Bhatnagar%E2%80%93Gross%E2%80%93Krook_operator\n";
4d0c080abSJoseph Pusztay 
5d0c080abSJoseph Pusztay #include <petscdmplex.h>
6d0c080abSJoseph Pusztay #include <petscdmswarm.h>
7d0c080abSJoseph Pusztay #include <petscts.h>
8d0c080abSJoseph Pusztay #include <petscdraw.h>
9d0c080abSJoseph Pusztay #include <petscviewer.h>
10d0c080abSJoseph Pusztay 
11d0c080abSJoseph Pusztay typedef struct {
12d0c080abSJoseph Pusztay   PetscInt    particlesPerCell; /* The number of partices per cell */
13d0c080abSJoseph Pusztay   PetscReal   momentTol;        /* Tolerance for checking moment conservation */
14d0c080abSJoseph Pusztay   PetscBool   monitorhg;        /* Flag for using the TS histogram monitor */
15d0c080abSJoseph Pusztay   PetscBool   monitorsp;        /* Flag for using the TS scatter monitor */
16d0c080abSJoseph Pusztay   PetscBool   monitorks;        /* Monitor to perform KS test to the maxwellian */
17d0c080abSJoseph Pusztay   PetscBool   error;            /* Flag for printing the error */
18d0c080abSJoseph Pusztay   PetscInt    ostep;            /* print the energy at each ostep time steps */
19d0c080abSJoseph Pusztay   PetscDraw   draw;             /* The draw object for histogram monitoring */
20d0c080abSJoseph Pusztay   PetscDrawHG drawhg;           /* The histogram draw context for monitoring */
21d0c080abSJoseph Pusztay   PetscDrawSP drawsp;           /* The scatter plot draw context for the monitor */
22d0c080abSJoseph Pusztay   PetscDrawSP drawks;           /* Scatterplot draw context for KS test */
23d0c080abSJoseph Pusztay } AppCtx;
24d0c080abSJoseph Pusztay 
25d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
26d71ae5a4SJacob Faibussowitsch {
27d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
28d0c080abSJoseph Pusztay   options->monitorhg        = PETSC_FALSE;
29d0c080abSJoseph Pusztay   options->monitorsp        = PETSC_FALSE;
30d0c080abSJoseph Pusztay   options->monitorks        = PETSC_FALSE;
31d0c080abSJoseph Pusztay   options->particlesPerCell = 1;
32d0c080abSJoseph Pusztay   options->momentTol        = 100.0 * PETSC_MACHINE_EPSILON;
33d0c080abSJoseph Pusztay   options->ostep            = 100;
34d0c080abSJoseph Pusztay 
35d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");
369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitorhg", "Flag to use the TS histogram monitor", "ex28.c", options->monitorhg, &options->monitorhg, NULL));
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitorsp", "Flag to use the TS scatter plot monitor", "ex28.c", options->monitorsp, &options->monitorsp, NULL));
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitorks", "Flag to plot KS test results", "ex28.c", options->monitorks, &options->monitorks, NULL));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex28.c", options->particlesPerCell, &options->particlesPerCell, NULL));
40f3fa974cSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex28.c", options->ostep, &options->ostep, NULL));
41d0609cedSBarry Smith   PetscOptionsEnd();
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43d0c080abSJoseph Pusztay }
44d0c080abSJoseph Pusztay 
45d0c080abSJoseph Pusztay /* Create the mesh for velocity space */
46d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
47d71ae5a4SJacob Faibussowitsch {
48d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
499566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
509566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
519566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
529566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54d0c080abSJoseph Pusztay }
55d0c080abSJoseph Pusztay 
56d0c080abSJoseph Pusztay /* Since we are putting the same number of particles in each cell, this amounts to a uniform distribution of v */
57d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialCoordinates(DM sw)
58d71ae5a4SJacob Faibussowitsch {
59d0c080abSJoseph Pusztay   AppCtx        *user;
60d0c080abSJoseph Pusztay   PetscRandom    rnd;
61d0c080abSJoseph Pusztay   DM             dm;
62d0c080abSJoseph Pusztay   DMPolytopeType ct;
63d0c080abSJoseph Pusztay   PetscBool      simplex;
64d0c080abSJoseph Pusztay   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *vals;
65d0c080abSJoseph Pusztay   PetscInt       dim, d, cStart, cEnd, c, Np, p;
66d0c080abSJoseph Pusztay 
67d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
699566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
709566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
71d0c080abSJoseph Pusztay 
729566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(sw, &user));
73d0c080abSJoseph Pusztay   Np = user->particlesPerCell;
749566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
759566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dm));
769f4ada15SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
79d0c080abSJoseph Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
81d0c080abSJoseph Pusztay   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
829566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
839566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals));
84d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
85d0c080abSJoseph Pusztay     if (Np == 1) {
869566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
87d0c080abSJoseph Pusztay       for (d = 0; d < dim; ++d) {
88d0c080abSJoseph Pusztay         coords[c * dim + d] = centroid[d];
89d0c080abSJoseph Pusztay         if ((coords[c * dim + d] >= -1) && (coords[c * dim + d] <= 1)) {
90d0c080abSJoseph Pusztay           vals[c] = 1.0;
91d0c080abSJoseph Pusztay         } else {
92d0c080abSJoseph Pusztay           vals[c] = 0.;
93d0c080abSJoseph Pusztay         }
94d0c080abSJoseph Pusztay       }
95d0c080abSJoseph Pusztay     } else {
969566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
97d0c080abSJoseph Pusztay       for (p = 0; p < Np; ++p) {
98d0c080abSJoseph Pusztay         const PetscInt n   = c * Np + p;
99d0c080abSJoseph Pusztay         PetscReal      sum = 0.0, refcoords[3];
100d0c080abSJoseph Pusztay 
101d0c080abSJoseph Pusztay         for (d = 0; d < dim; ++d) {
1029566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
103d0c080abSJoseph Pusztay           sum += refcoords[d];
104d0c080abSJoseph Pusztay         }
1059371c9d4SSatish Balay         if (simplex && sum > 0.0)
1069371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
107d0c080abSJoseph Pusztay         vals[n] = 1.0;
1089566063dSJacob Faibussowitsch         PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim]));
109d0c080abSJoseph Pusztay       }
110d0c080abSJoseph Pusztay     }
111d0c080abSJoseph Pusztay   }
1129566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1139566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals));
1149566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
1159566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117d0c080abSJoseph Pusztay }
118d0c080abSJoseph Pusztay 
119da81f932SPierre Jolivet /* The initial conditions are just the initial particle weights */
120d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
121d71ae5a4SJacob Faibussowitsch {
122d0c080abSJoseph Pusztay   DM           dm;
123d0c080abSJoseph Pusztay   AppCtx      *user;
124d0c080abSJoseph Pusztay   PetscReal   *vals;
125d0c080abSJoseph Pusztay   PetscScalar *initialConditions;
126d0c080abSJoseph Pusztay   PetscInt     dim, d, cStart, cEnd, c, Np, p, n;
127d0c080abSJoseph Pusztay 
128d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
1299566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(u, &n));
1309566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dmSw, &user));
131d0c080abSJoseph Pusztay   Np = user->particlesPerCell;
1329566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dmSw, &dm));
1339566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
13563a3b9bcSJacob Faibussowitsch   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);
1369566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dmSw, "w_q", NULL, NULL, (void **)&vals));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(u, &initialConditions));
138d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
139d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
140d0c080abSJoseph Pusztay       const PetscInt n = c * Np + p;
141ad540459SPierre Jolivet       for (d = 0; d < dim; d++) initialConditions[n] = vals[n];
142d0c080abSJoseph Pusztay     }
143d0c080abSJoseph Pusztay   }
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(u, &initialConditions));
1459566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dmSw, "w_q", NULL, NULL, (void **)&vals));
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147d0c080abSJoseph Pusztay }
148d0c080abSJoseph Pusztay 
149d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
150d71ae5a4SJacob Faibussowitsch {
151d0c080abSJoseph Pusztay   PetscInt *cellid;
152d0c080abSJoseph Pusztay   PetscInt  dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
153d0c080abSJoseph Pusztay 
154d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
1559566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1569566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1579566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1589566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
1599566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1609566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
1619566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL));
1629566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1639566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
1649566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1659566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
1669566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
1679566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
168d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
169d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
170d0c080abSJoseph Pusztay       const PetscInt n = c * Np + p;
171d0c080abSJoseph Pusztay       cellid[n]        = c;
172d0c080abSJoseph Pusztay     }
173d0c080abSJoseph Pusztay   }
1749566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
1759566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
1769566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178d0c080abSJoseph Pusztay }
179d0c080abSJoseph Pusztay 
180d0c080abSJoseph Pusztay /*
181d0c080abSJoseph Pusztay   f_t = 1/\tau \left( f_eq - f \right)
182d0c080abSJoseph Pusztay   n_t = 1/\tau \left( \int f_eq - \int f \right) = 1/\tau (n - n) = 0
183d0c080abSJoseph Pusztay   v_t = 1/\tau \left( \int v f_eq - \int v f \right) = 1/\tau (v - v) = 0
184d0c080abSJoseph Pusztay   E_t = 1/\tau \left( \int v^2 f_eq - \int v^2 f \right) = 1/\tau (T - T) = 0
185d0c080abSJoseph Pusztay 
186d0c080abSJoseph Pusztay   Let's look at a single cell:
187d0c080abSJoseph Pusztay 
188d0c080abSJoseph Pusztay     \int_C f_t             = 1/\tau \left( \int_C f_eq - \int_C f \right)
189d0c080abSJoseph Pusztay     \sum_{x_i \in C} w^i_t = 1/\tau \left( F_eq(C) - \sum_{x_i \in C} w_i \right)
190d0c080abSJoseph Pusztay */
191d0c080abSJoseph Pusztay 
192d0c080abSJoseph Pusztay /* This computes the 1D Maxwellian distribution for given mass n, velocity v, and temperature T */
193d71ae5a4SJacob Faibussowitsch static PetscReal ComputePDF(PetscReal m, PetscReal n, PetscReal T, PetscReal v[])
194d71ae5a4SJacob Faibussowitsch {
195d0c080abSJoseph Pusztay   return (n / PetscSqrtReal(2.0 * PETSC_PI * T / m)) * PetscExpReal(-0.5 * m * PetscSqr(v[0]) / T);
196d0c080abSJoseph Pusztay }
197d0c080abSJoseph Pusztay 
198d0c080abSJoseph Pusztay /*
199d0c080abSJoseph Pusztay   erf z = \frac{2}{\sqrt\pi} \int^z_0 dt e^{-t^2} and erf \infty = 1
200d0c080abSJoseph Pusztay 
201d0c080abSJoseph Pusztay   We begin with our distribution
202d0c080abSJoseph Pusztay 
203d0c080abSJoseph Pusztay     \sqrt{\frac{m}{2 \pi T}} e^{-m v^2/2T}
204d0c080abSJoseph Pusztay 
205d0c080abSJoseph Pusztay   Let t = \sqrt{m/2T} v, z = \sqrt{m/2T} w, so that we now have
206d0c080abSJoseph Pusztay 
207d0c080abSJoseph Pusztay       \sqrt{\frac{m}{2 \pi T}} \int^w_0 dv e^{-m v^2/2T}
208d0c080abSJoseph Pusztay     = \sqrt{\frac{m}{2 \pi T}} \int^{\sqrt{m/2T} w}_0 \sqrt{2T/m} dt e^{-t^2}
209d0c080abSJoseph Pusztay     = 1/\sqrt{\pi} \int^{\sqrt{m/2T} w}_0 dt e^{-t^2}
210d0c080abSJoseph Pusztay     = 1/2 erf(\sqrt{m/2T} w)
211d0c080abSJoseph Pusztay */
212d71ae5a4SJacob Faibussowitsch static PetscReal ComputeCDF(PetscReal m, PetscReal n, PetscReal T, PetscReal va, PetscReal vb)
213d71ae5a4SJacob Faibussowitsch {
214d0c080abSJoseph Pusztay   PetscReal alpha = PetscSqrtReal(0.5 * m / T);
215d0c080abSJoseph Pusztay   PetscReal za    = alpha * va;
216d0c080abSJoseph Pusztay   PetscReal zb    = alpha * vb;
217d0c080abSJoseph Pusztay   PetscReal sum   = 0.0;
218d0c080abSJoseph Pusztay 
219d0c080abSJoseph Pusztay   sum += zb >= 0 ? erf(zb) : -erf(-zb);
220d0c080abSJoseph Pusztay   sum -= za >= 0 ? erf(za) : -erf(-za);
221d0c080abSJoseph Pusztay   return 0.5 * n * sum;
222d0c080abSJoseph Pusztay }
223d0c080abSJoseph Pusztay 
224d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckDistribution(DM dm, PetscReal m, PetscReal n, PetscReal T, PetscReal v[])
225d71ae5a4SJacob Faibussowitsch {
226d0c080abSJoseph Pusztay   PetscSection coordSection;
227d0c080abSJoseph Pusztay   Vec          coordsLocal;
228d0c080abSJoseph Pusztay   PetscReal   *xq, *wq;
229d0c080abSJoseph Pusztay   PetscReal    vmin, vmax, neq, veq, Teq;
230d0c080abSJoseph Pusztay   PetscInt     Nq = 100, q, cStart, cEnd, c;
231d0c080abSJoseph Pusztay 
2327510d9b0SBarry Smith   PetscFunctionBeginUser;
2339566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dm, &vmin, &vmax));
234d0c080abSJoseph Pusztay   /* Check analytic over entire line */
235d0c080abSJoseph Pusztay   neq = ComputeCDF(m, n, T, vmin, vmax);
23663a3b9bcSJacob Faibussowitsch   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));
237d0c080abSJoseph Pusztay   /* Check analytic over cells */
2389566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
2399566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
2409566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
241d0c080abSJoseph Pusztay   neq = 0.0;
242d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
243d0c080abSJoseph Pusztay     PetscScalar *vcoords = NULL;
244d0c080abSJoseph Pusztay 
2459566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords));
246d0c080abSJoseph Pusztay     neq += ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
2479566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords));
248d0c080abSJoseph Pusztay   }
24963a3b9bcSJacob Faibussowitsch   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));
250d0c080abSJoseph Pusztay   /* Check quadrature over entire line */
2519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nq, &xq, Nq, &wq));
2529566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussQuadrature(100, vmin, vmax, xq, wq));
253d0c080abSJoseph Pusztay   neq = 0.0;
254ad540459SPierre Jolivet   for (q = 0; q < Nq; ++q) neq += ComputePDF(m, n, T, &xq[q]) * wq[q];
25563a3b9bcSJacob Faibussowitsch   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));
256d0c080abSJoseph Pusztay   /* Check omemnts with quadrature */
257d0c080abSJoseph Pusztay   veq = 0.0;
258ad540459SPierre Jolivet   for (q = 0; q < Nq; ++q) veq += xq[q] * ComputePDF(m, n, T, &xq[q]) * wq[q];
259d0c080abSJoseph Pusztay   veq /= neq;
26063a3b9bcSJacob Faibussowitsch   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]));
261d0c080abSJoseph Pusztay   Teq = 0.0;
262ad540459SPierre Jolivet   for (q = 0; q < Nq; ++q) Teq += PetscSqr(xq[q]) * ComputePDF(m, n, T, &xq[q]) * wq[q];
263d0c080abSJoseph Pusztay   Teq = Teq * m / neq - PetscSqr(veq);
26463a3b9bcSJacob Faibussowitsch   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));
2659566063dSJacob Faibussowitsch   PetscCall(PetscFree2(xq, wq));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267d0c080abSJoseph Pusztay }
268d0c080abSJoseph Pusztay 
269d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
270d71ae5a4SJacob Faibussowitsch {
271d0c080abSJoseph Pusztay   const PetscScalar *u;
272d0c080abSJoseph Pusztay   PetscSection       coordSection;
273d0c080abSJoseph Pusztay   Vec                coordsLocal;
274d0c080abSJoseph Pusztay   PetscScalar       *r, *coords;
275d0c080abSJoseph Pusztay   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;
276d0c080abSJoseph Pusztay   PetscInt           dim, d, Np, Ncp, p, cStart, cEnd, c;
277d0c080abSJoseph Pusztay   DM                 dmSw, plex;
278d0c080abSJoseph Pusztay 
279d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
2809566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
2819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2829566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
2839566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dmSw));
2849566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(dmSw, &plex));
2859566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dmSw, &dim));
2869566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(plex, &coordsLocal));
2879566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(plex, &coordSection));
2889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
289d0c080abSJoseph Pusztay   Np /= dim;
290d0c080abSJoseph Pusztay   Ncp = Np / (cEnd - cStart);
291d0c080abSJoseph Pusztay   /* Calculate moments of particle distribution, note that velocity is in the coordinate */
2929566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
293d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
294d0c080abSJoseph Pusztay     PetscReal m1 = 0.0, m2 = 0.0;
295d0c080abSJoseph Pusztay 
2969371c9d4SSatish Balay     for (d = 0; d < dim; ++d) {
2979371c9d4SSatish Balay       m1 += PetscRealPart(coords[p * dim + d]);
2989371c9d4SSatish Balay       m2 += PetscSqr(PetscRealPart(coords[p * dim + d]));
2999371c9d4SSatish Balay     }
300d0c080abSJoseph Pusztay     n += u[p];
301d0c080abSJoseph Pusztay     v += u[p] * m1;
302d0c080abSJoseph Pusztay     E += u[p] * m2;
303d0c080abSJoseph Pusztay   }
304d0c080abSJoseph Pusztay   v /= n;
305d0c080abSJoseph Pusztay   T = E * m / n - v * v;
30663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(ts, "Time %.2f: mass %.4f velocity: %+.4f temperature: %.4f\n", (double)t, (double)n, (double)v, (double)T));
3079566063dSJacob Faibussowitsch   PetscCall(CheckDistribution(plex, m, n, T, &v));
308d0c080abSJoseph Pusztay   /*
309d0c080abSJoseph Pusztay      Begin cellwise evaluation of the collision operator. Essentially, penalize the weights of the particles
310d0c080abSJoseph Pusztay      in that cell against the maxwellian for the number of particles expected to be in that cell
311d0c080abSJoseph Pusztay   */
312d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
313d0c080abSJoseph Pusztay     PetscScalar *vcoords    = NULL;
314d0c080abSJoseph Pusztay     PetscReal    relaxation = 1.0, neq;
315d0c080abSJoseph Pusztay     PetscInt     sp         = c * Ncp, q;
316d0c080abSJoseph Pusztay 
317d0c080abSJoseph Pusztay     /* Calculate equilibrium occupation for this velocity cell */
3189566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords));
319d0c080abSJoseph Pusztay     neq = ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
3209566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords));
321d0c080abSJoseph Pusztay     for (q = 0; q < Ncp; ++q) r[sp + q] = (1.0 / relaxation) * (neq - u[sp + q]);
322d0c080abSJoseph Pusztay   }
323d0c080abSJoseph Pusztay   /* Check update */
324d0c080abSJoseph Pusztay   for (p = 0; p < Np; ++p) {
325d0c080abSJoseph Pusztay     PetscReal    m1 = 0.0, m2 = 0.0;
326d0c080abSJoseph Pusztay     PetscScalar *vcoords = NULL;
327d0c080abSJoseph Pusztay 
3289371c9d4SSatish Balay     for (d = 0; d < dim; ++d) {
3299371c9d4SSatish Balay       m1 += PetscRealPart(coords[p * dim + d]);
3309371c9d4SSatish Balay       m2 += PetscSqr(PetscRealPart(coords[p * dim + d]));
3319371c9d4SSatish Balay     }
332d0c080abSJoseph Pusztay     cn += r[p];
333d0c080abSJoseph Pusztay     cv += r[p] * m1;
334d0c080abSJoseph Pusztay     cE += r[p] * m2;
335d0c080abSJoseph Pusztay     pE += u[p] * m2;
3369566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords));
337d0c080abSJoseph Pusztay     eqE += ComputeCDF(m, n, T, vcoords[0], vcoords[1]) * m2;
3389566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords));
339d0c080abSJoseph Pusztay   }
34063a3b9bcSJacob Faibussowitsch   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));
3419566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
3429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
3439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
345d0c080abSJoseph Pusztay }
346d0c080abSJoseph Pusztay 
347d71ae5a4SJacob Faibussowitsch static PetscErrorCode HGMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
348d71ae5a4SJacob Faibussowitsch {
349d0c080abSJoseph Pusztay   AppCtx            *user = (AppCtx *)ctx;
350d0c080abSJoseph Pusztay   const PetscScalar *u;
351d0c080abSJoseph Pusztay   DM                 sw, dm;
352d0c080abSJoseph Pusztay   PetscInt           dim, Np, p;
353d0c080abSJoseph Pusztay 
354d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
3553ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
356*f4f49eeaSPierre Jolivet   if ((user->ostep > 0) && (!(step % user->ostep))) {
357d0c080abSJoseph Pusztay     PetscDrawAxis axis;
358d0c080abSJoseph Pusztay 
3599566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &sw));
3609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(sw, &dm));
3619566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
3629566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGReset(user->drawhg));
3639566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGGetAxis(user->drawhg, &axis));
3649566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "f(V)"));
3659566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLimits(axis, -3, 3, 0, 100));
3669566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGSetLimits(user->drawhg, -3.0, 3.0, 0, 10));
3679566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(U, &Np));
368d0c080abSJoseph Pusztay     Np /= dim;
3699566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(U, &u));
370d0c080abSJoseph Pusztay     /* get points from solution vector */
3719566063dSJacob Faibussowitsch     for (p = 0; p < Np; ++p) PetscCall(PetscDrawHGAddValue(user->drawhg, u[p]));
3729566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGDraw(user->drawhg));
3739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U, &u));
374d0c080abSJoseph Pusztay   }
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
376d0c080abSJoseph Pusztay }
377d0c080abSJoseph Pusztay 
378d71ae5a4SJacob Faibussowitsch static PetscErrorCode SPMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
379d71ae5a4SJacob Faibussowitsch {
380d0c080abSJoseph Pusztay   AppCtx            *user = (AppCtx *)ctx;
381d0c080abSJoseph Pusztay   const PetscScalar *u;
382d0c080abSJoseph Pusztay   PetscReal         *v, *coords;
383d0c080abSJoseph Pusztay   PetscInt           Np, p;
384d0c080abSJoseph Pusztay   DM                 dmSw;
385d0c080abSJoseph Pusztay 
386d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
3873ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
388*f4f49eeaSPierre Jolivet   if ((user->ostep > 0) && (!(step % user->ostep))) {
389d0c080abSJoseph Pusztay     PetscDrawAxis axis;
390d0c080abSJoseph Pusztay 
3919566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dmSw));
3929566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPReset(user->drawsp));
3939566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(user->drawsp, &axis));
3949566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "w"));
3959566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(U, &Np));
3969566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(U, &u));
397d0c080abSJoseph Pusztay     /* get points from solution vector */
3989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Np, &v));
399d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
4009566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
4019566063dSJacob Faibussowitsch     for (p = 0; p < Np - 1; ++p) PetscCall(PetscDrawSPAddPoint(user->drawsp, &coords[p], &v[p]));
4029566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(user->drawsp, PETSC_TRUE));
4039566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U, &u));
4049566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
4059566063dSJacob Faibussowitsch     PetscCall(PetscFree(v));
406d0c080abSJoseph Pusztay   }
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408d0c080abSJoseph Pusztay }
409d0c080abSJoseph Pusztay 
410d71ae5a4SJacob Faibussowitsch static PetscErrorCode KSConv(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
411d71ae5a4SJacob Faibussowitsch {
412d0c080abSJoseph Pusztay   AppCtx            *user = (AppCtx *)ctx;
413d0c080abSJoseph Pusztay   const PetscScalar *u;
414d0c080abSJoseph Pusztay   PetscScalar        sup;
415d0c080abSJoseph Pusztay   PetscReal         *v, *coords, T = 0., vel = 0., step_cast, w_sum;
416d0c080abSJoseph Pusztay   PetscInt           dim, Np, p, cStart, cEnd;
417d0c080abSJoseph Pusztay   DM                 sw, plex;
418d0c080abSJoseph Pusztay 
419d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
4203ba16761SJacob Faibussowitsch   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
421*f4f49eeaSPierre Jolivet   if ((user->ostep > 0) && (!(step % user->ostep))) {
422d0c080abSJoseph Pusztay     PetscDrawAxis axis;
4239566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPGetAxis(user->drawks, &axis));
4249566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "t", "D_n"));
4259566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPSetLimits(user->drawks, 0., 100, 0., 3.5));
4269566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &sw));
4279566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(U, &Np));
4289566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(U, &u));
429d0c080abSJoseph Pusztay     /* get points from solution vector */
4309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Np, &v));
4319566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetCellDM(sw, &plex));
4329566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(sw, &dim));
4339566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
434d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
4359566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
436d0c080abSJoseph Pusztay     w_sum = 0.;
437d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
438d0c080abSJoseph Pusztay       w_sum += u[p];
439d0c080abSJoseph Pusztay       T += u[p] * coords[p] * coords[p];
440d0c080abSJoseph Pusztay       vel += u[p] * coords[p];
441d0c080abSJoseph Pusztay     }
442d0c080abSJoseph Pusztay     vel /= w_sum;
443d0c080abSJoseph Pusztay     T   = T / w_sum - vel * vel;
444d0c080abSJoseph Pusztay     sup = 0.0;
445d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
446d0c080abSJoseph Pusztay       PetscReal tmp = 0.;
447d0c080abSJoseph Pusztay       tmp           = PetscAbs(u[p] - ComputePDF(1.0, w_sum, T, &coords[p * dim]));
448d0c080abSJoseph Pusztay       if (tmp > sup) sup = tmp;
449d0c080abSJoseph Pusztay     }
450d0c080abSJoseph Pusztay     step_cast = (PetscReal)step;
4519566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPAddPoint(user->drawks, &step_cast, &sup));
4529566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDraw(user->drawks, PETSC_TRUE));
4539566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U, &u));
4549566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
4559566063dSJacob Faibussowitsch     PetscCall(PetscFree(v));
456d0c080abSJoseph Pusztay   }
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458d0c080abSJoseph Pusztay }
459d0c080abSJoseph Pusztay 
460d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitializeSolve(TS ts, Vec u)
461d71ae5a4SJacob Faibussowitsch {
462d0c080abSJoseph Pusztay   DM      dm;
463d0c080abSJoseph Pusztay   AppCtx *user;
464d0c080abSJoseph Pusztay 
465d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
4669566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
4679566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &user));
4689566063dSJacob Faibussowitsch   PetscCall(SetInitialCoordinates(dm));
4699566063dSJacob Faibussowitsch   PetscCall(SetInitialConditions(dm, u));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471d0c080abSJoseph Pusztay }
472d0c080abSJoseph Pusztay /*
473d0c080abSJoseph Pusztay      A single particle is generated for each velocity space cell using the dmswarmpicfield_coor and is used to evaluate collisions in that cell.
474d0c080abSJoseph Pusztay      0 weight ghost particles are initialized outside of a small velocity domain to ensure the tails of the amxwellian are resolved.
475d0c080abSJoseph Pusztay    */
476d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
477d71ae5a4SJacob Faibussowitsch {
478d0c080abSJoseph Pusztay   TS       ts;     /* nonlinear solver */
479d0c080abSJoseph Pusztay   DM       dm, sw; /* Velocity space mesh and Particle Swarm */
480d0c080abSJoseph Pusztay   Vec      u, w;   /* swarm vector */
481d0c080abSJoseph Pusztay   MPI_Comm comm;
482d0c080abSJoseph Pusztay   AppCtx   user;
483d0c080abSJoseph Pusztay 
484327415f7SBarry Smith   PetscFunctionBeginUser;
4859566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
486d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
4879566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
4889566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
4899566063dSJacob Faibussowitsch   PetscCall(CreateParticles(dm, &sw, &user));
4909566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(sw, &user));
4919566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
4929566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, sw));
4939566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10.0));
4949566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.01));
4959566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 100000));
4969566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
497d0c080abSJoseph Pusztay   if (user.monitorhg) {
4989566063dSJacob Faibussowitsch     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
4999566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetFromOptions(user.draw));
5009566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGCreate(user.draw, 20, &user.drawhg));
5019566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGSetColor(user.drawhg, 3));
5029566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, HGMonitor, &user, NULL));
5039371c9d4SSatish Balay   } else if (user.monitorsp) {
5049566063dSJacob Faibussowitsch     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
5059566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetFromOptions(user.draw));
5069566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawsp));
5079566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, SPMonitor, &user, NULL));
5089371c9d4SSatish Balay   } else if (user.monitorks) {
5099566063dSJacob Faibussowitsch     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
5109566063dSJacob Faibussowitsch     PetscCall(PetscDrawSetFromOptions(user.draw));
5119566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawks));
5129566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts, KSConv, &user, NULL));
513d0c080abSJoseph Pusztay   }
5149566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user));
5159566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
5169566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
5179566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w));
5189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(w, &u));
5199566063dSJacob Faibussowitsch   PetscCall(VecCopy(w, u));
5209566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w));
5219566063dSJacob Faibussowitsch   PetscCall(TSComputeInitialCondition(ts, u));
5229566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
523d0c080abSJoseph Pusztay   if (user.monitorhg) {
5249566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(user.draw));
5259566063dSJacob Faibussowitsch     PetscCall(PetscDrawHGDestroy(&user.drawhg));
5269566063dSJacob Faibussowitsch     PetscCall(PetscDrawDestroy(&user.draw));
527d0c080abSJoseph Pusztay   }
528d0c080abSJoseph Pusztay   if (user.monitorsp) {
5299566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(user.draw));
5309566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDestroy(&user.drawsp));
5319566063dSJacob Faibussowitsch     PetscCall(PetscDrawDestroy(&user.draw));
532d0c080abSJoseph Pusztay   }
533d0c080abSJoseph Pusztay   if (user.monitorks) {
5349566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(user.draw));
5359566063dSJacob Faibussowitsch     PetscCall(PetscDrawSPDestroy(&user.drawks));
5369566063dSJacob Faibussowitsch     PetscCall(PetscDrawDestroy(&user.draw));
537d0c080abSJoseph Pusztay   }
5389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
5399566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
5409566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
5419566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
5429566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
543b122ec5aSJacob Faibussowitsch   return 0;
544d0c080abSJoseph Pusztay }
545d0c080abSJoseph Pusztay 
546d0c080abSJoseph Pusztay /*TEST
547d0c080abSJoseph Pusztay    build:
54830602db0SMatthew G. Knepley      requires: double !complex
549d0c080abSJoseph Pusztay    test:
550d0c080abSJoseph Pusztay      suffix: 1
55130602db0SMatthew G. Knepley      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
552d0c080abSJoseph Pusztay    test:
553d0c080abSJoseph Pusztay      suffix: 2
55430602db0SMatthew G. Knepley      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
555d0c080abSJoseph Pusztay TEST*/
556