Actual source code: ex28.c

  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";

  5: #include <petscdmplex.h>
  6: #include <petscdmswarm.h>
  7: #include <petscts.h>
  8: #include <petscdraw.h>
  9: #include <petscviewer.h>

 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;

 25: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 26: {
 27:   PetscFunctionBeginUser;
 28:   options->monitorhg        = PETSC_FALSE;
 29:   options->monitorsp        = PETSC_FALSE;
 30:   options->monitorks        = PETSC_FALSE;
 31:   options->particlesPerCell = 1;
 32:   options->momentTol        = 100.0 * PETSC_MACHINE_EPSILON;
 33:   options->ostep            = 100;

 35:   PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");
 36:   PetscCall(PetscOptionsBool("-monitorhg", "Flag to use the TS histogram monitor", "ex28.c", options->monitorhg, &options->monitorhg, NULL));
 37:   PetscCall(PetscOptionsBool("-monitorsp", "Flag to use the TS scatter plot monitor", "ex28.c", options->monitorsp, &options->monitorsp, NULL));
 38:   PetscCall(PetscOptionsBool("-monitorks", "Flag to plot KS test results", "ex28.c", options->monitorks, &options->monitorks, NULL));
 39:   PetscCall(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex28.c", options->particlesPerCell, &options->particlesPerCell, NULL));
 40:   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex28.c", options->ostep, &options->ostep, NULL));
 41:   PetscOptionsEnd();
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: /* Create the mesh for velocity space */
 46: static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
 47: {
 48:   PetscFunctionBeginUser;
 49:   PetscCall(DMCreate(comm, dm));
 50:   PetscCall(DMSetType(*dm, DMPLEX));
 51:   PetscCall(DMSetFromOptions(*dm));
 52:   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
 53:   PetscFunctionReturn(PETSC_SUCCESS);
 54: }

 56: /* Since we are putting the same number of particles in each cell, this amounts to a uniform distribution of v */
 57: static PetscErrorCode SetInitialCoordinates(DM sw)
 58: {
 59:   AppCtx        *user;
 60:   PetscRandom    rnd;
 61:   DM             dm;
 62:   DMPolytopeType ct;
 63:   PetscBool      simplex;
 64:   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *vals;
 65:   PetscInt       dim, d, cStart, cEnd, c, Np, p;

 67:   PetscFunctionBeginUser;
 68:   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd));
 69:   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
 70:   PetscCall(PetscRandomSetFromOptions(rnd));

 72:   PetscCall(DMGetApplicationContext(sw, &user));
 73:   Np = user->particlesPerCell;
 74:   PetscCall(DMGetDimension(sw, &dim));
 75:   PetscCall(DMSwarmGetCellDM(sw, &dm));
 76:   PetscCall(DMGetCoordinatesLocalSetUp(dm));
 77:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
 78:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
 79:   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
 80:   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
 81:   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
 82:   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
 83:   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&vals));
 84:   for (c = cStart; c < cEnd; ++c) {
 85:     if (Np == 1) {
 86:       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
 87:       for (d = 0; d < dim; ++d) {
 88:         coords[c * dim + d] = centroid[d];
 89:         if ((coords[c * dim + d] >= -1) && (coords[c * dim + d] <= 1)) {
 90:           vals[c] = 1.0;
 91:         } else {
 92:           vals[c] = 0.;
 93:         }
 94:       }
 95:     } else {
 96:       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
 97:       for (p = 0; p < Np; ++p) {
 98:         const PetscInt n   = c * Np + p;
 99:         PetscReal      sum = 0.0, refcoords[3];

101:         for (d = 0; d < dim; ++d) {
102:           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
103:           sum += refcoords[d];
104:         }
105:         if (simplex && sum > 0.0)
106:           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
107:         vals[n] = 1.0;
108:         PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n * dim]));
109:       }
110:     }
111:   }
112:   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
113:   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&vals));
114:   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
115:   PetscCall(PetscRandomDestroy(&rnd));
116:   PetscFunctionReturn(PETSC_SUCCESS);
117: }

