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