xref: /petsc/src/dm/impls/swarm/tests/ex5.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Vlasov example of central orbits\n";
2 
3 /*
4   To visualize the orbit, we can used
5 
6     -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -1 -ts_monitor_sp_swarm_phase 0 -draw_size 500,500
7 
8   and we probably want it to run fast and not check convergence
9 
10     -convest_num_refine 0 -ts_dt 0.01 -ts_max_steps 100 -ts_max_time 100 -output_step 10
11 */
12 
13 #include <petscts.h>
14 #include <petscdmplex.h>
15 #include <petscdmswarm.h>
16 #include <petsc/private/dmpleximpl.h> /* For norm and dot */
17 
18 PETSC_EXTERN PetscErrorCode circleSingleX(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
19 PETSC_EXTERN PetscErrorCode circleSingleV(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
20 PETSC_EXTERN PetscErrorCode circleMultipleX(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
21 PETSC_EXTERN PetscErrorCode circleMultipleV(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
22 
23 typedef struct {
24   PetscBool error; /* Flag for printing the error */
25   PetscInt  ostep; /* print the energy at each ostep time steps */
26 } AppCtx;
27 
28 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
29   PetscFunctionBeginUser;
30   options->error = PETSC_FALSE;
31   options->ostep = 100;
32 
33   PetscOptionsBegin(comm, "", "Central Orbit Options", "DMSWARM");
34   PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex5.c", options->error, &options->error, NULL));
35   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex5.c", options->ostep, &options->ostep, PETSC_NULL));
36   PetscOptionsEnd();
37   PetscFunctionReturn(0);
38 }
39 
40 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
41   PetscFunctionBeginUser;
42   PetscCall(DMCreate(comm, dm));
43   PetscCall(DMSetType(*dm, DMPLEX));
44   PetscCall(DMSetFromOptions(*dm));
45   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
46   PetscFunctionReturn(0);
47 }
48 
49 static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw) {
50   PetscReal v0[1] = {1.};
51   PetscInt  dim;
52 
53   PetscFunctionBeginUser;
54   PetscCall(DMGetDimension(dm, &dim));
55   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
56   PetscCall(DMSetType(*sw, DMSWARM));
57   PetscCall(DMSetDimension(*sw, dim));
58   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
59   PetscCall(DMSwarmSetCellDM(*sw, dm));
60   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
61   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
62   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
63   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL));
64   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initVelocity", dim, PETSC_REAL));
65   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
66   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
67   PetscCall(DMSwarmInitializeCoordinates(*sw));
68   PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0));
69   PetscCall(DMSetFromOptions(*sw));
70   PetscCall(DMSetApplicationContext(*sw, user));
71   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
72   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
73   {
74     Vec gc, gc0, gv, gv0;
75 
76     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
77     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0));
78     PetscCall(VecCopy(gc, gc0));
79     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
80     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0));
81     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "velocity", &gv));
82     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initVelocity", &gv0));
83     PetscCall(VecCopy(gv, gv0));
84     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "velocity", &gv));
85     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initVelocity", &gv0));
86   }
87   PetscFunctionReturn(0);
88 }
89 
90 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx) {
91   DM                 sw;
92   const PetscReal   *coords, *vel;
93   const PetscScalar *u;
94   PetscScalar       *g;
95   PetscInt           dim, d, Np, p;
96 
97   PetscFunctionBeginUser;
98   PetscCall(TSGetDM(ts, &sw));
99   PetscCall(DMGetDimension(sw, &dim));
100   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
101   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
102   PetscCall(VecGetLocalSize(U, &Np));
103   PetscCall(VecGetArrayRead(U, &u));
104   PetscCall(VecGetArray(G, &g));
105   Np /= 2 * dim;
106   for (p = 0; p < Np; ++p) {
107     const PetscReal x0    = coords[p * dim + 0];
108     const PetscReal vy0   = vel[p * dim + 1];
109     const PetscReal omega = vy0 / x0;
110 
111     for (d = 0; d < dim; ++d) {
112       g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d];
113       g[(p * 2 + 1) * dim + d] = -PetscSqr(omega) * u[(p * 2 + 0) * dim + d];
114     }
115   }
116   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
117   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
118   PetscCall(VecRestoreArrayRead(U, &u));
119   PetscCall(VecRestoreArray(G, &g));
120   PetscFunctionReturn(0);
121 }
122 
123 /* J_{ij} = dF_i/dx_j
124    J_p = (  0   1)
125          (-w^2  0)
126 */
127 static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, void *ctx) {
128   DM               sw;
129   const PetscReal *coords, *vel;
130   PetscInt         dim, d, Np, p, rStart;
131 
132   PetscFunctionBeginUser;
133   PetscCall(TSGetDM(ts, &sw));
134   PetscCall(DMGetDimension(sw, &dim));
135   PetscCall(VecGetLocalSize(U, &Np));
136   PetscCall(MatGetOwnershipRange(J, &rStart, NULL));
137   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
138   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
139   Np /= 2 * dim;
140   for (p = 0; p < Np; ++p) {
141     const PetscReal x0      = coords[p * dim + 0];
142     const PetscReal vy0     = vel[p * dim + 1];
143     const PetscReal omega   = vy0 / x0;
144     PetscScalar     vals[4] = {0., 1., -PetscSqr(omega), 0.};
145 
146     for (d = 0; d < dim; ++d) {
147       const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
148       PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
149     }
150   }
151   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
152   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
153   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
154   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
155   PetscFunctionReturn(0);
156 }
157 
158 static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) {
159   const PetscScalar *v;
160   PetscScalar       *xres;
161   PetscInt           Np, p;
162 
163   PetscFunctionBeginUser;
164   PetscCall(VecGetLocalSize(Xres, &Np));
165   PetscCall(VecGetArrayRead(V, &v));
166   PetscCall(VecGetArray(Xres, &xres));
167   for (p = 0; p < Np; ++p) xres[p] = v[p];
168   PetscCall(VecRestoreArrayRead(V, &v));
169   PetscCall(VecRestoreArray(Xres, &xres));
170   PetscFunctionReturn(0);
171 }
172 
173 static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) {
174   DM                 sw;
175   const PetscScalar *x;
176   const PetscReal   *coords, *vel;
177   PetscScalar       *vres;
178   PetscInt           Np, p, dim, d;
179 
180   PetscFunctionBeginUser;
181   PetscCall(TSGetDM(ts, &sw));
182   PetscCall(DMGetDimension(sw, &dim));
183   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
184   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
185   PetscCall(VecGetLocalSize(Vres, &Np));
186   PetscCall(VecGetArrayRead(X, &x));
187   PetscCall(VecGetArray(Vres, &vres));
188   PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension must be 2");
189   Np /= dim;
190   for (p = 0; p < Np; ++p) {
191     const PetscReal x0    = coords[p * dim + 0];
192     const PetscReal vy0   = vel[p * dim + 1];
193     const PetscReal omega = vy0 / x0;
194 
195     for (d = 0; d < dim; ++d) vres[p * dim + d] = -PetscSqr(omega) * x[p * dim + d];
196   }
197   PetscCall(VecRestoreArrayRead(X, &x));
198   PetscCall(VecRestoreArray(Vres, &vres));
199   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
200   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
201   PetscFunctionReturn(0);
202 }
203 
204 static PetscErrorCode CreateSolution(TS ts) {
205   DM       sw;
206   Vec      u;
207   PetscInt dim, Np;
208 
209   PetscFunctionBegin;
210   PetscCall(TSGetDM(ts, &sw));
211   PetscCall(DMGetDimension(sw, &dim));
212   PetscCall(DMSwarmGetLocalSize(sw, &Np));
213   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
214   PetscCall(VecSetBlockSize(u, dim));
215   PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE));
216   PetscCall(VecSetUp(u));
217   PetscCall(TSSetSolution(ts, u));
218   PetscCall(VecDestroy(&u));
219   PetscFunctionReturn(0);
220 }
221 
222 static PetscErrorCode SetProblem(TS ts) {
223   AppCtx *user;
224   DM      sw;
225 
226   PetscFunctionBegin;
227   PetscCall(TSGetDM(ts, &sw));
228   PetscCall(DMGetApplicationContext(sw, (void **)&user));
229   // Define unified system for (X, V)
230   {
231     Mat      J;
232     PetscInt dim, Np;
233 
234     PetscCall(DMGetDimension(sw, &dim));
235     PetscCall(DMSwarmGetLocalSize(sw, &Np));
236     PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
237     PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE));
238     PetscCall(MatSetBlockSize(J, 2 * dim));
239     PetscCall(MatSetFromOptions(J));
240     PetscCall(MatSetUp(J));
241     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
242     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &user));
243     PetscCall(MatDestroy(&J));
244   }
245   // Define split system for X and V
246   {
247     Vec             u;
248     IS              isx, isv, istmp;
249     const PetscInt *idx;
250     PetscInt        dim, Np, rstart;
251 
252     PetscCall(TSGetSolution(ts, &u));
253     PetscCall(DMGetDimension(sw, &dim));
254     PetscCall(DMSwarmGetLocalSize(sw, &Np));
255     PetscCall(VecGetOwnershipRange(u, &rstart, NULL));
256     PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 0, 2, &istmp));
257     PetscCall(ISGetIndices(istmp, &idx));
258     PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx));
259     PetscCall(ISRestoreIndices(istmp, &idx));
260     PetscCall(ISDestroy(&istmp));
261     PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 1, 2, &istmp));
262     PetscCall(ISGetIndices(istmp, &idx));
263     PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv));
264     PetscCall(ISRestoreIndices(istmp, &idx));
265     PetscCall(ISDestroy(&istmp));
266     PetscCall(TSRHSSplitSetIS(ts, "position", isx));
267     PetscCall(TSRHSSplitSetIS(ts, "momentum", isv));
268     PetscCall(ISDestroy(&isx));
269     PetscCall(ISDestroy(&isv));
270     PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, &user));
271     PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, &user));
272   }
273   // Define symplectic formulation U_t = S . G, where G = grad F
274   {
275     //PetscCall(TSDiscGradSetFormulation(ts, RHSJacobianS, RHSObjectiveF, RHSFunctionG, &user));
276   }
277   PetscFunctionReturn(0);
278 }
279 
280 static PetscErrorCode DMSwarmTSRedistribute(TS ts) {
281   DM        sw;
282   Vec       u;
283   PetscReal t, maxt, dt;
284   PetscInt  n, maxn;
285 
286   PetscFunctionBegin;
287   PetscCall(TSGetDM(ts, &sw));
288   PetscCall(TSGetTime(ts, &t));
289   PetscCall(TSGetMaxTime(ts, &maxt));
290   PetscCall(TSGetTimeStep(ts, &dt));
291   PetscCall(TSGetStepNumber(ts, &n));
292   PetscCall(TSGetMaxSteps(ts, &maxn));
293 
294   PetscCall(TSReset(ts));
295   PetscCall(TSSetDM(ts, sw));
296   /* TODO Check whether TS was set from options */
297   PetscCall(TSSetFromOptions(ts));
298   PetscCall(TSSetTime(ts, t));
299   PetscCall(TSSetMaxTime(ts, maxt));
300   PetscCall(TSSetTimeStep(ts, dt));
301   PetscCall(TSSetStepNumber(ts, n));
302   PetscCall(TSSetMaxSteps(ts, maxn));
303 
304   PetscCall(CreateSolution(ts));
305   PetscCall(SetProblem(ts));
306   PetscCall(TSGetSolution(ts, &u));
307   PetscFunctionReturn(0);
308 }
309 
310 PetscErrorCode circleSingleX(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx) {
311   x[0] = p + 1.;
312   x[1] = 0.;
313   return 0;
314 }
315 
316 PetscErrorCode circleSingleV(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar v[], void *ctx) {
317   v[0] = 0.;
318   v[1] = PetscSqrtReal(1000. / (p + 1.));
319   return 0;
320 }
321 
322 /* Put 5 particles into each circle */
323 PetscErrorCode circleMultipleX(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx) {
324   const PetscInt  n   = 5;
325   const PetscReal r0  = (p / n) + 1.;
326   const PetscReal th0 = (2. * PETSC_PI * (p % n)) / n;
327 
328   x[0] = r0 * PetscCosReal(th0);
329   x[1] = r0 * PetscSinReal(th0);
330   return 0;
331 }
332 
333 /* Put 5 particles into each circle */
334 PetscErrorCode circleMultipleV(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar v[], void *ctx) {
335   const PetscInt  n     = 5;
336   const PetscReal r0    = (p / n) + 1.;
337   const PetscReal th0   = (2. * PETSC_PI * (p % n)) / n;
338   const PetscReal omega = PetscSqrtReal(1000. / r0) / r0;
339 
340   v[0] = -r0 * omega * PetscSinReal(th0);
341   v[1] = r0 * omega * PetscCosReal(th0);
342   return 0;
343 }
344 
345 /*
346   InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values.
347 
348   Input Parameters:
349 + ts         - The TS
350 - useInitial - Flag to also set the initial conditions to the current coodinates and velocities and setup the problem
351 
352   Output Parameters:
353 . u - The initialized solution vector
354 
355   Level: advanced
356 
357 .seealso: InitializeSolve()
358 */
359 static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial) {
360   DM      sw;
361   Vec     u, gc, gv, gc0, gv0;
362   IS      isx, isv;
363   AppCtx *user;
364 
365   PetscFunctionBeginUser;
366   PetscCall(TSGetDM(ts, &sw));
367   PetscCall(DMGetApplicationContext(sw, &user));
368   if (useInitial) {
369     PetscReal v0[1] = {1.};
370 
371     PetscCall(DMSwarmInitializeCoordinates(sw));
372     PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0));
373     PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
374     PetscCall(DMSwarmTSRedistribute(ts));
375   }
376   PetscCall(TSGetSolution(ts, &u));
377   PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
378   PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
379   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
380   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0));
381   if (useInitial) PetscCall(VecCopy(gc, gc0));
382   PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc));
383   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
384   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0));
385   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
386   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initVelocity", &gv0));
387   if (useInitial) PetscCall(VecCopy(gv, gv0));
388   PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv));
389   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
390   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initVelocity", &gv0));
391   PetscFunctionReturn(0);
392 }
393 
394 static PetscErrorCode InitializeSolve(TS ts, Vec u) {
395   PetscFunctionBegin;
396   PetscCall(TSSetSolution(ts, u));
397   PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE));
398   PetscFunctionReturn(0);
399 }
400 
401 static PetscErrorCode ComputeError(TS ts, Vec U, Vec E) {
402   MPI_Comm           comm;
403   DM                 sw;
404   AppCtx            *user;
405   const PetscScalar *u;
406   const PetscReal   *coords, *vel;
407   PetscScalar       *e;
408   PetscReal          t;
409   PetscInt           dim, Np, p;
410 
411   PetscFunctionBeginUser;
412   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
413   PetscCall(TSGetDM(ts, &sw));
414   PetscCall(DMGetApplicationContext(sw, &user));
415   PetscCall(DMGetDimension(sw, &dim));
416   PetscCall(TSGetSolveTime(ts, &t));
417   PetscCall(VecGetArray(E, &e));
418   PetscCall(VecGetArrayRead(U, &u));
419   PetscCall(VecGetLocalSize(U, &Np));
420   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
421   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
422   Np /= 2 * dim;
423   for (p = 0; p < Np; ++p) {
424     /* TODO generalize initial conditions and project into plane instead of assuming x-y */
425     const PetscReal    r0    = DMPlex_NormD_Internal(dim, &coords[p * dim]);
426     const PetscReal    th0   = PetscAtan2Real(coords[p * dim + 1], coords[p * dim + 0]);
427     const PetscReal    v0    = DMPlex_NormD_Internal(dim, &vel[p * dim]);
428     const PetscReal    omega = v0 / r0;
429     const PetscReal    ct    = PetscCosReal(omega * t + th0);
430     const PetscReal    st    = PetscSinReal(omega * t + th0);
431     const PetscScalar *x     = &u[(p * 2 + 0) * dim];
432     const PetscScalar *v     = &u[(p * 2 + 1) * dim];
433     const PetscReal    xe[3] = {r0 * ct, r0 * st, 0.0};
434     const PetscReal    ve[3] = {-v0 * st, v0 * ct, 0.0};
435     PetscInt           d;
436 
437     for (d = 0; d < dim; ++d) {
438       e[(p * 2 + 0) * dim + d] = x[d] - xe[d];
439       e[(p * 2 + 1) * dim + d] = v[d] - ve[d];
440     }
441     if (user->error) {
442       const PetscReal en   = 0.5 * DMPlex_DotRealD_Internal(dim, v, v);
443       const PetscReal exen = 0.5 * PetscSqr(v0);
444       PetscCall(PetscPrintf(comm, "t %.4g: p%" PetscInt_FMT " error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g (%.10lf%%)\n", (double)t, p, (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 0) * dim]), (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 1) * dim]), (double)x[0], (double)x[1], (double)v[0], (double)v[1], (double)xe[0], (double)xe[1], (double)ve[0], (double)ve[1], (double)en, (double)exen, (double)(PetscAbsReal(exen - en) * 100. / exen)));
445     }
446   }
447   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
448   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
449   PetscCall(VecRestoreArrayRead(U, &u));
450   PetscCall(VecRestoreArray(E, &e));
451   PetscFunctionReturn(0);
452 }
453 
454 static PetscErrorCode EnergyMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) {
455   const PetscInt     ostep = ((AppCtx *)ctx)->ostep;
456   DM                 sw;
457   const PetscScalar *u;
458   PetscInt           dim, Np, p;
459 
460   PetscFunctionBeginUser;
461   if (step % ostep == 0) {
462     PetscCall(TSGetDM(ts, &sw));
463     PetscCall(DMGetDimension(sw, &dim));
464     PetscCall(VecGetArrayRead(U, &u));
465     PetscCall(VecGetLocalSize(U, &Np));
466     Np /= 2 * dim;
467     if (!step) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Time     Step Part     Energy\n"));
468     for (p = 0; p < Np; ++p) {
469       const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]);
470 
471       PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " %4" PetscInt_FMT " %10.4lf\n", (double)t, step, p, (double)(0.5 * v2)));
472     }
473     PetscCall(PetscSynchronizedFlush(PetscObjectComm((PetscObject)ts), NULL));
474     PetscCall(VecRestoreArrayRead(U, &u));
475   }
476   PetscFunctionReturn(0);
477 }
478 
479 static PetscErrorCode MigrateParticles(TS ts) {
480   DM sw;
481 
482   PetscFunctionBeginUser;
483   PetscCall(TSGetDM(ts, &sw));
484   PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre"));
485   {
486     Vec u, gc, gv;
487     IS  isx, isv;
488 
489     PetscCall(TSGetSolution(ts, &u));
490     PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
491     PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
492     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
493     PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc));
494     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
495     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
496     PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv));
497     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
498   }
499   PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
500   PetscCall(DMSwarmTSRedistribute(ts));
501   PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE));
502   //PetscCall(VecViewFromOptions(u, NULL, "-sol_migrate_view"));
503   PetscFunctionReturn(0);
504 }
505 
506 int main(int argc, char **argv) {
507   DM     dm, sw;
508   TS     ts;
509   Vec    u;
510   AppCtx user;
511 
512   PetscFunctionBeginUser;
513   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
514   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
515   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
516   PetscCall(CreateSwarm(dm, &user, &sw));
517   PetscCall(DMSetApplicationContext(sw, &user));
518 
519   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
520   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
521   PetscCall(TSSetDM(ts, sw));
522   PetscCall(TSSetMaxTime(ts, 0.1));
523   PetscCall(TSSetTimeStep(ts, 0.00001));
524   PetscCall(TSSetMaxSteps(ts, 100));
525   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
526   PetscCall(TSMonitorSet(ts, EnergyMonitor, &user, NULL));
527   PetscCall(TSSetFromOptions(ts));
528   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
529   PetscCall(TSSetComputeExactError(ts, ComputeError));
530   PetscCall(TSSetPostStep(ts, MigrateParticles));
531 
532   PetscCall(CreateSolution(ts));
533   PetscCall(TSGetSolution(ts, &u));
534   PetscCall(TSComputeInitialCondition(ts, u));
535   PetscCall(TSSolve(ts, NULL));
536 
537   PetscCall(TSDestroy(&ts));
538   PetscCall(DMDestroy(&sw));
539   PetscCall(DMDestroy(&dm));
540   PetscCall(PetscFinalize());
541   return 0;
542 }
543 
544 /*TEST
545 
546    build:
547      requires: double !complex
548 
549    testset:
550      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
551      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
552            -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
553            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
554            -dm_view -output_step 50 -error
555      test:
556        suffix: bsi_2d_1
557        args: -ts_basicsymplectic_type 1
558      test:
559        suffix: bsi_2d_2
560        args: -ts_basicsymplectic_type 2
561      test:
562        suffix: bsi_2d_3
563        args: -ts_basicsymplectic_type 3
564      test:
565        suffix: bsi_2d_4
566        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
567 
568    testset:
569      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
570      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
571            -ts_type theta -ts_theta_theta 0.5 -ts_convergence_estimate -convest_num_refine 2 \
572              -mat_type baij -ksp_error_if_not_converged -pc_type lu \
573            -dm_view -output_step 50 -error
574      test:
575        suffix: im_2d_0
576        args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5
577 
578    testset:
579      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
580      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 10,10 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 -petscpartitioner_type simple \
581            -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
582            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
583            -dm_view -output_step 50
584      test:
585        suffix: bsi_2d_mesh_1
586        args: -ts_basicsymplectic_type 1 -error
587      test:
588        suffix: bsi_2d_mesh_1_par_2
589        nsize: 2
590        args: -ts_basicsymplectic_type 1
591      test:
592        suffix: bsi_2d_mesh_1_par_3
593        nsize: 3
594        args: -ts_basicsymplectic_type 1
595      test:
596        suffix: bsi_2d_mesh_1_par_4
597        nsize: 4
598        args: -ts_basicsymplectic_type 1 -dm_swarm_num_particles 0,0,2,0
599 
600    testset:
601      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
602      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 10,10 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
603            -dm_swarm_num_particles 10 -dm_swarm_coordinate_function circleMultipleX -dm_swarm_velocity_function circleMultipleV \
604            -ts_convergence_estimate -convest_num_refine 2 \
605            -dm_view -output_step 50 -error
606      test:
607        suffix: bsi_2d_multiple_1
608        args: -ts_type basicsymplectic -ts_basicsymplectic_type 1
609      test:
610        suffix: bsi_2d_multiple_2
611        args: -ts_type basicsymplectic -ts_basicsymplectic_type 2
612      test:
613        suffix: bsi_2d_multiple_3
614        args: -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_dt 0.001
615      test:
616        suffix: im_2d_multiple_0
617        args: -ts_type theta -ts_theta_theta 0.5 \
618                -mat_type baij -ksp_error_if_not_converged -pc_type lu
619 
620 TEST*/
621