119: /* The initial conditions are just the initial particle weights */
120: static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
121: {
122:   DM           dm;
123:   AppCtx      *user;
124:   PetscReal   *vals;
125:   PetscScalar *initialConditions;
126:   PetscInt     dim, d, cStart, cEnd, c, Np, p, n;

128:   PetscFunctionBeginUser;
129:   PetscCall(VecGetLocalSize(u, &n));
130:   PetscCall(DMGetApplicationContext(dmSw, &user));
131:   Np = user->particlesPerCell;
132:   PetscCall(DMSwarmGetCellDM(dmSw, &dm));
133:   PetscCall(DMGetDimension(dm, &dim));
134:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
135:   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);
136:   PetscCall(DMSwarmGetField(dmSw, "w_q", NULL, NULL, (void **)&vals));
137:   PetscCall(VecGetArray(u, &initialConditions));
138:   for (c = cStart; c < cEnd; ++c) {
139:     for (p = 0; p < Np; ++p) {
140:       const PetscInt n = c * Np + p;
141:       for (d = 0; d < dim; d++) initialConditions[n] = vals[n];
142:     }
143:   }
144:   PetscCall(VecRestoreArray(u, &initialConditions));
145:   PetscCall(DMSwarmRestoreField(dmSw, "w_q", NULL, NULL, (void **)&vals));
146:   PetscFunctionReturn(PETSC_SUCCESS);
147: }

149: static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
150: {
151:   PetscInt *cellid;
152:   PetscInt  dim, cStart, cEnd, c, Np = user->particlesPerCell, p;

154:   PetscFunctionBeginUser;
155:   PetscCall(DMGetDimension(dm, &dim));
156:   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
157:   PetscCall(DMSetType(*sw, DMSWARM));
158:   PetscCall(DMSetDimension(*sw, dim));
159:   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
160:   PetscCall(DMSwarmSetCellDM(*sw, dm));
161:   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL));
162:   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
163:   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
164:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
165:   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
166:   PetscCall(DMSetFromOptions(*sw));
167:   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
168:   for (c = cStart; c < cEnd; ++c) {
169:     for (p = 0; p < Np; ++p) {
170:       const PetscInt n = c * Np + p;
171:       cellid[n]        = c;
172:     }
173:   }
174:   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
175:   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
176:   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
177:   PetscFunctionReturn(PETSC_SUCCESS);
178: }

180: /*
181:   f_t = 1/\tau \left( f_eq - f \right)
182:   n_t = 1/\tau \left( \int f_eq - \int f \right) = 1/\tau (n - n) = 0
183:   v_t = 1/\tau \left( \int v f_eq - \int v f \right) = 1/\tau (v - v) = 0
184:   E_t = 1/\tau \left( \int v^2 f_eq - \int v^2 f \right) = 1/\tau (T - T) = 0

186:   Let's look at a single cell:

188:     \int_C f_t             = 1/\tau \left( \int_C f_eq - \int_C f \right)
189:     \sum_{x_i \in C} w^i_t = 1/\tau \left( F_eq(C) - \sum_{x_i \in C} w_i \right)
190: */

192: /* This computes the 1D Maxwellian distribution for given mass n, velocity v, and temperature T */
193: static PetscReal ComputePDF(PetscReal m, PetscReal n, PetscReal T, PetscReal v[])
194: {
195:   return (n / PetscSqrtReal(2.0 * PETSC_PI * T / m)) * PetscExpReal(-0.5 * m * PetscSqr(v[0]) / T);
196: }

