xref: /petsc/src/ts/impls/symplectic/basicsymplectic/basicsymplectic.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
371cc06b55SBarry 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 
451cc06b55SBarry 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 
551cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegisterDestroy()`
56e49d4f37SHong Zhang @*/
TSBasicSymplecticRegisterAll(void)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 
891cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticRegister()`, `TSBasicSymplecticRegisterAll()`, `TSBASICSYMPLECTIC`
90e49d4f37SHong Zhang @*/
TSBasicSymplecticRegisterDestroy(void)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 
1131cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSBASICSYMPLECTIC`
114e49d4f37SHong Zhang @*/
TSBasicSymplecticInitializePackage(void)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 
1311cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSBASICSYMPLECTIC`
132e49d4f37SHong Zhang @*/
TSBasicSymplecticFinalizePackage(void)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 
1551d27aa22SBarry Smith   Note:
156e49d4f37SHong Zhang   Several symplectic methods are provided, this function is only needed to create new methods.
157e49d4f37SHong Zhang 
1581cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`
159e49d4f37SHong Zhang @*/
TSBasicSymplecticRegister(TSRosWType name,PetscInt order,PetscInt s,PetscReal c[],PetscReal d[])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;
1664f572ea9SToby Isaac   PetscAssertPointer(name, 1);
1674f572ea9SToby Isaac   PetscAssertPointer(c, 4);
1684f572ea9SToby Isaac   PetscAssertPointer(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
18825de42d4SHong Zhang  q_{i+1} = q_i + c_i*g(p_i)*h
18925de42d4SHong Zhang  p_{i+1} = p_i + d_i*f(q_{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
19625de42d4SHong Zhang - Update the position of the particle by adding to it its velocity multiplied by c_1
19725de42d4SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by d_1
19825de42d4SHong Zhang - Update the position of the particle by adding to it its (updated) velocity multiplied by c_2
19925de42d4SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by d_2
20037fdd005SBarry Smith .ve
201e49d4f37SHong Zhang 
202e49d4f37SHong Zhang */
TSStep_BasicSymplectic(TS ts)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;
210f9813ea2SStefano Zampini   PetscBool             stageok = PETSC_TRUE;
211f023dd21SMatthew G. Knepley   PetscReal             ptime, next_time_step = ts->time_step;
212f023dd21SMatthew G. Knepley   PetscInt              n;
213e49d4f37SHong Zhang 
214e49d4f37SHong Zhang   PetscFunctionBegin;
215f023dd21SMatthew G. Knepley   PetscCall(TSGetStepNumber(ts, &n));
216f023dd21SMatthew G. Knepley   PetscCall(TSSetStepNumber(subts_p, n));
217f023dd21SMatthew G. Knepley   PetscCall(TSSetStepNumber(subts_q, n));
218f023dd21SMatthew G. Knepley   PetscCall(TSGetTime(ts, &ptime));
219f023dd21SMatthew G. Knepley   PetscCall(TSSetTime(subts_p, ptime));
220f023dd21SMatthew G. Knepley   PetscCall(TSSetTime(subts_q, ptime));
2219566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(update, is_q, &q_update));
2229566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(update, is_p, &p_update));
223f023dd21SMatthew G. Knepley   for (PetscInt iter = 0; iter < scheme->s; iter++) {
224f9813ea2SStefano Zampini     PetscCall(TSPreStage(ts, ptime));
225f9813ea2SStefano Zampini     PetscCall(VecGetSubVector(solution, is_q, &q));
226f9813ea2SStefano Zampini     PetscCall(VecGetSubVector(solution, is_p, &p));
227e49d4f37SHong Zhang     /* update position q */
22825de42d4SHong Zhang     if (scheme->c[iter]) {
229f9813ea2SStefano Zampini       PetscCall(TSComputeRHSFunction(subts_q, ptime, p, q_update));
23025de42d4SHong Zhang       PetscCall(VecAXPY(q, scheme->c[iter] * ts->time_step, q_update));
23125de42d4SHong Zhang     }
23225de42d4SHong Zhang     /* update velocity p */
23325de42d4SHong Zhang     if (scheme->d[iter]) {
234f9813ea2SStefano Zampini       ptime = ptime + scheme->d[iter] * ts->time_step;
23525de42d4SHong Zhang       PetscCall(TSComputeRHSFunction(subts_p, ptime, q, p_update));
23625de42d4SHong Zhang       PetscCall(VecAXPY(p, scheme->d[iter] * ts->time_step, p_update));
237e49d4f37SHong Zhang     }
2389566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(solution, is_q, &q));
2399566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(solution, is_p, &p));
240f9813ea2SStefano Zampini     PetscCall(TSPostStage(ts, ptime, 0, &solution));
241f9813ea2SStefano Zampini     PetscCall(TSAdaptCheckStage(ts->adapt, ts, ptime, solution, &stageok));
242f9813ea2SStefano Zampini     if (!stageok) goto finally;
243f9813ea2SStefano Zampini     PetscCall(TSFunctionDomainError(ts, ptime, solution, &stageok));
244f9813ea2SStefano Zampini     if (!stageok) goto finally;
245f9813ea2SStefano Zampini   }
246f9813ea2SStefano Zampini 
247f9813ea2SStefano Zampini finally:
248f9813ea2SStefano Zampini   if (!stageok) ts->reason = TS_DIVERGED_STEP_REJECTED;
249f9813ea2SStefano Zampini   else ts->ptime += next_time_step;
2509566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(update, is_q, &q_update));
2519566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(update, is_p, &p_update));
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
253e49d4f37SHong Zhang }
254e49d4f37SHong Zhang 
DMCoarsenHook_BasicSymplectic(DM fine,DM coarse,PetscCtx ctx)255*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, PetscCtx ctx)
256d71ae5a4SJacob Faibussowitsch {
257e49d4f37SHong Zhang   PetscFunctionBegin;
2583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
259e49d4f37SHong Zhang }
260e49d4f37SHong Zhang 
DMRestrictHook_BasicSymplectic(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,PetscCtx ctx)261*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx)
262d71ae5a4SJacob Faibussowitsch {
263e49d4f37SHong Zhang   PetscFunctionBegin;
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
265e49d4f37SHong Zhang }
266e49d4f37SHong Zhang 
DMSubDomainHook_BasicSymplectic(DM dm,DM subdm,PetscCtx ctx)267*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, PetscCtx ctx)
268d71ae5a4SJacob Faibussowitsch {
269e49d4f37SHong Zhang   PetscFunctionBegin;
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
271e49d4f37SHong Zhang }
272e49d4f37SHong Zhang 
DMSubDomainRestrictHook_BasicSymplectic(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,PetscCtx ctx)273*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx)
274d71ae5a4SJacob Faibussowitsch {
275e49d4f37SHong Zhang   PetscFunctionBegin;
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
277e49d4f37SHong Zhang }
278e49d4f37SHong Zhang 
TSSetUp_BasicSymplectic(TS ts)279d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BasicSymplectic(TS ts)
280d71ae5a4SJacob Faibussowitsch {
281e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
282e49d4f37SHong Zhang   DM                  dm;
283e49d4f37SHong Zhang 
284e49d4f37SHong Zhang   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q));
2869566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p));
287aaa8cc7dSPierre 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");
2889566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q));
2899566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p));
2903c633725SBarry 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");
291e49d4f37SHong Zhang 
2929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update));
293e49d4f37SHong Zhang 
2949566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
2959566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */
2969566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
297e49d4f37SHong Zhang   if (dm) {
2989566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts));
2999566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts));
300e49d4f37SHong Zhang   }
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
302e49d4f37SHong Zhang }
303e49d4f37SHong Zhang 
TSReset_BasicSymplectic(TS ts)304d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BasicSymplectic(TS ts)
305d71ae5a4SJacob Faibussowitsch {
306e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
307e49d4f37SHong Zhang 
308e49d4f37SHong Zhang   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bsymp->update));
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311e49d4f37SHong Zhang }
312e49d4f37SHong Zhang 
TSDestroy_BasicSymplectic(TS ts)313d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BasicSymplectic(TS ts)
314d71ae5a4SJacob Faibussowitsch {
315e49d4f37SHong Zhang   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(TSReset_BasicSymplectic(ts));
3172e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL));
3182e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL));
3199566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
321e49d4f37SHong Zhang }
322e49d4f37SHong Zhang 
TSSetFromOptions_BasicSymplectic(TS ts,PetscOptionItems PetscOptionsObject)323ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems PetscOptionsObject)
324d71ae5a4SJacob Faibussowitsch {
325e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
326e49d4f37SHong Zhang 
327e49d4f37SHong Zhang   PetscFunctionBegin;
328d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options");
329e49d4f37SHong Zhang   {
330e49d4f37SHong Zhang     BasicSymplecticSchemeLink link;
331e49d4f37SHong Zhang     PetscInt                  count, choice;
332e49d4f37SHong Zhang     PetscBool                 flg;
333e49d4f37SHong Zhang     const char              **namelist;
334e49d4f37SHong Zhang 
335fbccb6d4SPierre Jolivet     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++);
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 
TSInterpolate_BasicSymplectic(TS ts,PetscReal t,Vec X)346d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X)
347d71ae5a4SJacob Faibussowitsch {
348e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp  = (TS_BasicSymplectic *)ts->data;
349e49d4f37SHong Zhang   Vec                 update = bsymp->update;
350e49d4f37SHong Zhang   PetscReal           alpha  = (ts->ptime - t) / ts->time_step;
351e49d4f37SHong Zhang 
352e49d4f37SHong Zhang   PetscFunctionBegin;
3539566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol));
3549566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol));
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
356e49d4f37SHong Zhang }
357e49d4f37SHong Zhang 
TSComputeLinearStability_BasicSymplectic(TS ts,PetscReal xr,PetscReal xi,PetscReal * yr,PetscReal * yi)358d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
359d71ae5a4SJacob Faibussowitsch {
360e49d4f37SHong Zhang   PetscFunctionBegin;
361e49d4f37SHong Zhang   *yr = 1.0 + xr;
362e49d4f37SHong Zhang   *yi = xi;
3633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
364e49d4f37SHong Zhang }
365e49d4f37SHong Zhang 
366cc4c1da9SBarry Smith /*@
367e49d4f37SHong Zhang   TSBasicSymplecticSetType - Set the type of the basic symplectic method
368e49d4f37SHong Zhang 
369c3339decSBarry Smith   Logically Collective
370e49d4f37SHong Zhang 
371d8d19677SJose E. Roman   Input Parameters:
372e49d4f37SHong Zhang + ts        - timestepping context
373e49d4f37SHong Zhang - bsymptype - type of the symplectic scheme
374e49d4f37SHong Zhang 
375bcf0153eSBarry Smith   Options Database Key:
376147403d9SBarry Smith . -ts_basicsymplectic_type <scheme> - select the scheme
377e49d4f37SHong Zhang 
378bcf0153eSBarry Smith   Level: intermediate
379bcf0153eSBarry Smith 
380bcf0153eSBarry Smith   Note:
3811d27aa22SBarry Smith   The symplectic solver always expects a two-way splitting with the split names being "position" and "momentum".
3821d27aa22SBarry Smith   Each split is associated with an `IS` object and a sub-`TS`
383147403d9SBarry Smith   that is intended to store the user-provided RHS function.
384e49d4f37SHong Zhang 
38542747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`
386e49d4f37SHong Zhang @*/
TSBasicSymplecticSetType(TS ts,TSBasicSymplecticType bsymptype)387d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype)
388d71ae5a4SJacob Faibussowitsch {
389e49d4f37SHong Zhang   PetscFunctionBegin;
390e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
391cac4c232SBarry Smith   PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype));
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
393e49d4f37SHong Zhang }
394e49d4f37SHong Zhang 
395cc4c1da9SBarry Smith /*@
396e49d4f37SHong Zhang   TSBasicSymplecticGetType - Get the type of the basic symplectic method
397e49d4f37SHong Zhang 
398c3339decSBarry Smith   Logically Collective
399e49d4f37SHong Zhang 
400d8d19677SJose E. Roman   Input Parameters:
401e49d4f37SHong Zhang + ts        - timestepping context
402e49d4f37SHong Zhang - bsymptype - type of the basic symplectic scheme
403e49d4f37SHong Zhang 
404e49d4f37SHong Zhang   Level: intermediate
405bcf0153eSBarry Smith 
4061cc06b55SBarry Smith .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
407e49d4f37SHong Zhang @*/
TSBasicSymplecticGetType(TS ts,TSBasicSymplecticType * bsymptype)408d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype)
409d71ae5a4SJacob Faibussowitsch {
410e49d4f37SHong Zhang   PetscFunctionBegin;
411e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
412cac4c232SBarry Smith   PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype));
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
414e49d4f37SHong Zhang }
415e49d4f37SHong Zhang 
TSBasicSymplecticSetType_BasicSymplectic(TS ts,TSBasicSymplecticType bsymptype)416d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype)
417d71ae5a4SJacob Faibussowitsch {
418e49d4f37SHong Zhang   TS_BasicSymplectic       *bsymp = (TS_BasicSymplectic *)ts->data;
419e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
420e49d4f37SHong Zhang   PetscBool                 match;
421e49d4f37SHong Zhang 
422e49d4f37SHong Zhang   PetscFunctionBegin;
423e49d4f37SHong Zhang   if (bsymp->scheme) {
4249566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match));
4253ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
426e49d4f37SHong Zhang   }
427e49d4f37SHong Zhang   for (link = BasicSymplecticSchemeList; link; link = link->next) {
4289566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match));
429e49d4f37SHong Zhang     if (match) {
430e49d4f37SHong Zhang       bsymp->scheme = &link->sch;
4313ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
432e49d4f37SHong Zhang     }
433e49d4f37SHong Zhang   }
43498921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", bsymptype);
435e49d4f37SHong Zhang }
436e49d4f37SHong Zhang 
TSBasicSymplecticGetType_BasicSymplectic(TS ts,TSBasicSymplecticType * bsymptype)437d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts, TSBasicSymplecticType *bsymptype)
438d71ae5a4SJacob Faibussowitsch {
439e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
440e49d4f37SHong Zhang 
441e49d4f37SHong Zhang   PetscFunctionBegin;
442e49d4f37SHong Zhang   *bsymptype = bsymp->scheme->name;
4433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
444e49d4f37SHong Zhang }
445e49d4f37SHong Zhang 
446e49d4f37SHong Zhang /*MC
4471d27aa22SBarry Smith   TSBASICSYMPLECTIC - ODE solver using basic symplectic integration schemes <https://en.wikipedia.org/wiki/Symplectic_integrator>
448e49d4f37SHong Zhang 
449bcf0153eSBarry Smith   These methods are intended for separable Hamiltonian systems
4501d27aa22SBarry Smith 
4511d27aa22SBarry Smith   $$
4521d27aa22SBarry Smith   \begin{align*}
45362e6f1d1SMatthew Knepley   \dot q &= \frac{dH(q,p,t)}{dp}   \\
45462e6f1d1SMatthew Knepley   \dot p &= -\frac{dH(q,p,t)}{dq}
4551d27aa22SBarry Smith   \end{align*}
4561d27aa22SBarry Smith   $$
457e49d4f37SHong Zhang 
458e49d4f37SHong Zhang   where the Hamiltonian can be split into the sum of kinetic energy and potential energy
4591d27aa22SBarry Smith 
4601d27aa22SBarry Smith   $$
461bcf0153eSBarry Smith   H(q,p,t) = T(p,t) + V(q,t).
4621d27aa22SBarry Smith   $$
463e49d4f37SHong Zhang 
464145b44c9SPierre Jolivet   As a result, the system can be generally represented by
4651d27aa22SBarry Smith 
4661d27aa22SBarry Smith   $$
4671d27aa22SBarry Smith   \begin{align*}
46862e6f1d1SMatthew Knepley   \dot q &= f(p,t) = \frac{dT(p,t)}{dp} \\
46962e6f1d1SMatthew Knepley   \dot p &= g(q,t) = -\frac{dV(q,t)}{dq}
4701d27aa22SBarry Smith   \end{align*}
4711d27aa22SBarry Smith   $$
472e49d4f37SHong Zhang 
47362e6f1d1SMatthew Knepley   and solved iteratively with $i \in [0, n]$
4741d27aa22SBarry Smith 
4751d27aa22SBarry Smith   $$
4761d27aa22SBarry Smith   \begin{align*}
47762e6f1d1SMatthew Knepley   q_{new} &= q_{old} + h d_i f(p_{old}, t_{old}) \\
47862e6f1d1SMatthew Knepley   t_{new} &= t_{old} + h d_i \\
47962e6f1d1SMatthew Knepley   p_{new} &= p_{old} + h c_i g(q_{new}, t_{new})
4801d27aa22SBarry Smith   \end{align*}
4811d27aa22SBarry Smith   $$
482e49d4f37SHong Zhang 
483bcf0153eSBarry Smith   The solution vector should contain both q and p, which correspond to (generalized) position and momentum respectively. Note that the momentum component
484bcf0153eSBarry 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".
485bcf0153eSBarry Smith   Each split is associated with an `IS` object and a sub-`TS` that is intended to store the user-provided RHS function.
486e49d4f37SHong Zhang 
487e49d4f37SHong Zhang   Level: beginner
488e49d4f37SHong Zhang 
4891cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()`, `TSType`
490e49d4f37SHong Zhang M*/
TSCreate_BasicSymplectic(TS ts)491d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts)
492d71ae5a4SJacob Faibussowitsch {
493e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp;
494e49d4f37SHong Zhang 
495e49d4f37SHong Zhang   PetscFunctionBegin;
4969566063dSJacob Faibussowitsch   PetscCall(TSBasicSymplecticInitializePackage());
4974dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&bsymp));
498e49d4f37SHong Zhang   ts->data = (void *)bsymp;
499e49d4f37SHong Zhang 
500e49d4f37SHong Zhang   ts->ops->setup           = TSSetUp_BasicSymplectic;
501e49d4f37SHong Zhang   ts->ops->step            = TSStep_BasicSymplectic;
502e49d4f37SHong Zhang   ts->ops->reset           = TSReset_BasicSymplectic;
503e49d4f37SHong Zhang   ts->ops->destroy         = TSDestroy_BasicSymplectic;
504e49d4f37SHong Zhang   ts->ops->setfromoptions  = TSSetFromOptions_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