xref: /petsc/src/ts/impls/symplectic/basicsymplectic/basicsymplectic.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1e49d4f37SHong Zhang /*
2e49d4f37SHong Zhang   Code for Timestepping with basic symplectic integrators for separable Hamiltonian systems
3e49d4f37SHong Zhang */
4e49d4f37SHong Zhang #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
5e49d4f37SHong Zhang #include <petscdm.h>
6e49d4f37SHong Zhang 
7e49d4f37SHong Zhang static TSBasicSymplecticType TSBasicSymplecticDefault = TSBASICSYMPLECTICSIEULER;
8e49d4f37SHong Zhang static PetscBool             TSBasicSymplecticRegisterAllCalled;
9e49d4f37SHong Zhang static PetscBool             TSBasicSymplecticPackageInitialized;
10e49d4f37SHong Zhang 
11e49d4f37SHong Zhang typedef struct _BasicSymplecticScheme     *BasicSymplecticScheme;
12e49d4f37SHong Zhang typedef struct _BasicSymplecticSchemeLink *BasicSymplecticSchemeLink;
13e49d4f37SHong Zhang 
14e49d4f37SHong Zhang struct _BasicSymplecticScheme {
15e49d4f37SHong Zhang   char      *name;
16e49d4f37SHong Zhang   PetscInt   order;
17e49d4f37SHong Zhang   PetscInt   s; /* number of stages */
18e49d4f37SHong Zhang   PetscReal *c, *d;
19e49d4f37SHong Zhang };
20e49d4f37SHong Zhang struct _BasicSymplecticSchemeLink {
21e49d4f37SHong Zhang   struct _BasicSymplecticScheme sch;
22e49d4f37SHong Zhang   BasicSymplecticSchemeLink     next;
23e49d4f37SHong Zhang };
24e49d4f37SHong Zhang static BasicSymplecticSchemeLink BasicSymplecticSchemeList;
25e49d4f37SHong Zhang typedef struct {
26e49d4f37SHong Zhang   TS                    subts_p, subts_q; /* sub TS contexts that holds the RHSFunction pointers */
27e49d4f37SHong Zhang   IS                    is_p, is_q;       /* IS sets for position and momentum respectively */
28e49d4f37SHong Zhang   Vec                   update;           /* a nest work vector for generalized coordinates */
29e49d4f37SHong Zhang   BasicSymplecticScheme scheme;
30e49d4f37SHong Zhang } TS_BasicSymplectic;
31e49d4f37SHong Zhang 
32e49d4f37SHong Zhang /*MC
33e49d4f37SHong Zhang   TSBASICSYMPLECTICSIEULER - first order semi-implicit Euler method
34bcf0153eSBarry Smith 
35e49d4f37SHong Zhang   Level: intermediate
36bcf0153eSBarry Smith 
37*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`
38e49d4f37SHong Zhang M*/
39e49d4f37SHong Zhang 
40e49d4f37SHong Zhang /*MC
41aaa8cc7dSPierre Jolivet   TSBASICSYMPLECTICVELVERLET - second order Velocity Verlet method (leapfrog method with starting process and determining velocity and position at the same time)
42bcf0153eSBarry Smith 
43e49d4f37SHong Zhang Level: intermediate
44bcf0153eSBarry Smith 
45*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`
46e49d4f37SHong Zhang M*/
47e49d4f37SHong Zhang 
48e49d4f37SHong Zhang /*@C
49bcf0153eSBarry Smith   TSBasicSymplecticRegisterAll - Registers all of the basic symplectic integration methods in `TSBASICSYMPLECTIC`
50e49d4f37SHong Zhang 
51e49d4f37SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
52e49d4f37SHong Zhang 
53e49d4f37SHong Zhang   Level: advanced
54e49d4f37SHong Zhang 
55*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegisterDestroy()`
56e49d4f37SHong Zhang @*/
57d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticRegisterAll(void)
58d71ae5a4SJacob Faibussowitsch {
59e49d4f37SHong Zhang   PetscFunctionBegin;
603ba16761SJacob Faibussowitsch   if (TSBasicSymplecticRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
61e49d4f37SHong Zhang   TSBasicSymplecticRegisterAllCalled = PETSC_TRUE;
62e49d4f37SHong Zhang   {
63e49d4f37SHong Zhang     PetscReal c[1] = {1.0}, d[1] = {1.0};
649566063dSJacob Faibussowitsch     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICSIEULER, 1, 1, c, d));
65e49d4f37SHong Zhang   }
66e49d4f37SHong Zhang   {
67e49d4f37SHong Zhang     PetscReal c[2] = {0, 1.0}, d[2] = {0.5, 0.5};
689566063dSJacob Faibussowitsch     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICVELVERLET, 2, 2, c, d));
69e49d4f37SHong Zhang   }
70e49d4f37SHong Zhang   {
71e49d4f37SHong Zhang     PetscReal c[3] = {1, -2.0 / 3.0, 2.0 / 3.0}, d[3] = {-1.0 / 24.0, 3.0 / 4.0, 7.0 / 24.0};
729566063dSJacob Faibussowitsch     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC3, 3, 3, c, d));
73e49d4f37SHong Zhang   }
74e49d4f37SHong Zhang   {
75e49d4f37SHong Zhang #define CUBEROOTOFTWO 1.2599210498948731647672106
76e49d4f37SHong Zhang     PetscReal c[4] = {1.0 / 2.0 / (2.0 - CUBEROOTOFTWO), (1.0 - CUBEROOTOFTWO) / 2.0 / (2.0 - CUBEROOTOFTWO), (1.0 - CUBEROOTOFTWO) / 2.0 / (2.0 - CUBEROOTOFTWO), 1.0 / 2.0 / (2.0 - CUBEROOTOFTWO)}, d[4] = {1.0 / (2.0 - CUBEROOTOFTWO), -CUBEROOTOFTWO / (2.0 - CUBEROOTOFTWO), 1.0 / (2.0 - CUBEROOTOFTWO), 0};
779566063dSJacob Faibussowitsch     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC4, 4, 4, c, d));
78e49d4f37SHong Zhang   }
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80e49d4f37SHong Zhang }
81e49d4f37SHong Zhang 
82e49d4f37SHong Zhang /*@C
83bcf0153eSBarry Smith    TSBasicSymplecticRegisterDestroy - Frees the list of schemes that were registered by `TSBasicSymplecticRegister()`.
84e49d4f37SHong Zhang 
85e49d4f37SHong Zhang    Not Collective
86e49d4f37SHong Zhang 
87e49d4f37SHong Zhang    Level: advanced
88e49d4f37SHong Zhang 
89*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticRegister()`, `TSBasicSymplecticRegisterAll()`, `TSBASICSYMPLECTIC`
90e49d4f37SHong Zhang @*/
91d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticRegisterDestroy(void)
92d71ae5a4SJacob Faibussowitsch {
93e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
94e49d4f37SHong Zhang 
95e49d4f37SHong Zhang   PetscFunctionBegin;
96e49d4f37SHong Zhang   while ((link = BasicSymplecticSchemeList)) {
97e49d4f37SHong Zhang     BasicSymplecticScheme scheme = &link->sch;
98e49d4f37SHong Zhang     BasicSymplecticSchemeList    = link->next;
999566063dSJacob Faibussowitsch     PetscCall(PetscFree2(scheme->c, scheme->d));
1009566063dSJacob Faibussowitsch     PetscCall(PetscFree(scheme->name));
1019566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
102e49d4f37SHong Zhang   }
103e49d4f37SHong Zhang   TSBasicSymplecticRegisterAllCalled = PETSC_FALSE;
1043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105e49d4f37SHong Zhang }
106e49d4f37SHong Zhang 
107e49d4f37SHong Zhang /*@C
108bcf0153eSBarry Smith   TSBasicSymplecticInitializePackage - This function initializes everything in the `TSBASICSYMPLECTIC` package. It is called
109bcf0153eSBarry Smith   from `TSInitializePackage()`.
110e49d4f37SHong Zhang 
111e49d4f37SHong Zhang   Level: developer
112e49d4f37SHong Zhang 
113*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSBASICSYMPLECTIC`
114e49d4f37SHong Zhang @*/
115d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticInitializePackage(void)
116d71ae5a4SJacob Faibussowitsch {
117e49d4f37SHong Zhang   PetscFunctionBegin;
1183ba16761SJacob Faibussowitsch   if (TSBasicSymplecticPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
119e49d4f37SHong Zhang   TSBasicSymplecticPackageInitialized = PETSC_TRUE;
1209566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticRegisterAll());
1219566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSBasicSymplecticFinalizePackage));
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123e49d4f37SHong Zhang }
124e49d4f37SHong Zhang 
125e49d4f37SHong Zhang /*@C
126bcf0153eSBarry Smith   TSBasicSymplecticFinalizePackage - This function destroys everything in the `TSBASICSYMPLECTIC` package. It is
127bcf0153eSBarry Smith   called from `PetscFinalize()`.
128e49d4f37SHong Zhang 
129e49d4f37SHong Zhang   Level: developer
130e49d4f37SHong Zhang 
131*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSBASICSYMPLECTIC`
132e49d4f37SHong Zhang @*/
133d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticFinalizePackage(void)
134d71ae5a4SJacob Faibussowitsch {
135e49d4f37SHong Zhang   PetscFunctionBegin;
136e49d4f37SHong Zhang   TSBasicSymplecticPackageInitialized = PETSC_FALSE;
1379566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticRegisterDestroy());
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139e49d4f37SHong Zhang }
140e49d4f37SHong Zhang 
141e49d4f37SHong Zhang /*@C
142e49d4f37SHong Zhang    TSBasicSymplecticRegister - register a basic symplectic integration scheme by providing the coefficients.
143e49d4f37SHong Zhang 
144e49d4f37SHong Zhang    Not Collective, but the same schemes should be registered on all processes on which they will be used
145e49d4f37SHong Zhang 
146e49d4f37SHong Zhang    Input Parameters:
147e49d4f37SHong Zhang +  name - identifier for method
148e49d4f37SHong Zhang .  order - approximation order of method
149e49d4f37SHong Zhang .  s - number of stages, this is the dimension of the matrices below
150e49d4f37SHong Zhang .  c - coefficients for updating generalized position (dimension s)
151e49d4f37SHong Zhang -  d - coefficients for updating generalized momentum (dimension s)
152e49d4f37SHong Zhang 
153bcf0153eSBarry Smith    Level: advanced
154bcf0153eSBarry Smith 
155e49d4f37SHong Zhang    Notes:
156e49d4f37SHong Zhang    Several symplectic methods are provided, this function is only needed to create new methods.
157e49d4f37SHong Zhang 
158*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`
159e49d4f37SHong Zhang @*/
160d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticRegister(TSRosWType name, PetscInt order, PetscInt s, PetscReal c[], PetscReal d[])
161d71ae5a4SJacob Faibussowitsch {
162e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
163e49d4f37SHong Zhang   BasicSymplecticScheme     scheme;
164e49d4f37SHong Zhang 
165e49d4f37SHong Zhang   PetscFunctionBegin;
166e49d4f37SHong Zhang   PetscValidCharPointer(name, 1);
167dadcf809SJacob Faibussowitsch   PetscValidRealPointer(c, 4);
168dadcf809SJacob Faibussowitsch   PetscValidRealPointer(d, 5);
169e49d4f37SHong Zhang 
1709566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticInitializePackage());
1719566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
172e49d4f37SHong Zhang   scheme = &link->sch;
1739566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &scheme->name));
174e49d4f37SHong Zhang   scheme->order = order;
175e49d4f37SHong Zhang   scheme->s     = s;
1769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(s, &scheme->c, s, &scheme->d));
1779566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(scheme->c, c, s));
1789566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(scheme->d, d, s));
179e49d4f37SHong Zhang   link->next                = BasicSymplecticSchemeList;
180e49d4f37SHong Zhang   BasicSymplecticSchemeList = link;
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182e49d4f37SHong Zhang }
183e49d4f37SHong Zhang 
184e49d4f37SHong Zhang /*
185e49d4f37SHong Zhang The simplified form of the equations are:
186e49d4f37SHong Zhang 
18737fdd005SBarry Smith .vb
18837fdd005SBarry Smith  p_{i+1} = p_i + c_i*g(q_i)*h
18937fdd005SBarry Smith  q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h
19037fdd005SBarry Smith .ve
191e49d4f37SHong Zhang 
192e49d4f37SHong Zhang Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p.
193e49d4f37SHong Zhang 
194e49d4f37SHong Zhang To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps:
19537fdd005SBarry Smith .vb
196e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration multiplied by c_1
197e49d4f37SHong Zhang - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1
198e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2
199e49d4f37SHong Zhang - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2
20037fdd005SBarry Smith .ve
201e49d4f37SHong Zhang 
202e49d4f37SHong Zhang */
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_BasicSymplectic(TS ts)
204d71ae5a4SJacob Faibussowitsch {
205e49d4f37SHong Zhang   TS_BasicSymplectic   *bsymp    = (TS_BasicSymplectic *)ts->data;
206e49d4f37SHong Zhang   BasicSymplecticScheme scheme   = bsymp->scheme;
207e49d4f37SHong Zhang   Vec                   solution = ts->vec_sol, update = bsymp->update, q, p, q_update, p_update;
208e49d4f37SHong Zhang   IS                    is_q = bsymp->is_q, is_p = bsymp->is_p;
209e49d4f37SHong Zhang   TS                    subts_q = bsymp->subts_q, subts_p = bsymp->subts_p;
210e49d4f37SHong Zhang   PetscBool             stageok;
211e49d4f37SHong Zhang   PetscReal             next_time_step = ts->time_step;
212e49d4f37SHong Zhang   PetscInt              iter;
213e49d4f37SHong Zhang 
214e49d4f37SHong Zhang   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(solution, is_q, &q));
2169566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(solution, is_p, &p));
2179566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(update, is_q, &q_update));
2189566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(update, is_p, &p_update));
219e49d4f37SHong Zhang 
220e49d4f37SHong Zhang   for (iter = 0; iter < scheme->s; iter++) {
2219566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, ts->ptime));
222e49d4f37SHong Zhang     /* update velocity p */
223e49d4f37SHong Zhang     if (scheme->c[iter]) {
2249566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(subts_p, ts->ptime, q, p_update));
2259566063dSJacob Faibussowitsch       PetscCall(VecAXPY(p, scheme->c[iter] * ts->time_step, p_update));
226e49d4f37SHong Zhang     }
227e49d4f37SHong Zhang     /* update position q */
228e49d4f37SHong Zhang     if (scheme->d[iter]) {
2299566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(subts_q, ts->ptime, p, q_update));
2309566063dSJacob Faibussowitsch       PetscCall(VecAXPY(q, scheme->d[iter] * ts->time_step, q_update));
231e49d4f37SHong Zhang       ts->ptime = ts->ptime + scheme->d[iter] * ts->time_step;
232e49d4f37SHong Zhang     }
2339566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, ts->ptime, 0, &solution));
2349566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok));
2359371c9d4SSatish Balay     if (!stageok) {
2369371c9d4SSatish Balay       ts->reason = TS_DIVERGED_STEP_REJECTED;
2373ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
2389371c9d4SSatish Balay     }
2399566063dSJacob Faibussowitsch     PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok));
2409371c9d4SSatish Balay     if (!stageok) {
2419371c9d4SSatish Balay       ts->reason = TS_DIVERGED_STEP_REJECTED;
2423ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
2439371c9d4SSatish Balay     }
244e49d4f37SHong Zhang   }
245e49d4f37SHong Zhang 
246e49d4f37SHong Zhang   ts->time_step = next_time_step;
2479566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(solution, is_q, &q));
2489566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(solution, is_p, &p));
2499566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(update, is_q, &q_update));
2509566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(update, is_p, &p_update));
2513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
252e49d4f37SHong Zhang }
253e49d4f37SHong Zhang 
254d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, void *ctx)
255d71ae5a4SJacob Faibussowitsch {
256e49d4f37SHong Zhang   PetscFunctionBegin;
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258e49d4f37SHong Zhang }
259e49d4f37SHong Zhang 
260d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
261d71ae5a4SJacob Faibussowitsch {
262e49d4f37SHong Zhang   PetscFunctionBegin;
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
264e49d4f37SHong Zhang }
265e49d4f37SHong Zhang 
266d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, void *ctx)
267d71ae5a4SJacob Faibussowitsch {
268e49d4f37SHong Zhang   PetscFunctionBegin;
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270e49d4f37SHong Zhang }
271e49d4f37SHong Zhang 
272d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
273d71ae5a4SJacob Faibussowitsch {
274e49d4f37SHong Zhang   PetscFunctionBegin;
2753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
276e49d4f37SHong Zhang }
277e49d4f37SHong Zhang 
278d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BasicSymplectic(TS ts)
279d71ae5a4SJacob Faibussowitsch {
280e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
281e49d4f37SHong Zhang   DM                  dm;
282e49d4f37SHong Zhang 
283e49d4f37SHong Zhang   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q));
2859566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p));
286aaa8cc7dSPierre Jolivet   PetscCheck(bsymp->is_q && bsymp->is_p, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up RHSSplits with TSRHSSplitSetIS() using split names position and momentum respectively in order to use -ts_type basicsymplectic");
2879566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q));
2889566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p));
2893c633725SBarry Smith   PetscCheck(bsymp->subts_q && bsymp->subts_p, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up the RHSFunctions for position and momentum using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS");
290e49d4f37SHong Zhang 
2919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update));
292e49d4f37SHong Zhang 
2939566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
2949566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */
2959566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
296e49d4f37SHong Zhang   if (dm) {
2979566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts));
2989566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts));
299e49d4f37SHong Zhang   }
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301e49d4f37SHong Zhang }
302e49d4f37SHong Zhang 
303d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BasicSymplectic(TS ts)
304d71ae5a4SJacob Faibussowitsch {
305e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
306e49d4f37SHong Zhang 
307e49d4f37SHong Zhang   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bsymp->update));
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
310e49d4f37SHong Zhang }
311e49d4f37SHong Zhang 
312d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BasicSymplectic(TS ts)
313d71ae5a4SJacob Faibussowitsch {
314e49d4f37SHong Zhang   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(TSReset_BasicSymplectic(ts));
3162e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL));
3172e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL));
3189566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
3193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
320e49d4f37SHong Zhang }
321e49d4f37SHong Zhang 
322d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems *PetscOptionsObject)
323d71ae5a4SJacob Faibussowitsch {
324e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
325e49d4f37SHong Zhang 
326e49d4f37SHong Zhang   PetscFunctionBegin;
327d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options");
328e49d4f37SHong Zhang   {
329e49d4f37SHong Zhang     BasicSymplecticSchemeLink link;
330e49d4f37SHong Zhang     PetscInt                  count, choice;
331e49d4f37SHong Zhang     PetscBool                 flg;
332e49d4f37SHong Zhang     const char              **namelist;
333e49d4f37SHong Zhang 
3349371c9d4SSatish Balay     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++)
3359371c9d4SSatish Balay       ;
3369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
337e49d4f37SHong Zhang     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) namelist[count] = link->sch.name;
3389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_basicsymplectic_type", "Family of basic symplectic integration method", "TSBasicSymplecticSetType", (const char *const *)namelist, count, bsymp->scheme->name, &choice, &flg));
3399566063dSJacob Faibussowitsch     if (flg) PetscCall(TSBasicSymplecticSetType(ts, namelist[choice]));
3409566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
341e49d4f37SHong Zhang   }
342d0609cedSBarry Smith   PetscOptionsHeadEnd();
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
344e49d4f37SHong Zhang }
345e49d4f37SHong Zhang 
346d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BasicSymplectic(TS ts, PetscViewer viewer)
347d71ae5a4SJacob Faibussowitsch {
348e49d4f37SHong Zhang   PetscFunctionBegin;
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
350e49d4f37SHong Zhang }
351e49d4f37SHong Zhang 
352d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X)
353d71ae5a4SJacob Faibussowitsch {
354e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp  = (TS_BasicSymplectic *)ts->data;
355e49d4f37SHong Zhang   Vec                 update = bsymp->update;
356e49d4f37SHong Zhang   PetscReal           alpha  = (ts->ptime - t) / ts->time_step;
357e49d4f37SHong Zhang 
358e49d4f37SHong Zhang   PetscFunctionBegin;
3599566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol));
3609566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol));
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
362e49d4f37SHong Zhang }
363e49d4f37SHong Zhang 
364d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
365d71ae5a4SJacob Faibussowitsch {
366e49d4f37SHong Zhang   PetscFunctionBegin;
367e49d4f37SHong Zhang   *yr = 1.0 + xr;
368e49d4f37SHong Zhang   *yi = xi;
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
370e49d4f37SHong Zhang }
371e49d4f37SHong Zhang 
372e49d4f37SHong Zhang /*@C
373e49d4f37SHong Zhang   TSBasicSymplecticSetType - Set the type of the basic symplectic method
374e49d4f37SHong Zhang 
375c3339decSBarry Smith   Logically Collective
376e49d4f37SHong Zhang 
377d8d19677SJose E. Roman   Input Parameters:
378e49d4f37SHong Zhang +  ts - timestepping context
379e49d4f37SHong Zhang -  bsymptype - type of the symplectic scheme
380e49d4f37SHong Zhang 
381bcf0153eSBarry Smith   Options Database Key:
382147403d9SBarry Smith .  -ts_basicsymplectic_type <scheme> - select the scheme
383e49d4f37SHong Zhang 
384bcf0153eSBarry Smith   Level: intermediate
385bcf0153eSBarry Smith 
386bcf0153eSBarry Smith   Note:
387bcf0153eSBarry Smith     The symplectic solver always expects a two-way splitting with the split names being "position" and "momentum". Each split is associated with an `IS` object and a sub-`TS`
388147403d9SBarry Smith     that is intended to store the user-provided RHS function.
389e49d4f37SHong Zhang 
390*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
391e49d4f37SHong Zhang @*/
392d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype)
393d71ae5a4SJacob Faibussowitsch {
394e49d4f37SHong Zhang   PetscFunctionBegin;
395e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
396cac4c232SBarry Smith   PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype));
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398e49d4f37SHong Zhang }
399e49d4f37SHong Zhang 
400e49d4f37SHong Zhang /*@C
401e49d4f37SHong Zhang   TSBasicSymplecticGetType - Get the type of the basic symplectic method
402e49d4f37SHong Zhang 
403c3339decSBarry Smith   Logically Collective
404e49d4f37SHong Zhang 
405d8d19677SJose E. Roman   Input Parameters:
406e49d4f37SHong Zhang +  ts - timestepping context
407e49d4f37SHong Zhang -  bsymptype - type of the basic symplectic scheme
408e49d4f37SHong Zhang 
409e49d4f37SHong Zhang   Level: intermediate
410bcf0153eSBarry Smith 
411*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
412e49d4f37SHong Zhang @*/
413d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype)
414d71ae5a4SJacob Faibussowitsch {
415e49d4f37SHong Zhang   PetscFunctionBegin;
416e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
417cac4c232SBarry Smith   PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype));
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
419e49d4f37SHong Zhang }
420e49d4f37SHong Zhang 
421d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype)
422d71ae5a4SJacob Faibussowitsch {
423e49d4f37SHong Zhang   TS_BasicSymplectic       *bsymp = (TS_BasicSymplectic *)ts->data;
424e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
425e49d4f37SHong Zhang   PetscBool                 match;
426e49d4f37SHong Zhang 
427e49d4f37SHong Zhang   PetscFunctionBegin;
428e49d4f37SHong Zhang   if (bsymp->scheme) {
4299566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match));
4303ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
431e49d4f37SHong Zhang   }
432e49d4f37SHong Zhang   for (link = BasicSymplecticSchemeList; link; link = link->next) {
4339566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match));
434e49d4f37SHong Zhang     if (match) {
435e49d4f37SHong Zhang       bsymp->scheme = &link->sch;
4363ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
437e49d4f37SHong Zhang     }
438e49d4f37SHong Zhang   }
43998921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", bsymptype);
440e49d4f37SHong Zhang }
441e49d4f37SHong Zhang 
442d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts, TSBasicSymplecticType *bsymptype)
443d71ae5a4SJacob Faibussowitsch {
444e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
445e49d4f37SHong Zhang 
446e49d4f37SHong Zhang   PetscFunctionBegin;
447e49d4f37SHong Zhang   *bsymptype = bsymp->scheme->name;
4483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
449e49d4f37SHong Zhang }
450e49d4f37SHong Zhang 
451e49d4f37SHong Zhang /*MC
452bcf0153eSBarry Smith   TSBASICSYMPLECTIC - ODE solver using basic symplectic integration schemes
453e49d4f37SHong Zhang 
454bcf0153eSBarry Smith   These methods are intended for separable Hamiltonian systems
455bcf0153eSBarry Smith .vb
456bcf0153eSBarry Smith   qdot = dH(q,p,t)/dp
457bcf0153eSBarry Smith   pdot = -dH(q,p,t)/dq
458bcf0153eSBarry Smith .ve
459e49d4f37SHong Zhang 
460e49d4f37SHong Zhang   where the Hamiltonian can be split into the sum of kinetic energy and potential energy
461bcf0153eSBarry Smith .vb
462bcf0153eSBarry Smith   H(q,p,t) = T(p,t) + V(q,t).
463bcf0153eSBarry Smith .ve
464e49d4f37SHong Zhang 
465e49d4f37SHong Zhang   As a result, the system can be genearlly represented by
466bcf0153eSBarry Smith .vb
467bcf0153eSBarry Smith   qdot = f(p,t) = dT(p,t)/dp
468bcf0153eSBarry Smith   pdot = g(q,t) = -dV(q,t)/dq
469bcf0153eSBarry Smith .ve
470e49d4f37SHong Zhang 
471e49d4f37SHong Zhang   and solved iteratively with
472bcf0153eSBarry Smith .vb
473bcf0153eSBarry Smith   q_new = q_old + d_i*h*f(p_old,t_old)
474bcf0153eSBarry Smith   t_new = t_old + d_i*h
475bcf0153eSBarry Smith   p_new = p_old + c_i*h*g(p_new,t_new)
476bcf0153eSBarry Smith   i=0,1,...,n.
477bcf0153eSBarry Smith .ve
478e49d4f37SHong Zhang 
479bcf0153eSBarry Smith   The solution vector should contain both q and p, which correspond to (generalized) position and momentum respectively. Note that the momentum component
480bcf0153eSBarry Smith   could simply be velocity in some representations. The symplectic solver always expects a two-way splitting with the split names being "position" and "momentum".
481bcf0153eSBarry Smith   Each split is associated with an `IS` object and a sub-`TS` that is intended to store the user-provided RHS function.
482e49d4f37SHong Zhang 
483e49d4f37SHong Zhang   Level: beginner
484e49d4f37SHong Zhang 
485bcf0153eSBarry Smith   Reference:
486bcf0153eSBarry Smith . * -  wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator)
487e49d4f37SHong Zhang 
488*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()`, `TSType`
489e49d4f37SHong Zhang M*/
490d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts)
491d71ae5a4SJacob Faibussowitsch {
492e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp;
493e49d4f37SHong Zhang 
494e49d4f37SHong Zhang   PetscFunctionBegin;
4959566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticInitializePackage());
4964dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&bsymp));
497e49d4f37SHong Zhang   ts->data = (void *)bsymp;
498e49d4f37SHong Zhang 
499e49d4f37SHong Zhang   ts->ops->setup           = TSSetUp_BasicSymplectic;
500e49d4f37SHong Zhang   ts->ops->step            = TSStep_BasicSymplectic;
501e49d4f37SHong Zhang   ts->ops->reset           = TSReset_BasicSymplectic;
502e49d4f37SHong Zhang   ts->ops->destroy         = TSDestroy_BasicSymplectic;
503e49d4f37SHong Zhang   ts->ops->setfromoptions  = TSSetFromOptions_BasicSymplectic;
504e49d4f37SHong Zhang   ts->ops->view            = TSView_BasicSymplectic;
505e49d4f37SHong Zhang   ts->ops->interpolate     = TSInterpolate_BasicSymplectic;
506e49d4f37SHong Zhang   ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic;
507e49d4f37SHong Zhang 
5089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", TSBasicSymplecticSetType_BasicSymplectic));
5099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", TSBasicSymplecticGetType_BasicSymplectic));
510e49d4f37SHong Zhang 
5119566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticSetType(ts, TSBasicSymplecticDefault));
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
513e49d4f37SHong Zhang }
514