198: /*
199:   erf z = \frac{2}{\sqrt\pi} \int^z_0 dt e^{-t^2} and erf \infty = 1

201:   We begin with our distribution

203:     \sqrt{\frac{m}{2 \pi T}} e^{-m v^2/2T}

205:   Let t = \sqrt{m/2T} v, z = \sqrt{m/2T} w, so that we now have

207:       \sqrt{\frac{m}{2 \pi T}} \int^w_0 dv e^{-m v^2/2T}
208:     = \sqrt{\frac{m}{2 \pi T}} \int^{\sqrt{m/2T} w}_0 \sqrt{2T/m} dt e^{-t^2}
209:     = 1/\sqrt{\pi} \int^{\sqrt{m/2T} w}_0 dt e^{-t^2}
210:     = 1/2 erf(\sqrt{m/2T} w)
211: */
212: static PetscReal ComputeCDF(PetscReal m, PetscReal n, PetscReal T, PetscReal va, PetscReal vb)
213: {
214:   PetscReal alpha = PetscSqrtReal(0.5 * m / T);
215:   PetscReal za    = alpha * va;
216:   PetscReal zb    = alpha * vb;
217:   PetscReal sum   = 0.0;

219:   sum += zb >= 0 ? erf(zb) : -erf(-zb);
220:   sum -= za >= 0 ? erf(za) : -erf(-za);
221:   return 0.5 * n * sum;
222: }

224: static PetscErrorCode CheckDistribution(DM dm, PetscReal m, PetscReal n, PetscReal T, PetscReal v[])
225: {
226:   PetscSection coordSection;
227:   Vec          coordsLocal;
228:   PetscReal   *xq, *wq;
229:   PetscReal    vmin, vmax, neq, veq, Teq;
230:   PetscInt     Nq = 100, q, cStart, cEnd, c;

232:   PetscFunctionBeginUser;
233:   PetscCall(DMGetBoundingBox(dm, &vmin, &vmax));
234:   /* Check analytic over entire line */
235:   neq = ComputeCDF(m, n, T, vmin, vmax);
236:   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));
237:   /* Check analytic over cells */
238:   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
239:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
240:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
241:   neq = 0.0;
242:   for (c = cStart; c < cEnd; ++c) {
243:     PetscScalar *vcoords = NULL;

245:     PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords));
246:     neq += ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
247:     PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, NULL, &vcoords));
248:   }
249:   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));
250:   /* Check quadrature over entire line */
251:   PetscCall(PetscMalloc2(Nq, &xq, Nq, &wq));
252:   PetscCall(PetscDTGaussQuadrature(100, vmin, vmax, xq, wq));
253:   neq = 0.0;
254:   for (q = 0; q < Nq; ++q) neq += ComputePDF(m, n, T, &xq[q]) * wq[q];
255:   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));
256:   /* Check omemnts with quadrature */
257:   veq = 0.0;
258:   for (q = 0; q < Nq; ++q) veq += xq[q] * ComputePDF(m, n, T, &xq[q]) * wq[q];
259:   veq /= neq;
260:   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]));
261:   Teq = 0.0;
262:   for (q = 0; q < Nq; ++q) Teq += PetscSqr(xq[q]) * ComputePDF(m, n, T, &xq[q]) * wq[q];
263:   Teq = Teq * m / neq - PetscSqr(veq);
264:   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));
265:   PetscCall(PetscFree2(xq, wq));
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
270: {
271:   const PetscScalar *u;
272:   PetscSection       coordSection;
273:   Vec                coordsLocal;
274:   PetscScalar       *r, *coords;
275:   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;
276:   PetscInt           dim, d, Np, Ncp, p, cStart, cEnd, c;
277:   DM                 dmSw, plex;

279:   PetscFunctionBeginUser;
280:   PetscCall(VecGetLocalSize(U, &Np));
281:   PetscCall(VecGetArrayRead(U, &u));
282:   PetscCall(VecGetArray(R, &r));
283:   PetscCall(TSGetDM(ts, &dmSw));
284:   PetscCall(DMSwarmGetCellDM(dmSw, &plex));
285:   PetscCall(DMGetDimension(dmSw, &dim));
286:   PetscCall(DMGetCoordinatesLocal(plex, &coordsLocal));
287:   PetscCall(DMGetCoordinateSection(plex, &coordSection));
288:   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
289:   Np /= dim;
290:   Ncp = Np / (cEnd - cStart);
291:   /* Calculate moments of particle distribution, note that velocity is in the coordinate */
292:   PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
293:   for (p = 0; p < Np; ++p) {
294:     PetscReal m1 = 0.0, m2 = 0.0;

296:     for (d = 0; d < dim; ++d) {
297:       m1 += PetscRealPart(coords[p * dim + d]);
298:       m2 += PetscSqr(PetscRealPart(coords[p * dim + d]));
299:     }
300:     n += u[p];
301:     v += u[p] * m1;
302:     E += u[p] * m2;
303:   }
304:   v /= n;
305:   T = E * m / n - v * v;
306:   PetscCall(PetscInfo(ts, "Time %.2f: mass %.4f velocity: %+.4f temperature: %.4f\n", (double)t, (double)n, (double)v, (double)T));
307:   PetscCall(CheckDistribution(plex, m, n, T, &v));
308:   /*
309:      Begin cellwise evaluation of the collision operator. Essentially, penalize the weights of the particles
310:      in that cell against the maxwellian for the number of particles expected to be in that cell
311:   */
312:   for (c = cStart; c < cEnd; ++c) {
313:     PetscScalar *vcoords    = NULL;
314:     PetscReal    relaxation = 1.0, neq;
315:     PetscInt     sp         = c * Ncp, q;

317:     /* Calculate equilibrium occupation for this velocity cell */
318:     PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords));
319:     neq = ComputeCDF(m, n, T, vcoords[0], vcoords[1]);
320:     PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, c, NULL, &vcoords));
321:     for (q = 0; q < Ncp; ++q) r[sp + q] = (1.0 / relaxation) * (neq - u[sp + q]);
322:   }
323:   /* Check update */
324:   for (p = 0; p < Np; ++p) {
325:     PetscReal    m1 = 0.0, m2 = 0.0;
326:     PetscScalar *vcoords = NULL;

328:     for (d = 0; d < dim; ++d) {
329:       m1 += PetscRealPart(coords[p * dim + d]);
330:       m2 += PetscSqr(PetscRealPart(coords[p * dim + d]));
331:     }
332:     cn += r[p];
333:     cv += r[p] * m1;
334:     cE += r[p] * m2;
335:     pE += u[p] * m2;
336:     PetscCall(DMPlexVecGetClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords));
337:     eqE += ComputeCDF(m, n, T, vcoords[0], vcoords[1]) * m2;
338:     PetscCall(DMPlexVecRestoreClosure(plex, coordSection, coordsLocal, p, NULL, &vcoords));
339:   }
340:   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));
341:   PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
342:   PetscCall(VecRestoreArrayRead(U, &u));
343:   PetscCall(VecRestoreArray(R, &r));
344:   PetscFunctionReturn(PETSC_SUCCESS);
345: }

