xref: /petsc/src/ts/impls/symplectic/basicsymplectic/basicsymplectic.c (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
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 
37bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`
38e49d4f37SHong Zhang M*/
39e49d4f37SHong Zhang 
40e49d4f37SHong Zhang /*MC
41*aaa8cc7dSPierre 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 
45bcf0153eSBarry Smith .seealso: [](chapter_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 
55bcf0153eSBarry Smith .seealso: [](chapter_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 
89bcf0153eSBarry Smith .seealso: [](chapter_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 
113bcf0153eSBarry Smith .seealso: [](chapter_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 
131bcf0153eSBarry Smith .seealso: [](chapter_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 
158bcf0153eSBarry Smith .seealso: [](chapter_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 
187e49d4f37SHong Zhang $ p_{i+1} = p_i + c_i*g(q_i)*h
188e49d4f37SHong Zhang $ q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h
189e49d4f37SHong Zhang 
190e49d4f37SHong Zhang Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p.
191e49d4f37SHong Zhang 
192e49d4f37SHong Zhang To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps:
193e49d4f37SHong Zhang 
194e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration multiplied by c_1
195e49d4f37SHong Zhang - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1
196e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2
197e49d4f37SHong Zhang - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2
198e49d4f37SHong Zhang 
199e49d4f37SHong Zhang */
200d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_BasicSymplectic(TS ts)
201d71ae5a4SJacob Faibussowitsch {
202e49d4f37SHong Zhang   TS_BasicSymplectic   *bsymp    = (TS_BasicSymplectic *)ts->data;
203e49d4f37SHong Zhang   BasicSymplecticScheme scheme   = bsymp->scheme;
204e49d4f37SHong Zhang   Vec                   solution = ts->vec_sol, update = bsymp->update, q, p, q_update, p_update;
205e49d4f37SHong Zhang   IS                    is_q = bsymp->is_q, is_p = bsymp->is_p;
206e49d4f37SHong Zhang   TS                    subts_q = bsymp->subts_q, subts_p = bsymp->subts_p;
207e49d4f37SHong Zhang   PetscBool             stageok;
208e49d4f37SHong Zhang   PetscReal             next_time_step = ts->time_step;
209e49d4f37SHong Zhang   PetscInt              iter;
210e49d4f37SHong Zhang 
211e49d4f37SHong Zhang   PetscFunctionBegin;
2129566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(solution, is_q, &q));
2139566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(solution, is_p, &p));
2149566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(update, is_q, &q_update));
2159566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(update, is_p, &p_update));
216e49d4f37SHong Zhang 
217e49d4f37SHong Zhang   for (iter = 0; iter < scheme->s; iter++) {
2189566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, ts->ptime));
219e49d4f37SHong Zhang     /* update velocity p */
220e49d4f37SHong Zhang     if (scheme->c[iter]) {
2219566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(subts_p, ts->ptime, q, p_update));
2229566063dSJacob Faibussowitsch       PetscCall(VecAXPY(p, scheme->c[iter] * ts->time_step, p_update));
223e49d4f37SHong Zhang     }
224e49d4f37SHong Zhang     /* update position q */
225e49d4f37SHong Zhang     if (scheme->d[iter]) {
2269566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(subts_q, ts->ptime, p, q_update));
2279566063dSJacob Faibussowitsch       PetscCall(VecAXPY(q, scheme->d[iter] * ts->time_step, q_update));
228e49d4f37SHong Zhang       ts->ptime = ts->ptime + scheme->d[iter] * ts->time_step;
229e49d4f37SHong Zhang     }
2309566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, ts->ptime, 0, &solution));
2319566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok));
2329371c9d4SSatish Balay     if (!stageok) {
2339371c9d4SSatish Balay       ts->reason = TS_DIVERGED_STEP_REJECTED;
2343ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
2359371c9d4SSatish Balay     }
2369566063dSJacob Faibussowitsch     PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok));
2379371c9d4SSatish Balay     if (!stageok) {
2389371c9d4SSatish Balay       ts->reason = TS_DIVERGED_STEP_REJECTED;
2393ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
2409371c9d4SSatish Balay     }
241e49d4f37SHong Zhang   }
242e49d4f37SHong Zhang 
243e49d4f37SHong Zhang   ts->time_step = next_time_step;
2449566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(solution, is_q, &q));
2459566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(solution, is_p, &p));
2469566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(update, is_q, &q_update));
2479566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(update, is_p, &p_update));
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
249e49d4f37SHong Zhang }
250e49d4f37SHong Zhang 
251d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, void *ctx)
252d71ae5a4SJacob Faibussowitsch {
253e49d4f37SHong Zhang   PetscFunctionBegin;
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255e49d4f37SHong Zhang }
256e49d4f37SHong Zhang 
257d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
258d71ae5a4SJacob Faibussowitsch {
259e49d4f37SHong Zhang   PetscFunctionBegin;
2603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
261e49d4f37SHong Zhang }
262e49d4f37SHong Zhang 
263d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, void *ctx)
264d71ae5a4SJacob Faibussowitsch {
265e49d4f37SHong Zhang   PetscFunctionBegin;
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267e49d4f37SHong Zhang }
268e49d4f37SHong Zhang 
269d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
270d71ae5a4SJacob Faibussowitsch {
271e49d4f37SHong Zhang   PetscFunctionBegin;
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273e49d4f37SHong Zhang }
274e49d4f37SHong Zhang 
275d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BasicSymplectic(TS ts)
276d71ae5a4SJacob Faibussowitsch {
277e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
278e49d4f37SHong Zhang   DM                  dm;
279e49d4f37SHong Zhang 
280e49d4f37SHong Zhang   PetscFunctionBegin;
2819566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q));
2829566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p));
283*aaa8cc7dSPierre 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");
2849566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q));
2859566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p));
2863c633725SBarry 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");
287e49d4f37SHong Zhang 
2889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update));
289e49d4f37SHong Zhang 
2909566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
2919566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */
2929566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
293e49d4f37SHong Zhang   if (dm) {
2949566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts));
2959566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts));
296e49d4f37SHong Zhang   }
2973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
298e49d4f37SHong Zhang }
299e49d4f37SHong Zhang 
300d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BasicSymplectic(TS ts)
301d71ae5a4SJacob Faibussowitsch {
302e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
303e49d4f37SHong Zhang 
304e49d4f37SHong Zhang   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bsymp->update));
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307e49d4f37SHong Zhang }
308e49d4f37SHong Zhang 
309d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BasicSymplectic(TS ts)
310d71ae5a4SJacob Faibussowitsch {
311e49d4f37SHong Zhang   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(TSReset_BasicSymplectic(ts));
3132e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL));
3142e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL));
3159566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
317e49d4f37SHong Zhang }
318e49d4f37SHong Zhang 
319d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems *PetscOptionsObject)
320d71ae5a4SJacob Faibussowitsch {
321e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
322e49d4f37SHong Zhang 
323e49d4f37SHong Zhang   PetscFunctionBegin;
324d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options");
325e49d4f37SHong Zhang   {
326e49d4f37SHong Zhang     BasicSymplecticSchemeLink link;
327e49d4f37SHong Zhang     PetscInt                  count, choice;
328e49d4f37SHong Zhang     PetscBool                 flg;
329e49d4f37SHong Zhang     const char              **namelist;
330e49d4f37SHong Zhang 
3319371c9d4SSatish Balay     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++)
3329371c9d4SSatish Balay       ;
3339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
334e49d4f37SHong Zhang     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) namelist[count] = link->sch.name;
3359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_basicsymplectic_type", "Family of basic symplectic integration method", "TSBasicSymplecticSetType", (const char *const *)namelist, count, bsymp->scheme->name, &choice, &flg));
3369566063dSJacob Faibussowitsch     if (flg) PetscCall(TSBasicSymplecticSetType(ts, namelist[choice]));
3379566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
338e49d4f37SHong Zhang   }
339d0609cedSBarry Smith   PetscOptionsHeadEnd();
3403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
341e49d4f37SHong Zhang }
342e49d4f37SHong Zhang 
343d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BasicSymplectic(TS ts, PetscViewer viewer)
344d71ae5a4SJacob Faibussowitsch {
345e49d4f37SHong Zhang   PetscFunctionBegin;
3463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
347e49d4f37SHong Zhang }
348e49d4f37SHong Zhang 
349d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X)
350d71ae5a4SJacob Faibussowitsch {
351e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp  = (TS_BasicSymplectic *)ts->data;
352e49d4f37SHong Zhang   Vec                 update = bsymp->update;
353e49d4f37SHong Zhang   PetscReal           alpha  = (ts->ptime - t) / ts->time_step;
354e49d4f37SHong Zhang 
355e49d4f37SHong Zhang   PetscFunctionBegin;
3569566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol));
3579566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol));
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
359e49d4f37SHong Zhang }
360e49d4f37SHong Zhang 
361d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
362d71ae5a4SJacob Faibussowitsch {
363e49d4f37SHong Zhang   PetscFunctionBegin;
364e49d4f37SHong Zhang   *yr = 1.0 + xr;
365e49d4f37SHong Zhang   *yi = xi;
3663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
367e49d4f37SHong Zhang }
368e49d4f37SHong Zhang 
369e49d4f37SHong Zhang /*@C
370e49d4f37SHong Zhang   TSBasicSymplecticSetType - Set the type of the basic symplectic method
371e49d4f37SHong Zhang 
372c3339decSBarry Smith   Logically Collective
373e49d4f37SHong Zhang 
374d8d19677SJose E. Roman   Input Parameters:
375e49d4f37SHong Zhang +  ts - timestepping context
376e49d4f37SHong Zhang -  bsymptype - type of the symplectic scheme
377e49d4f37SHong Zhang 
378bcf0153eSBarry Smith   Options Database Key:
379147403d9SBarry Smith .  -ts_basicsymplectic_type <scheme> - select the scheme
380e49d4f37SHong Zhang 
381bcf0153eSBarry Smith   Level: intermediate
382bcf0153eSBarry Smith 
383bcf0153eSBarry Smith   Note:
384bcf0153eSBarry 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`
385147403d9SBarry Smith     that is intended to store the user-provided RHS function.
386e49d4f37SHong Zhang 
387bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
388e49d4f37SHong Zhang @*/
389d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype)
390d71ae5a4SJacob Faibussowitsch {
391e49d4f37SHong Zhang   PetscFunctionBegin;
392e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
393cac4c232SBarry Smith   PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype));
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
395e49d4f37SHong Zhang }
396e49d4f37SHong Zhang 
397e49d4f37SHong Zhang /*@C
398e49d4f37SHong Zhang   TSBasicSymplecticGetType - Get the type of the basic symplectic method
399e49d4f37SHong Zhang 
400c3339decSBarry Smith   Logically Collective
401e49d4f37SHong Zhang 
402d8d19677SJose E. Roman   Input Parameters:
403e49d4f37SHong Zhang +  ts - timestepping context
404e49d4f37SHong Zhang -  bsymptype - type of the basic symplectic scheme
405e49d4f37SHong Zhang 
406e49d4f37SHong Zhang   Level: intermediate
407bcf0153eSBarry Smith 
408bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
409e49d4f37SHong Zhang @*/
410d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype)
411d71ae5a4SJacob Faibussowitsch {
412e49d4f37SHong Zhang   PetscFunctionBegin;
413e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
414cac4c232SBarry Smith   PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype));
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
416e49d4f37SHong Zhang }
417e49d4f37SHong Zhang 
418d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype)
419d71ae5a4SJacob Faibussowitsch {
420e49d4f37SHong Zhang   TS_BasicSymplectic       *bsymp = (TS_BasicSymplectic *)ts->data;
421e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
422e49d4f37SHong Zhang   PetscBool                 match;
423e49d4f37SHong Zhang 
424e49d4f37SHong Zhang   PetscFunctionBegin;
425e49d4f37SHong Zhang   if (bsymp->scheme) {
4269566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match));
4273ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
428e49d4f37SHong Zhang   }
429e49d4f37SHong Zhang   for (link = BasicSymplecticSchemeList; link; link = link->next) {
4309566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match));
431e49d4f37SHong Zhang     if (match) {
432e49d4f37SHong Zhang       bsymp->scheme = &link->sch;
4333ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
434e49d4f37SHong Zhang     }
435e49d4f37SHong Zhang   }
43698921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", bsymptype);
437e49d4f37SHong Zhang }
438e49d4f37SHong Zhang 
439d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts, TSBasicSymplecticType *bsymptype)
440d71ae5a4SJacob Faibussowitsch {
441e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
442e49d4f37SHong Zhang 
443e49d4f37SHong Zhang   PetscFunctionBegin;
444e49d4f37SHong Zhang   *bsymptype = bsymp->scheme->name;
4453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
446e49d4f37SHong Zhang }
447e49d4f37SHong Zhang 
448e49d4f37SHong Zhang /*MC
449bcf0153eSBarry Smith   TSBASICSYMPLECTIC - ODE solver using basic symplectic integration schemes
450e49d4f37SHong Zhang 
451bcf0153eSBarry Smith   These methods are intended for separable Hamiltonian systems
452bcf0153eSBarry Smith .vb
453bcf0153eSBarry Smith   qdot = dH(q,p,t)/dp
454bcf0153eSBarry Smith   pdot = -dH(q,p,t)/dq
455bcf0153eSBarry Smith .ve
456e49d4f37SHong Zhang 
457e49d4f37SHong Zhang   where the Hamiltonian can be split into the sum of kinetic energy and potential energy
458bcf0153eSBarry Smith .vb
459bcf0153eSBarry Smith   H(q,p,t) = T(p,t) + V(q,t).
460bcf0153eSBarry Smith .ve
461e49d4f37SHong Zhang 
462e49d4f37SHong Zhang   As a result, the system can be genearlly represented by
463bcf0153eSBarry Smith .vb
464bcf0153eSBarry Smith   qdot = f(p,t) = dT(p,t)/dp
465bcf0153eSBarry Smith   pdot = g(q,t) = -dV(q,t)/dq
466bcf0153eSBarry Smith .ve
467e49d4f37SHong Zhang 
468e49d4f37SHong Zhang   and solved iteratively with
469bcf0153eSBarry Smith .vb
470bcf0153eSBarry Smith   q_new = q_old + d_i*h*f(p_old,t_old)
471bcf0153eSBarry Smith   t_new = t_old + d_i*h
472bcf0153eSBarry Smith   p_new = p_old + c_i*h*g(p_new,t_new)
473bcf0153eSBarry Smith   i=0,1,...,n.
474bcf0153eSBarry Smith .ve
475e49d4f37SHong Zhang 
476bcf0153eSBarry Smith   The solution vector should contain both q and p, which correspond to (generalized) position and momentum respectively. Note that the momentum component
477bcf0153eSBarry 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".
478bcf0153eSBarry Smith   Each split is associated with an `IS` object and a sub-`TS` that is intended to store the user-provided RHS function.
479e49d4f37SHong Zhang 
480e49d4f37SHong Zhang   Level: beginner
481e49d4f37SHong Zhang 
482bcf0153eSBarry Smith   Reference:
483bcf0153eSBarry Smith . * -  wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator)
484e49d4f37SHong Zhang 
485bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()`, `TSType`
486e49d4f37SHong Zhang M*/
487d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts)
488d71ae5a4SJacob Faibussowitsch {
489e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp;
490e49d4f37SHong Zhang 
491e49d4f37SHong Zhang   PetscFunctionBegin;
4929566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticInitializePackage());
4934dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&bsymp));
494e49d4f37SHong Zhang   ts->data = (void *)bsymp;
495e49d4f37SHong Zhang 
496e49d4f37SHong Zhang   ts->ops->setup           = TSSetUp_BasicSymplectic;
497e49d4f37SHong Zhang   ts->ops->step            = TSStep_BasicSymplectic;
498e49d4f37SHong Zhang   ts->ops->reset           = TSReset_BasicSymplectic;
499e49d4f37SHong Zhang   ts->ops->destroy         = TSDestroy_BasicSymplectic;
500e49d4f37SHong Zhang   ts->ops->setfromoptions  = TSSetFromOptions_BasicSymplectic;
501e49d4f37SHong Zhang   ts->ops->view            = TSView_BasicSymplectic;
502e49d4f37SHong Zhang   ts->ops->interpolate     = TSInterpolate_BasicSymplectic;
503e49d4f37SHong Zhang   ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic;
504e49d4f37SHong Zhang 
5059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", TSBasicSymplecticSetType_BasicSymplectic));
5069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", TSBasicSymplecticGetType_BasicSymplectic));
507e49d4f37SHong Zhang 
5089566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticSetType(ts, TSBasicSymplecticDefault));
5093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
510e49d4f37SHong Zhang }
511