xref: /petsc/src/ts/tests/ex27.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 static char help[] = "Particle Basis Landau Example using nonlinear solve + Implicit Midpoint-like time stepping.";
2 
3 /* TODO
4 
5 1) SNES is sensitive to epsilon (but not to h). Should we do continuation in it?
6 
7 2) Put this timestepper in library, maybe by changing DG
8 
9 3) Add monitor to visualize distributions
10 
11 */
12 
13 /* References
14   [1] https://arxiv.org/abs/1910.03080v2
15 */
16 
17 #include <petscdmplex.h>
18 #include <petscdmswarm.h>
19 #include <petscts.h>
20 #include <petscviewer.h>
21 #include <petscmath.h>
22 
23 typedef struct {
24   /* Velocity space grid and functions */
25   PetscInt  N;         /* The number of partices per spatial cell */
26   PetscReal L;         /* Velocity space is [-L, L]^d */
27   PetscReal h;         /* Spacing for grid 2L / N^{1/d} */
28   PetscReal epsilon;   /* gaussian regularization parameter */
29   PetscReal momentTol; /* Tolerance for checking moment conservation */
30 } AppCtx;
31 
32 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33 {
34   PetscErrorCode ierr;
35 
36   PetscFunctionBeginUser;
37   options->N         = 1;
38   options->momentTol = 100.0*PETSC_MACHINE_EPSILON;
39   options->L         = 1.0;
40   options->h         = -1.0;
41   options->epsilon   = -1.0;
42 
43   ierr = PetscOptionsBegin(comm, "", "Collision Options", "DMPLEX");PetscCall(ierr);
44   PetscCall(PetscOptionsInt("-N", "Number of particles per spatial cell", "ex27.c", options->N, &options->N, NULL));
45   PetscCall(PetscOptionsReal("-L", "Velocity-space extent", "ex27.c", options->L, &options->L, NULL));
46   PetscCall(PetscOptionsReal("-h", "Velocity-space resolution", "ex27.c", options->h, &options->h, NULL));
47   PetscCall(PetscOptionsReal("-epsilon", "Mollifier regularization parameter", "ex27.c", options->epsilon, &options->epsilon, NULL));
48   ierr = PetscOptionsEnd();PetscCall(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
53 {
54   PetscFunctionBeginUser;
55   PetscCall(DMCreate(comm, dm));
56   PetscCall(DMSetType(*dm, DMPLEX));
57   PetscCall(DMSetFromOptions(*dm));
58   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
59   PetscFunctionReturn(0);
60 }
61 
62 static PetscErrorCode SetInitialCoordinates(DM sw)
63 {
64   AppCtx        *user;
65   PetscRandom    rnd, rndv;
66   DM             dm;
67   DMPolytopeType ct;
68   PetscBool      simplex;
69   PetscReal     *centroid, *coords, *velocity, *xi0, *v0, *J, *invJ, detJ, *vals;
70   PetscInt       dim, d, cStart, cEnd, c, Np, p;
71 
72   PetscFunctionBeginUser;
73   /* Randomization for coordinates */
74   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rnd));
75   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
76   PetscCall(PetscRandomSetFromOptions(rnd));
77   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) sw), &rndv));
78   PetscCall(PetscRandomSetInterval(rndv, -1., 1.));
79   PetscCall(PetscRandomSetFromOptions(rndv));
80   PetscCall(DMGetApplicationContext(sw, &user));
81   Np   = user->N;
82   PetscCall(DMGetDimension(sw, &dim));
83   PetscCall(DMSwarmGetCellDM(sw, &dm));
84   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
85   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
86   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
87   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
88   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
89   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
90   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &velocity));
91   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &vals));
92   for (c = cStart; c < cEnd; ++c) {
93     if (Np == 1) {
94       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
95       for (d = 0; d < dim; ++d) {
96         coords[c*dim+d] = centroid[d];
97       }
98       vals[c] = 1.0;
99     } else {
100       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
101       for (p = 0; p < Np; ++p) {
102         const PetscInt n   = c*Np + p;
103         PetscReal      sum = 0.0, refcoords[3];
104 
105         for (d = 0; d < dim; ++d) {
106           PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d]));
107           sum += refcoords[d];
108         }
109         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
110         vals[n] = 1.0;
111         PetscCall(DMPlexReferenceToCoordinates(dm, c, 1, refcoords, &coords[n*dim]));
112       }
113     }
114   }
115   /* Random velocity IC */
116   for (c = cStart; c < cEnd; ++c) {
117     for (p = 0; p < Np; ++p) {
118       for (d = 0; d < dim; ++d) {
119         PetscReal v_val;
120 
121         PetscCall(PetscRandomGetValueReal(rndv, &v_val));
122         velocity[p*dim+d] = v_val;
123       }
124     }
125   }
126   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
127   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &velocity));
128   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &vals));
129   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
130   PetscCall(PetscRandomDestroy(&rnd));
131   PetscCall(PetscRandomDestroy(&rndv));
132   PetscFunctionReturn(0);
133 }
134 
135 /* Get velocities from swarm and place in solution vector */
136 static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
137 {
138   DM             dm;
139   AppCtx        *user;
140   PetscReal     *velocity;
141   PetscScalar   *initialConditions;
142   PetscInt       dim, d, cStart, cEnd, c, Np, p, n;
143 
144   PetscFunctionBeginUser;
145   PetscCall(VecGetLocalSize(u, &n));
146   PetscCall(DMGetApplicationContext(dmSw, &user));
147   Np   = user->N;
148   PetscCall(DMSwarmGetCellDM(dmSw, &dm));
149   PetscCall(DMGetDimension(dm, &dim));
150   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
151   PetscCall(DMSwarmGetField(dmSw, "velocity", NULL, NULL, (void **) &velocity));
152   PetscCall(VecGetArray(u, &initialConditions));
153   for (c = cStart; c < cEnd; ++c) {
154     for (p = 0; p < Np; ++p) {
155       const PetscInt n = c*Np + p;
156       for (d = 0; d < dim; d++) {
157         initialConditions[n*dim+d] = velocity[n*dim+d];
158       }
159     }
160   }
161   PetscCall(VecRestoreArray(u, &initialConditions));
162   PetscCall(DMSwarmRestoreField(dmSw, "velocity", NULL, NULL, (void **) &velocity));
163   PetscFunctionReturn(0);
164 }
165 
166 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
167 {
168   PetscInt      *cellid;
169   PetscInt       dim, cStart, cEnd, c, Np = user->N, p;
170   PetscBool      view = PETSC_FALSE;
171 
172   PetscFunctionBeginUser;
173   PetscCall(DMGetDimension(dm, &dim));
174   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
175   PetscCall(DMSetType(*sw, DMSWARM));
176   PetscCall(DMSetDimension(*sw, dim));
177   /* h = 2L/n and N = n^d */
178   if (user->h < 0.) user->h = 2.*user->L / PetscPowReal(user->N, 1./dim);
179   /* From Section 4 in [1], \epsilon = 0.64 h^.98 */
180   if (user->epsilon < 0.) user->epsilon = 0.64*pow(user->h, 1.98);
181   PetscCall(PetscOptionsGetBool(NULL, NULL, "-param_view", &view, NULL));
182   if (view) {
183     PetscCall(PetscPrintf(PETSC_COMM_SELF, "N: %D L: %g h: %g eps: %g\n", user->N, user->L, user->h, user->epsilon));
184   }
185   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
186   PetscCall(DMSwarmSetCellDM(*sw, dm));
187   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
188   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
189   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
190   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
191   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
192   PetscCall(DMSetFromOptions(*sw));
193   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
194   for (c = cStart; c < cEnd; ++c) {
195     for (p = 0; p < Np; ++p) {
196       const PetscInt n = c*Np + p;
197       cellid[n] = c;
198     }
199   }
200   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
201   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
202   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
203   PetscFunctionReturn(0);
204 }
205 
206 /* Internal dmplex function, same as found in dmpleximpl.h */
207 static void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
208 {
209   PetscInt d;
210 
211   for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d];
212 }
213 
214 /* Internal dmplex function, same as found in dmpleximpl.h */
215 static PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
216 {
217   PetscReal sum = 0.0;
218   PetscInt d;
219 
220   for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d])*y[d];
221   return sum;
222 }
223 
224 /* Internal dmplex function, same as found in dmpleximpl.h */
225 static void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
226 {
227   PetscScalar z[2];
228   z[0] = x[0]; z[1] = x[ldx];
229   y[0]   += A[0]*z[0] + A[1]*z[1];
230   y[ldx] += A[2]*z[0] + A[3]*z[1];
231   (void)PetscLogFlops(6.0);
232 }
233 
234 /* Internal dmplex function, same as found in dmpleximpl.h to avoid private includes. */
235 static void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
236 {
237   PetscScalar z[3];
238   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
239   y[0]     += A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
240   y[ldx]   += A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
241   y[ldx*2] += A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
242   (void)PetscLogFlops(15.0);
243 }
244 
245 /*
246   Gaussian - The Gaussian function G(x)
247 
248   Input Parameters:
249 +  dim   - The number of dimensions, or size of x
250 .  mu    - The mean, or center
251 .  sigma - The standard deviation, or width
252 -  x     - The evaluation point of the function
253 
254   Output Parameter:
255 . ret - The value G(x)
256 */
257 static PetscReal Gaussian(PetscInt dim, const PetscReal mu[], PetscReal sigma, const PetscReal x[])
258 {
259   PetscReal arg = 0.0;
260   PetscInt  d;
261 
262   for (d = 0; d < dim; ++d) arg += PetscSqr(x[d] - mu[d]);
263   return PetscPowReal(2.0*PETSC_PI*sigma, -dim/2.0) * PetscExpReal(-arg/(2.0*sigma));
264 }
265 
266 /*
267   ComputeGradS - Compute grad_v dS_eps/df
268 
269   Input Parameters:
270 + dim      - The dimension
271 . Np       - The number of particles
272 . vp       - The velocity v_p of the particle at which we evaluate
273 . velocity - The velocity field for all particles
274 . epsilon  - The regularization strength
275 
276   Output Parameter:
277 . integral - The output grad_v dS_eps/df (v_p)
278 
279   Note:
280   This comes from (3.6) in [1], and we are computing
281 $   \nabla_v S_p = \grad \psi_\epsilon(v_p - v) log \sum_q \psi_\epsilon(v - v_q)
282   which is discretized by using a one-point quadrature in each box l at its center v^c_l
283 $   \sum_l h^d \nabla\psi_\epsilon(v_p - v^c_l) \log\left( \sum_q w_q \psi_\epsilon(v^c_l - v_q) \right)
284   where h^d is the volume of each box.
285 */
286 static PetscErrorCode ComputeGradS(PetscInt dim, PetscInt Np, const PetscReal vp[], const PetscReal velocity[], PetscReal integral[], AppCtx *ctx)
287 {
288   PetscReal vc_l[3], L = ctx->L, h = ctx->h, epsilon = ctx->epsilon, init = 0.5*h - L;
289   PetscInt  nx = roundf(2.*L / h);
290   PetscInt  ny = dim > 1 ? nx : 1;
291   PetscInt  nz = dim > 2 ? nx : 1;
292   PetscInt  i, j, k, d, q, dbg = 0;
293 
294   PetscFunctionBeginHot;
295   for (d = 0; d < dim; ++d) integral[d] = 0.0;
296   for (k = 0, vc_l[2] = init; k < nz; ++k, vc_l[2] += h) {
297     for (j = 0, vc_l[1] = init; j < ny; ++j, vc_l[1] += h) {
298       for (i = 0, vc_l[0] = init; i < nx; ++i, vc_l[0] += h) {
299         PetscReal sum = 0.0;
300 
301         if (dbg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(%D %D) vc_l: %g %g\n", i, j, vc_l[0], vc_l[1]));
302         /* \log \sum_k \psi(v - v_k)  */
303         for (q = 0; q < Np; ++q) sum += Gaussian(dim, &velocity[q*dim], epsilon, vc_l);
304         sum = PetscLogReal(sum);
305         for (d = 0; d < dim; ++d) integral[d] += (-1./(epsilon))*PetscAbsReal(vp[d] - vc_l[d])*(Gaussian(dim, vp, epsilon, vc_l)) * sum;
306       }
307     }
308   }
309   PetscFunctionReturn(0);
310 }
311 
312 /* Q = 1/|xi| (I - xi xi^T / |xi|^2), xi = vp - vq */
313 static PetscErrorCode QCompute(PetscInt dim, const PetscReal vp[], const PetscReal vq[], PetscReal Q[])
314 {
315   PetscReal xi[3], xi2, xi3, mag;
316   PetscInt  d, e;
317 
318   PetscFunctionBeginHot;
319   DMPlex_WaxpyD_Internal(dim, -1.0, vq, vp, xi);
320   xi2 = DMPlex_DotD_Internal(dim, xi, xi);
321   mag = PetscSqrtReal(xi2);
322   xi3 = xi2 * mag;
323   for (d = 0; d < dim; ++d) {
324     for (e = 0; e < dim; ++e) {
325       Q[d*dim+e] = -xi[d]*xi[e] / xi3;
326     }
327     Q[d*dim+d] += 1. / mag;
328   }
329   PetscFunctionReturn(0);
330 }
331 
332 static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
333 {
334   AppCtx            *user = (AppCtx *) ctx;
335   PetscInt           dbg  = 0;
336   DM                 sw;                  /* Particles */
337   Vec                sol;                 /* Solution vector at current time */
338   const PetscScalar *u;                   /* input solution vector */
339   PetscScalar       *r;
340   PetscReal         *velocity;
341   PetscInt           dim, Np, p, q;
342 
343   PetscFunctionBeginUser;
344   PetscCall(VecZeroEntries(R));
345   PetscCall(TSGetDM(ts, &sw));
346   PetscCall(DMGetDimension(sw, &dim));
347   PetscCall(VecGetLocalSize(U, &Np));
348   PetscCall(TSGetSolution(ts, &sol));
349   PetscCall(VecGetArray(sol, &velocity));
350   PetscCall(VecGetArray(R, &r));
351   PetscCall(VecGetArrayRead(U, &u));
352   Np  /= dim;
353   if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Part  ppr     x        y\n"));
354   for (p = 0; p < Np; ++p) {
355     PetscReal gradS_p[3] = {0., 0., 0.};
356 
357     PetscCall(ComputeGradS(dim, Np, &velocity[p*dim], velocity, gradS_p, user));
358     for (q = 0; q < Np; ++q) {
359       PetscReal gradS_q[3] = {0., 0., 0.}, GammaS[3] = {0., 0., 0.}, Q[9];
360 
361       if (q == p) continue;
362       PetscCall(ComputeGradS(dim, Np, &velocity[q*dim], velocity, gradS_q, user));
363       DMPlex_WaxpyD_Internal(dim, -1.0, gradS_q, gradS_p, GammaS);
364       PetscCall(QCompute(dim, &u[p*dim], &u[q*dim], Q));
365       switch (dim) {
366         case 2: DMPlex_MultAdd2DReal_Internal(Q, 1, GammaS, &r[p*dim]);break;
367         case 3: DMPlex_MultAdd3DReal_Internal(Q, 1, GammaS, &r[p*dim]);break;
368         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Do not support dimension %D", dim);
369       }
370     }
371     if (dbg) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final %4D %10.8lf %10.8lf\n", p, r[p*dim+0], r[p*dim+1]));
372   }
373   PetscCall(VecRestoreArrayRead(U, &u));
374   PetscCall(VecRestoreArray(R, &r));
375   PetscCall(VecRestoreArray(sol, &velocity));
376   PetscCall(VecViewFromOptions(R, NULL, "-residual_view"));
377   PetscFunctionReturn(0);
378 }
379 
380 /*
381  TS Post Step Function. Copy the solution back into the swarm for migration. We may also need to reform
382  the solution vector in cases of particle migration, but we forgo that here since there is no velocity space grid
383  to migrate between.
384 */
385 static PetscErrorCode UpdateSwarm(TS ts)
386 {
387   PetscInt idx, n;
388   const PetscScalar *u;
389   PetscScalar *velocity;
390   DM sw;
391   Vec sol;
392 
393   PetscFunctionBeginUser;
394   PetscCall(TSGetDM(ts, &sw));
395   PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **) &velocity));
396   PetscCall(TSGetSolution(ts, &sol));
397   PetscCall(VecGetArrayRead(sol, &u));
398   PetscCall(VecGetLocalSize(sol, &n));
399   for (idx = 0; idx < n; ++idx) velocity[idx] = u[idx];
400   PetscCall(VecRestoreArrayRead(sol, &u));
401   PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **) &velocity));
402   PetscFunctionReturn(0);
403 }
404 
405 static PetscErrorCode InitializeSolve(TS ts, Vec u)
406 {
407   DM             dm;
408   AppCtx        *user;
409 
410   PetscFunctionBeginUser;
411   PetscCall(TSGetDM(ts, &dm));
412   PetscCall(DMGetApplicationContext(dm, &user));
413   PetscCall(SetInitialCoordinates(dm));
414   PetscCall(SetInitialConditions(dm, u));
415   PetscFunctionReturn(0);
416 }
417 
418 int main(int argc,char **argv)
419 {
420   TS             ts;     /* nonlinear solver */
421   DM             dm, sw; /* Velocity space mesh and Particle Swarm */
422   Vec            u, v;   /* problem vector */
423   MPI_Comm       comm;
424   AppCtx         user;
425 
426   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
427   comm = PETSC_COMM_WORLD;
428   PetscCall(ProcessOptions(comm, &user));
429   /* Initialize objects and set initial conditions */
430   PetscCall(CreateMesh(comm, &dm, &user));
431   PetscCall(CreateParticles(dm, &sw, &user));
432   PetscCall(DMSetApplicationContext(sw, &user));
433   PetscCall(DMSwarmVectorDefineField(sw, "velocity"));
434   PetscCall(TSCreate(comm, &ts));
435   PetscCall(TSSetDM(ts, sw));
436   PetscCall(TSSetMaxTime(ts, 10.0));
437   PetscCall(TSSetTimeStep(ts, 0.1));
438   PetscCall(TSSetMaxSteps(ts, 1));
439   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
440   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user));
441   PetscCall(TSSetFromOptions(ts));
442   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
443   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &v));
444   PetscCall(VecDuplicate(v, &u));
445   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &v));
446   PetscCall(TSComputeInitialCondition(ts, u));
447   PetscCall(TSSetPostStep(ts, UpdateSwarm));
448   PetscCall(TSSolve(ts, u));
449   PetscCall(VecDestroy(&u));
450   PetscCall(TSDestroy(&ts));
451   PetscCall(DMDestroy(&sw));
452   PetscCall(DMDestroy(&dm));
453   PetscCall(PetscFinalize());
454   return 0;
455 }
456 
457 /*TEST
458    build:
459      requires: triangle !single !complex
460    test:
461      suffix: midpoint
462      args: -N 3 -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1 -dm_view \
463            -ts_type theta -ts_theta_theta 0.5 -ts_dmswarm_monitor_moments -ts_monitor_frequency 1 -snes_fd
464 TEST*/
465