347: static PetscErrorCode HGMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
348: {
349:   AppCtx            *user = (AppCtx *)ctx;
350:   const PetscScalar *u;
351:   DM                 sw, dm;
352:   PetscInt           dim, Np, p;

354:   PetscFunctionBeginUser;
355:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
356:   if ((user->ostep > 0) && (!(step % user->ostep))) {
357:     PetscDrawAxis axis;

359:     PetscCall(TSGetDM(ts, &sw));
360:     PetscCall(DMSwarmGetCellDM(sw, &dm));
361:     PetscCall(DMGetDimension(dm, &dim));
362:     PetscCall(PetscDrawHGReset(user->drawhg));
363:     PetscCall(PetscDrawHGGetAxis(user->drawhg, &axis));
364:     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "f(V)"));
365:     PetscCall(PetscDrawAxisSetLimits(axis, -3, 3, 0, 100));
366:     PetscCall(PetscDrawHGSetLimits(user->drawhg, -3.0, 3.0, 0, 10));
367:     PetscCall(VecGetLocalSize(U, &Np));
368:     Np /= dim;
369:     PetscCall(VecGetArrayRead(U, &u));
370:     /* get points from solution vector */
371:     for (p = 0; p < Np; ++p) PetscCall(PetscDrawHGAddValue(user->drawhg, u[p]));
372:     PetscCall(PetscDrawHGDraw(user->drawhg));
373:     PetscCall(VecRestoreArrayRead(U, &u));
374:   }
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: static PetscErrorCode SPMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
379: {
380:   AppCtx            *user = (AppCtx *)ctx;
381:   const PetscScalar *u;
382:   PetscReal         *v, *coords;
383:   PetscInt           Np, p;
384:   DM                 dmSw;

386:   PetscFunctionBeginUser;
387:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
388:   if ((user->ostep > 0) && (!(step % user->ostep))) {
389:     PetscDrawAxis axis;

391:     PetscCall(TSGetDM(ts, &dmSw));
392:     PetscCall(PetscDrawSPReset(user->drawsp));
393:     PetscCall(PetscDrawSPGetAxis(user->drawsp, &axis));
394:     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "V", "w"));
395:     PetscCall(VecGetLocalSize(U, &Np));
396:     PetscCall(VecGetArrayRead(U, &u));
397:     /* get points from solution vector */
398:     PetscCall(PetscMalloc1(Np, &v));
399:     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
400:     PetscCall(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
401:     for (p = 0; p < Np - 1; ++p) PetscCall(PetscDrawSPAddPoint(user->drawsp, &coords[p], &v[p]));
402:     PetscCall(PetscDrawSPDraw(user->drawsp, PETSC_TRUE));
403:     PetscCall(VecRestoreArrayRead(U, &u));
404:     PetscCall(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
405:     PetscCall(PetscFree(v));
406:   }
407:   PetscFunctionReturn(PETSC_SUCCESS);
408: }

410: static PetscErrorCode KSConv(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
411: {
412:   AppCtx            *user = (AppCtx *)ctx;
413:   const PetscScalar *u;
414:   PetscScalar        sup;
415:   PetscReal         *v, *coords, T = 0., vel = 0., step_cast, w_sum;
416:   PetscInt           dim, Np, p, cStart, cEnd;
417:   DM                 sw, plex;

419:   PetscFunctionBeginUser;
420:   if (step < 0) PetscFunctionReturn(PETSC_SUCCESS);
421:   if ((user->ostep > 0) && (!(step % user->ostep))) {
422:     PetscDrawAxis axis;
423:     PetscCall(PetscDrawSPGetAxis(user->drawks, &axis));
424:     PetscCall(PetscDrawAxisSetLabels(axis, "Particles", "t", "D_n"));
425:     PetscCall(PetscDrawSPSetLimits(user->drawks, 0., 100, 0., 3.5));
426:     PetscCall(TSGetDM(ts, &sw));
427:     PetscCall(VecGetLocalSize(U, &Np));
428:     PetscCall(VecGetArrayRead(U, &u));
429:     /* get points from solution vector */
430:     PetscCall(PetscMalloc1(Np, &v));
431:     PetscCall(DMSwarmGetCellDM(sw, &plex));
432:     PetscCall(DMGetDimension(sw, &dim));
433:     PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
434:     for (p = 0; p < Np; ++p) v[p] = PetscRealPart(u[p]);
435:     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
436:     w_sum = 0.;
437:     for (p = 0; p < Np; ++p) {
438:       w_sum += u[p];
439:       T += u[p] * coords[p] * coords[p];
440:       vel += u[p] * coords[p];
441:     }
442:     vel /= w_sum;
443:     T   = T / w_sum - vel * vel;
444:     sup = 0.0;
445:     for (p = 0; p < Np; ++p) {
446:       PetscReal tmp = 0.;
447:       tmp           = PetscAbs(u[p] - ComputePDF(1.0, w_sum, T, &coords[p * dim]));
448:       if (tmp > sup) sup = tmp;
449:     }
450:     step_cast = (PetscReal)step;
451:     PetscCall(PetscDrawSPAddPoint(user->drawks, &step_cast, &sup));
452:     PetscCall(PetscDrawSPDraw(user->drawks, PETSC_TRUE));
453:     PetscCall(VecRestoreArrayRead(U, &u));
454:     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
455:     PetscCall(PetscFree(v));
456:   }
457:   PetscFunctionReturn(PETSC_SUCCESS);
458: }

460: static PetscErrorCode InitializeSolve(TS ts, Vec u)
461: {
462:   DM      dm;
463:   AppCtx *user;

465:   PetscFunctionBeginUser;
466:   PetscCall(TSGetDM(ts, &dm));
467:   PetscCall(DMGetApplicationContext(dm, &user));
468:   PetscCall(SetInitialCoordinates(dm));
469:   PetscCall(SetInitialConditions(dm, u));
470:   PetscFunctionReturn(PETSC_SUCCESS);
471: }
472: /*
473:      A single particle is generated for each velocity space cell using the dmswarmpicfield_coor and is used to evaluate collisions in that cell.
474:      0 weight ghost particles are initialized outside of a small velocity domain to ensure the tails of the amxwellian are resolved.
475:    */
476: int main(int argc, char **argv)
477: {
478:   TS       ts;     /* nonlinear solver */
479:   DM       dm, sw; /* Velocity space mesh and Particle Swarm */
480:   Vec      u, w;   /* swarm vector */
481:   MPI_Comm comm;
482:   AppCtx   user;

484:   PetscFunctionBeginUser;
485:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
486:   comm = PETSC_COMM_WORLD;
487:   PetscCall(ProcessOptions(comm, &user));
488:   PetscCall(CreateMesh(comm, &dm, &user));
489:   PetscCall(CreateParticles(dm, &sw, &user));
490:   PetscCall(DMSetApplicationContext(sw, &user));
491:   PetscCall(TSCreate(comm, &ts));
492:   PetscCall(TSSetDM(ts, sw));
493:   PetscCall(TSSetMaxTime(ts, 10.0));
494:   PetscCall(TSSetTimeStep(ts, 0.01));
495:   PetscCall(TSSetMaxSteps(ts, 100000));
496:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
497:   if (user.monitorhg) {
498:     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
499:     PetscCall(PetscDrawSetFromOptions(user.draw));
500:     PetscCall(PetscDrawHGCreate(user.draw, 20, &user.drawhg));
501:     PetscCall(PetscDrawHGSetColor(user.drawhg, 3));
502:     PetscCall(TSMonitorSet(ts, HGMonitor, &user, NULL));
503:   } else if (user.monitorsp) {
504:     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
505:     PetscCall(PetscDrawSetFromOptions(user.draw));
506:     PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawsp));
507:     PetscCall(TSMonitorSet(ts, SPMonitor, &user, NULL));
508:   } else if (user.monitorks) {
509:     PetscCall(PetscDrawCreate(comm, NULL, "monitor", 0, 0, 400, 300, &user.draw));
510:     PetscCall(PetscDrawSetFromOptions(user.draw));
511:     PetscCall(PetscDrawSPCreate(user.draw, 1, &user.drawks));
512:     PetscCall(TSMonitorSet(ts, KSConv, &user, NULL));
513:   }
514:   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user));
515:   PetscCall(TSSetFromOptions(ts));
516:   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
517:   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &w));
518:   PetscCall(VecDuplicate(w, &u));
519:   PetscCall(VecCopy(w, u));
520:   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &w));
521:   PetscCall(TSComputeInitialCondition(ts, u));
522:   PetscCall(TSSolve(ts, u));
523:   if (user.monitorhg) {
524:     PetscCall(PetscDrawSave(user.draw));
525:     PetscCall(PetscDrawHGDestroy(&user.drawhg));
526:     PetscCall(PetscDrawDestroy(&user.draw));
527:   }
528:   if (user.monitorsp) {
529:     PetscCall(PetscDrawSave(user.draw));
530:     PetscCall(PetscDrawSPDestroy(&user.drawsp));
531:     PetscCall(PetscDrawDestroy(&user.draw));
532:   }
533:   if (user.monitorks) {
534:     PetscCall(PetscDrawSave(user.draw));
535:     PetscCall(PetscDrawSPDestroy(&user.drawks));
536:     PetscCall(PetscDrawDestroy(&user.draw));
537:   }
538:   PetscCall(VecDestroy(&u));
539:   PetscCall(TSDestroy(&ts));
540:   PetscCall(DMDestroy(&sw));
541:   PetscCall(DMDestroy(&dm));
542:   PetscCall(PetscFinalize());
543:   return 0;
544: }

546: /*TEST
547:    build:
548:      requires: double !complex
549:    test:
550:      suffix: 1
551:      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
552:    test:
553:      suffix: 2
554:      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
555: TEST*/