xref: /petsc/src/ts/impls/symplectic/basicsymplectic/basicsymplectic.c (revision e49d4f379a92c3d2f3574da21f0d6f03788744a4)
1*e49d4f37SHong Zhang /*
2*e49d4f37SHong Zhang   Code for Timestepping with basic symplectic integrators for separable Hamiltonian systems
3*e49d4f37SHong Zhang */
4*e49d4f37SHong Zhang #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5*e49d4f37SHong Zhang #include <petscdm.h>
6*e49d4f37SHong Zhang 
7*e49d4f37SHong Zhang static TSBasicSymplecticType TSBasicSymplecticDefault = TSBASICSYMPLECTICSIEULER;
8*e49d4f37SHong Zhang static PetscBool TSBasicSymplecticRegisterAllCalled;
9*e49d4f37SHong Zhang static PetscBool TSBasicSymplecticPackageInitialized;
10*e49d4f37SHong Zhang 
11*e49d4f37SHong Zhang typedef struct _BasicSymplecticScheme *BasicSymplecticScheme;
12*e49d4f37SHong Zhang typedef struct _BasicSymplecticSchemeLink *BasicSymplecticSchemeLink;
13*e49d4f37SHong Zhang 
14*e49d4f37SHong Zhang struct _BasicSymplecticScheme {
15*e49d4f37SHong Zhang   char      *name;
16*e49d4f37SHong Zhang   PetscInt  order;
17*e49d4f37SHong Zhang   PetscInt  s;       /* number of stages */
18*e49d4f37SHong Zhang   PetscReal *c,*d;
19*e49d4f37SHong Zhang };
20*e49d4f37SHong Zhang struct _BasicSymplecticSchemeLink {
21*e49d4f37SHong Zhang   struct _BasicSymplecticScheme sch;
22*e49d4f37SHong Zhang   BasicSymplecticSchemeLink     next;
23*e49d4f37SHong Zhang };
24*e49d4f37SHong Zhang static BasicSymplecticSchemeLink BasicSymplecticSchemeList;
25*e49d4f37SHong Zhang typedef struct {
26*e49d4f37SHong Zhang   TS                    subts_p,subts_q; /* sub TS contexts that holds the RHSFunction pointers */
27*e49d4f37SHong Zhang   IS                    is_p,is_q; /* IS sets for position and momentum respectively */
28*e49d4f37SHong Zhang   Vec                   update;    /* a nest work vector for generalized coordinates */
29*e49d4f37SHong Zhang   BasicSymplecticScheme scheme;
30*e49d4f37SHong Zhang } TS_BasicSymplectic;
31*e49d4f37SHong Zhang 
32*e49d4f37SHong Zhang /*MC
33*e49d4f37SHong Zhang   TSBASICSYMPLECTICSIEULER - first order semi-implicit Euler method
34*e49d4f37SHong Zhang   Level: intermediate
35*e49d4f37SHong Zhang .seealso: TSBASICSYMPLECTIC
36*e49d4f37SHong Zhang M*/
37*e49d4f37SHong Zhang 
38*e49d4f37SHong Zhang /*MC
39*e49d4f37SHong Zhang   TSBASICSYMPLECTICVELVERLET - second order Velocity Verlet method (leapfrog method with starting process and determing velocity and position at the same time)
40*e49d4f37SHong Zhang Level: intermediate
41*e49d4f37SHong Zhang .seealso: TSBASICSYMPLECTIC
42*e49d4f37SHong Zhang M*/
43*e49d4f37SHong Zhang 
44*e49d4f37SHong Zhang /*@C
45*e49d4f37SHong Zhang   TSBasicSymplecticRegisterAll - Registers all of the basic symplectic integration methods in TSBasicSymplectic
46*e49d4f37SHong Zhang 
47*e49d4f37SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
48*e49d4f37SHong Zhang 
49*e49d4f37SHong Zhang   Level: advanced
50*e49d4f37SHong Zhang 
51*e49d4f37SHong Zhang .keywords: TS, TSBASICSYMPLECTIC, register, all
52*e49d4f37SHong Zhang 
53*e49d4f37SHong Zhang .seealso:  TSBasicSymplecticRegisterDestroy()
54*e49d4f37SHong Zhang @*/
55*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticRegisterAll(void)
56*e49d4f37SHong Zhang {
57*e49d4f37SHong Zhang   PetscErrorCode ierr;
58*e49d4f37SHong Zhang 
59*e49d4f37SHong Zhang   PetscFunctionBegin;
60*e49d4f37SHong Zhang   if (TSBasicSymplecticRegisterAllCalled) PetscFunctionReturn(0);
61*e49d4f37SHong Zhang   TSBasicSymplecticRegisterAllCalled = PETSC_TRUE;
62*e49d4f37SHong Zhang   {
63*e49d4f37SHong Zhang     PetscReal c[1] = {1.0},d[1] = {1.0};
64*e49d4f37SHong Zhang     ierr = TSBasicSymplecticRegister(TSBASICSYMPLECTICSIEULER,1,1,c,d);CHKERRQ(ierr);
65*e49d4f37SHong Zhang   }
66*e49d4f37SHong Zhang   {
67*e49d4f37SHong Zhang     PetscReal c[2] = {0,1.0},d[2] = {0.5,0.5};
68*e49d4f37SHong Zhang     ierr = TSBasicSymplecticRegister(TSBASICSYMPLECTICVELVERLET,2,2,c,d);CHKERRQ(ierr);
69*e49d4f37SHong Zhang   }
70*e49d4f37SHong Zhang   {
71*e49d4f37SHong 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};
72*e49d4f37SHong Zhang     ierr = TSBasicSymplecticRegister(TSBASICSYMPLECTIC3,3,3,c,d);CHKERRQ(ierr);
73*e49d4f37SHong Zhang   }
74*e49d4f37SHong Zhang   {
75*e49d4f37SHong Zhang #define CUBEROOTOFTWO 1.2599210498948731647672106
76*e49d4f37SHong 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};
77*e49d4f37SHong Zhang     ierr = TSBasicSymplecticRegister(TSBASICSYMPLECTIC4,4,4,c,d);CHKERRQ(ierr);
78*e49d4f37SHong Zhang   }
79*e49d4f37SHong Zhang   PetscFunctionReturn(0);
80*e49d4f37SHong Zhang }
81*e49d4f37SHong Zhang 
82*e49d4f37SHong Zhang /*@C
83*e49d4f37SHong Zhang    TSBasicSymplecticRegisterDestroy - Frees the list of schemes that were registered by TSBasicSymplecticRegister().
84*e49d4f37SHong Zhang 
85*e49d4f37SHong Zhang    Not Collective
86*e49d4f37SHong Zhang 
87*e49d4f37SHong Zhang    Level: advanced
88*e49d4f37SHong Zhang 
89*e49d4f37SHong Zhang .keywords: TSBasicSymplectic, register, destroy
90*e49d4f37SHong Zhang .seealso: TSBasicSymplecticRegister(), TSBasicSymplecticRegisterAll()
91*e49d4f37SHong Zhang @*/
92*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticRegisterDestroy(void)
93*e49d4f37SHong Zhang {
94*e49d4f37SHong Zhang   PetscErrorCode            ierr;
95*e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
96*e49d4f37SHong Zhang 
97*e49d4f37SHong Zhang   PetscFunctionBegin;
98*e49d4f37SHong Zhang   while ((link = BasicSymplecticSchemeList)) {
99*e49d4f37SHong Zhang     BasicSymplecticScheme scheme = &link->sch;
100*e49d4f37SHong Zhang     BasicSymplecticSchemeList = link->next;
101*e49d4f37SHong Zhang     ierr = PetscFree2(scheme->c,scheme->d);CHKERRQ(ierr);
102*e49d4f37SHong Zhang     ierr = PetscFree(scheme->name);CHKERRQ(ierr);
103*e49d4f37SHong Zhang     ierr = PetscFree(link);CHKERRQ(ierr);
104*e49d4f37SHong Zhang   }
105*e49d4f37SHong Zhang   TSBasicSymplecticRegisterAllCalled = PETSC_FALSE;
106*e49d4f37SHong Zhang   PetscFunctionReturn(0);
107*e49d4f37SHong Zhang }
108*e49d4f37SHong Zhang 
109*e49d4f37SHong Zhang /*@C
110*e49d4f37SHong Zhang   TSBasicSymplecticInitializePackage - This function initializes everything in the TSBasicSymplectic package. It is called
111*e49d4f37SHong Zhang   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_BasicSymplectic()
112*e49d4f37SHong Zhang   when using static libraries.
113*e49d4f37SHong Zhang 
114*e49d4f37SHong Zhang   Level: developer
115*e49d4f37SHong Zhang 
116*e49d4f37SHong Zhang .keywords: TS, TSBasicSymplectic, initialize, package
117*e49d4f37SHong Zhang .seealso: PetscInitialize()
118*e49d4f37SHong Zhang @*/
119*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticInitializePackage(void)
120*e49d4f37SHong Zhang {
121*e49d4f37SHong Zhang   PetscErrorCode ierr;
122*e49d4f37SHong Zhang 
123*e49d4f37SHong Zhang   PetscFunctionBegin;
124*e49d4f37SHong Zhang   if (TSBasicSymplecticPackageInitialized) PetscFunctionReturn(0);
125*e49d4f37SHong Zhang   TSBasicSymplecticPackageInitialized = PETSC_TRUE;
126*e49d4f37SHong Zhang   ierr = TSBasicSymplecticRegisterAll();CHKERRQ(ierr);
127*e49d4f37SHong Zhang   ierr = PetscRegisterFinalize(TSBasicSymplecticFinalizePackage);CHKERRQ(ierr);
128*e49d4f37SHong Zhang   PetscFunctionReturn(0);
129*e49d4f37SHong Zhang }
130*e49d4f37SHong Zhang 
131*e49d4f37SHong Zhang /*@C
132*e49d4f37SHong Zhang   TSBasicSymplecticFinalizePackage - This function destroys everything in the TSBasicSymplectic package. It is
133*e49d4f37SHong Zhang   called from PetscFinalize().
134*e49d4f37SHong Zhang 
135*e49d4f37SHong Zhang   Level: developer
136*e49d4f37SHong Zhang 
137*e49d4f37SHong Zhang .keywords: Petsc, destroy, package
138*e49d4f37SHong Zhang .seealso: PetscFinalize()
139*e49d4f37SHong Zhang @*/
140*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticFinalizePackage(void)
141*e49d4f37SHong Zhang {
142*e49d4f37SHong Zhang   PetscErrorCode ierr;
143*e49d4f37SHong Zhang 
144*e49d4f37SHong Zhang   PetscFunctionBegin;
145*e49d4f37SHong Zhang   TSBasicSymplecticPackageInitialized = PETSC_FALSE;
146*e49d4f37SHong Zhang   ierr = TSBasicSymplecticRegisterDestroy();CHKERRQ(ierr);
147*e49d4f37SHong Zhang   PetscFunctionReturn(0);
148*e49d4f37SHong Zhang }
149*e49d4f37SHong Zhang 
150*e49d4f37SHong Zhang /*@C
151*e49d4f37SHong Zhang    TSBasicSymplecticRegister - register a basic symplectic integration scheme by providing the coefficients.
152*e49d4f37SHong Zhang 
153*e49d4f37SHong Zhang    Not Collective, but the same schemes should be registered on all processes on which they will be used
154*e49d4f37SHong Zhang 
155*e49d4f37SHong Zhang    Input Parameters:
156*e49d4f37SHong Zhang +  name - identifier for method
157*e49d4f37SHong Zhang .  order - approximation order of method
158*e49d4f37SHong Zhang .  s - number of stages, this is the dimension of the matrices below
159*e49d4f37SHong Zhang .  c - coefficients for updating generalized position (dimension s)
160*e49d4f37SHong Zhang -  d - coefficients for updating generalized momentum (dimension s)
161*e49d4f37SHong Zhang 
162*e49d4f37SHong Zhang    Notes:
163*e49d4f37SHong Zhang    Several symplectic methods are provided, this function is only needed to create new methods.
164*e49d4f37SHong Zhang 
165*e49d4f37SHong Zhang    Level: advanced
166*e49d4f37SHong Zhang 
167*e49d4f37SHong Zhang .keywords: TS, register
168*e49d4f37SHong Zhang 
169*e49d4f37SHong Zhang .seealso: TSBasicSymplectic
170*e49d4f37SHong Zhang @*/
171*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticRegister(TSRosWType name,PetscInt order,PetscInt s,PetscReal c[],PetscReal d[])
172*e49d4f37SHong Zhang {
173*e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
174*e49d4f37SHong Zhang   BasicSymplecticScheme     scheme;
175*e49d4f37SHong Zhang   PetscErrorCode  ierr;
176*e49d4f37SHong Zhang 
177*e49d4f37SHong Zhang   PetscFunctionBegin;
178*e49d4f37SHong Zhang   PetscValidCharPointer(name,1);
179*e49d4f37SHong Zhang   PetscValidPointer(c,4);
180*e49d4f37SHong Zhang   PetscValidPointer(d,4);
181*e49d4f37SHong Zhang 
182*e49d4f37SHong Zhang   ierr = PetscCalloc1(1,&link);CHKERRQ(ierr);
183*e49d4f37SHong Zhang   scheme = &link->sch;
184*e49d4f37SHong Zhang   ierr = PetscStrallocpy(name,&scheme->name);CHKERRQ(ierr);
185*e49d4f37SHong Zhang   scheme->order = order;
186*e49d4f37SHong Zhang   scheme->s = s;
187*e49d4f37SHong Zhang   ierr = PetscMalloc2(s,&scheme->c,s,&scheme->d);CHKERRQ(ierr);
188*e49d4f37SHong Zhang   ierr = PetscMemcpy(scheme->c,c,s*sizeof(c[0]));CHKERRQ(ierr);
189*e49d4f37SHong Zhang   ierr = PetscMemcpy(scheme->d,d,s*sizeof(d[0]));CHKERRQ(ierr);
190*e49d4f37SHong Zhang   link->next = BasicSymplecticSchemeList;
191*e49d4f37SHong Zhang   BasicSymplecticSchemeList = link;
192*e49d4f37SHong Zhang   PetscFunctionReturn(0);
193*e49d4f37SHong Zhang }
194*e49d4f37SHong Zhang 
195*e49d4f37SHong Zhang /*
196*e49d4f37SHong Zhang The simplified form of the equations are:
197*e49d4f37SHong Zhang 
198*e49d4f37SHong Zhang $ p_{i+1} = p_i + c_i*g(q_i)*h
199*e49d4f37SHong Zhang $ q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h
200*e49d4f37SHong Zhang 
201*e49d4f37SHong Zhang Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p.
202*e49d4f37SHong Zhang 
203*e49d4f37SHong Zhang To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps:
204*e49d4f37SHong Zhang 
205*e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration multiplied by c_1
206*e49d4f37SHong Zhang - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1
207*e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2
208*e49d4f37SHong Zhang - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2
209*e49d4f37SHong Zhang 
210*e49d4f37SHong Zhang */
211*e49d4f37SHong Zhang static PetscErrorCode TSStep_BasicSymplectic(TS ts)
212*e49d4f37SHong Zhang {
213*e49d4f37SHong Zhang   TS_BasicSymplectic    *bsymp = (TS_BasicSymplectic*)ts->data;
214*e49d4f37SHong Zhang   BasicSymplecticScheme scheme = bsymp->scheme;
215*e49d4f37SHong Zhang   Vec                   solution = ts->vec_sol,update = bsymp->update,q,p,q_update,p_update;
216*e49d4f37SHong Zhang   IS                    is_q = bsymp->is_q,is_p = bsymp->is_p;
217*e49d4f37SHong Zhang   TS                    subts_q = bsymp->subts_q,subts_p = bsymp->subts_p;
218*e49d4f37SHong Zhang   PetscBool             stageok;
219*e49d4f37SHong Zhang   PetscReal             next_time_step = ts->time_step;
220*e49d4f37SHong Zhang   PetscInt              iter;
221*e49d4f37SHong Zhang   PetscErrorCode        ierr;
222*e49d4f37SHong Zhang 
223*e49d4f37SHong Zhang   PetscFunctionBegin;
224*e49d4f37SHong Zhang   ierr = VecGetSubVector(solution,is_q,&q);CHKERRQ(ierr);
225*e49d4f37SHong Zhang   ierr = VecGetSubVector(solution,is_p,&p);CHKERRQ(ierr);
226*e49d4f37SHong Zhang   ierr = VecGetSubVector(update,is_q,&q_update);CHKERRQ(ierr);
227*e49d4f37SHong Zhang   ierr = VecGetSubVector(update,is_p,&p_update);CHKERRQ(ierr);
228*e49d4f37SHong Zhang 
229*e49d4f37SHong Zhang   for (iter = 0;iter<scheme->s;iter++) {
230*e49d4f37SHong Zhang     ierr = TSPreStage(ts,ts->ptime);CHKERRQ(ierr);
231*e49d4f37SHong Zhang     /* update velocity p */
232*e49d4f37SHong Zhang     if (scheme->c[iter]) {
233*e49d4f37SHong Zhang       ierr = TSComputeRHSFunction(subts_p,ts->ptime,q,p_update);CHKERRQ(ierr);
234*e49d4f37SHong Zhang       ierr = VecAXPY(p,scheme->c[iter]*ts->time_step,p_update);CHKERRQ(ierr);
235*e49d4f37SHong Zhang     }
236*e49d4f37SHong Zhang     /* update position q */
237*e49d4f37SHong Zhang     if (scheme->d[iter]) {
238*e49d4f37SHong Zhang       ierr = TSComputeRHSFunction(subts_q,ts->ptime,p,q_update);CHKERRQ(ierr);
239*e49d4f37SHong Zhang       ierr = VecAXPY(q,scheme->d[iter]*ts->time_step,q_update);CHKERRQ(ierr);
240*e49d4f37SHong Zhang       ts->ptime = ts->ptime+scheme->d[iter]*ts->time_step;
241*e49d4f37SHong Zhang     }
242*e49d4f37SHong Zhang     ierr = TSPostStage(ts,ts->ptime,0,&solution);CHKERRQ(ierr);
243*e49d4f37SHong Zhang     ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok);CHKERRQ(ierr);
244*e49d4f37SHong Zhang     if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
245*e49d4f37SHong Zhang     ierr = TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok);CHKERRQ(ierr);
246*e49d4f37SHong Zhang     if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);}
247*e49d4f37SHong Zhang   }
248*e49d4f37SHong Zhang 
249*e49d4f37SHong Zhang   ts->time_step = next_time_step;
250*e49d4f37SHong Zhang   ierr = VecRestoreSubVector(solution,is_q,&q);CHKERRQ(ierr);
251*e49d4f37SHong Zhang   ierr = VecRestoreSubVector(solution,is_p,&p);CHKERRQ(ierr);
252*e49d4f37SHong Zhang   ierr = VecRestoreSubVector(update,is_q,&q_update);CHKERRQ(ierr);
253*e49d4f37SHong Zhang   ierr = VecRestoreSubVector(update,is_p,&p_update);CHKERRQ(ierr);
254*e49d4f37SHong Zhang   PetscFunctionReturn(0);
255*e49d4f37SHong Zhang }
256*e49d4f37SHong Zhang 
257*e49d4f37SHong Zhang static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine,DM coarse,void *ctx)
258*e49d4f37SHong Zhang {
259*e49d4f37SHong Zhang   PetscFunctionBegin;
260*e49d4f37SHong Zhang   PetscFunctionReturn(0);
261*e49d4f37SHong Zhang }
262*e49d4f37SHong Zhang 
263*e49d4f37SHong Zhang static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
264*e49d4f37SHong Zhang {
265*e49d4f37SHong Zhang   PetscFunctionBegin;
266*e49d4f37SHong Zhang   PetscFunctionReturn(0);
267*e49d4f37SHong Zhang }
268*e49d4f37SHong Zhang 
269*e49d4f37SHong Zhang static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm,DM subdm,void *ctx)
270*e49d4f37SHong Zhang {
271*e49d4f37SHong Zhang   PetscFunctionBegin;
272*e49d4f37SHong Zhang   PetscFunctionReturn(0);
273*e49d4f37SHong Zhang }
274*e49d4f37SHong Zhang 
275*e49d4f37SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
276*e49d4f37SHong Zhang {
277*e49d4f37SHong Zhang 
278*e49d4f37SHong Zhang   PetscFunctionBegin;
279*e49d4f37SHong Zhang   PetscFunctionReturn(0);
280*e49d4f37SHong Zhang }
281*e49d4f37SHong Zhang 
282*e49d4f37SHong Zhang static PetscErrorCode TSSetUp_BasicSymplectic(TS ts)
283*e49d4f37SHong Zhang {
284*e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data;
285*e49d4f37SHong Zhang   DM                 dm;
286*e49d4f37SHong Zhang   PetscErrorCode     ierr;
287*e49d4f37SHong Zhang 
288*e49d4f37SHong Zhang   PetscFunctionBegin;
289*e49d4f37SHong Zhang   ierr = TSRHSSplitGetIS(ts,"position",&bsymp->is_q);CHKERRQ(ierr);
290*e49d4f37SHong Zhang   ierr = TSRHSSplitGetIS(ts,"momentum",&bsymp->is_p);CHKERRQ(ierr);
291*e49d4f37SHong Zhang   if (!bsymp->is_q || !bsymp->is_p) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names positon and momentum respectively in order to use -ts_type basicsymplectic");
292*e49d4f37SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"position",&bsymp->subts_q);CHKERRQ(ierr);
293*e49d4f37SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"momentum",&bsymp->subts_p);CHKERRQ(ierr);
294*e49d4f37SHong Zhang   if (!bsymp->subts_q || !bsymp->subts_p) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for position and momentum using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS");
295*e49d4f37SHong Zhang 
296*e49d4f37SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&bsymp->update);CHKERRQ(ierr);
297*e49d4f37SHong Zhang 
298*e49d4f37SHong Zhang   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
299*e49d4f37SHong Zhang   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); /* make sure to use fixed time stepping */
300*e49d4f37SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
301*e49d4f37SHong Zhang   if (dm) {
302*e49d4f37SHong Zhang     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_BasicSymplectic,DMRestrictHook_BasicSymplectic,ts);CHKERRQ(ierr);
303*e49d4f37SHong Zhang     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_BasicSymplectic,DMSubDomainRestrictHook_BasicSymplectic,ts);CHKERRQ(ierr);
304*e49d4f37SHong Zhang   }
305*e49d4f37SHong Zhang   PetscFunctionReturn(0);
306*e49d4f37SHong Zhang }
307*e49d4f37SHong Zhang 
308*e49d4f37SHong Zhang static PetscErrorCode TSReset_BasicSymplectic(TS ts)
309*e49d4f37SHong Zhang {
310*e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data;
311*e49d4f37SHong Zhang   PetscErrorCode     ierr;
312*e49d4f37SHong Zhang 
313*e49d4f37SHong Zhang   PetscFunctionBegin;
314*e49d4f37SHong Zhang   ierr = VecDestroy(&bsymp->update);CHKERRQ(ierr);
315*e49d4f37SHong Zhang   PetscFunctionReturn(0);
316*e49d4f37SHong Zhang }
317*e49d4f37SHong Zhang 
318*e49d4f37SHong Zhang static PetscErrorCode TSDestroy_BasicSymplectic(TS ts)
319*e49d4f37SHong Zhang {
320*e49d4f37SHong Zhang   PetscErrorCode ierr;
321*e49d4f37SHong Zhang 
322*e49d4f37SHong Zhang   PetscFunctionBegin;
323*e49d4f37SHong Zhang   ierr = TSReset_BasicSymplectic(ts);CHKERRQ(ierr);
324*e49d4f37SHong Zhang   ierr = PetscFree(ts->data);CHKERRQ(ierr);
325*e49d4f37SHong Zhang   PetscFunctionReturn(0);
326*e49d4f37SHong Zhang }
327*e49d4f37SHong Zhang 
328*e49d4f37SHong Zhang static PetscErrorCode TSSetFromOptions_BasicSymplectic(PetscOptionItems *PetscOptionsObject,TS ts)
329*e49d4f37SHong Zhang {
330*e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data;
331*e49d4f37SHong Zhang   PetscErrorCode     ierr;
332*e49d4f37SHong Zhang 
333*e49d4f37SHong Zhang   PetscFunctionBegin;
334*e49d4f37SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"Basic symplectic integrator options");CHKERRQ(ierr);
335*e49d4f37SHong Zhang   {
336*e49d4f37SHong Zhang     BasicSymplecticSchemeLink link;
337*e49d4f37SHong Zhang     PetscInt                  count,choice;
338*e49d4f37SHong Zhang     PetscBool                 flg;
339*e49d4f37SHong Zhang     const char                **namelist;
340*e49d4f37SHong Zhang 
341*e49d4f37SHong Zhang     for (link=BasicSymplecticSchemeList,count=0; link; link=link->next,count++) ;
342*e49d4f37SHong Zhang     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
343*e49d4f37SHong Zhang     for (link=BasicSymplecticSchemeList,count=0; link; link=link->next,count++) namelist[count] = link->sch.name;
344*e49d4f37SHong Zhang     ierr = PetscOptionsEList("-ts_basicsymplectic_type","Family of basic symplectic integration method","TSBasicSymplecticSetType",(const char*const*)namelist,count,bsymp->scheme->name,&choice,&flg);CHKERRQ(ierr);
345*e49d4f37SHong Zhang     if (flg) {ierr = TSBasicSymplecticSetType(ts,namelist[choice]);CHKERRQ(ierr);}
346*e49d4f37SHong Zhang     ierr = PetscFree(namelist);CHKERRQ(ierr);
347*e49d4f37SHong Zhang   }
348*e49d4f37SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
349*e49d4f37SHong Zhang   PetscFunctionReturn(0);
350*e49d4f37SHong Zhang }
351*e49d4f37SHong Zhang 
352*e49d4f37SHong Zhang static PetscErrorCode TSView_BasicSymplectic(TS ts,PetscViewer viewer)
353*e49d4f37SHong Zhang {
354*e49d4f37SHong Zhang   PetscFunctionBegin;
355*e49d4f37SHong Zhang   PetscFunctionReturn(0);
356*e49d4f37SHong Zhang }
357*e49d4f37SHong Zhang 
358*e49d4f37SHong Zhang static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts,PetscReal t,Vec X)
359*e49d4f37SHong Zhang {
360*e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data;
361*e49d4f37SHong Zhang   Vec                update = bsymp->update;
362*e49d4f37SHong Zhang   PetscReal          alpha = (ts->ptime - t)/ts->time_step;
363*e49d4f37SHong Zhang   PetscErrorCode     ierr;
364*e49d4f37SHong Zhang 
365*e49d4f37SHong Zhang   PetscFunctionBegin;
366*e49d4f37SHong Zhang   ierr = VecWAXPY(X,-ts->time_step,update,ts->vec_sol);CHKERRQ(ierr);
367*e49d4f37SHong Zhang   ierr = VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol);CHKERRQ(ierr);
368*e49d4f37SHong Zhang   PetscFunctionReturn(0);
369*e49d4f37SHong Zhang }
370*e49d4f37SHong Zhang 
371*e49d4f37SHong Zhang static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
372*e49d4f37SHong Zhang {
373*e49d4f37SHong Zhang   PetscFunctionBegin;
374*e49d4f37SHong Zhang   *yr = 1.0 + xr;
375*e49d4f37SHong Zhang   *yi = xi;
376*e49d4f37SHong Zhang   PetscFunctionReturn(0);
377*e49d4f37SHong Zhang }
378*e49d4f37SHong Zhang 
379*e49d4f37SHong Zhang /*@C
380*e49d4f37SHong Zhang   TSBasicSymplecticSetType - Set the type of the basic symplectic method
381*e49d4f37SHong Zhang 
382*e49d4f37SHong Zhang   Logically Collective on TS
383*e49d4f37SHong Zhang 
384*e49d4f37SHong Zhang   Input Parameter:
385*e49d4f37SHong Zhang +  ts - timestepping context
386*e49d4f37SHong Zhang -  bsymptype - type of the symplectic scheme
387*e49d4f37SHong Zhang 
388*e49d4f37SHong Zhang   Options Database:
389*e49d4f37SHong Zhang .  -ts_basicsymplectic_type <scheme>
390*e49d4f37SHong Zhang 
391*e49d4f37SHong Zhang   Notes:
392*e49d4f37SHong Zhang   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 that is intended to store the user-provided RHS function.
393*e49d4f37SHong Zhang 
394*e49d4f37SHong Zhang   Level: intermediate
395*e49d4f37SHong Zhang @*/
396*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticSetType(TS ts,TSBasicSymplecticType bsymptype)
397*e49d4f37SHong Zhang {
398*e49d4f37SHong Zhang   PetscErrorCode ierr;
399*e49d4f37SHong Zhang 
400*e49d4f37SHong Zhang   PetscFunctionBegin;
401*e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
402*e49d4f37SHong Zhang   ierr = PetscTryMethod(ts,"TSBasicSymplecticSetType_C",(TS,TSBasicSymplecticType),(ts,bsymptype));CHKERRQ(ierr);
403*e49d4f37SHong Zhang   PetscFunctionReturn(0);
404*e49d4f37SHong Zhang }
405*e49d4f37SHong Zhang 
406*e49d4f37SHong Zhang /*@C
407*e49d4f37SHong Zhang   TSBasicSymplecticGetType - Get the type of the basic symplectic method
408*e49d4f37SHong Zhang 
409*e49d4f37SHong Zhang   Logically Collective on TS
410*e49d4f37SHong Zhang 
411*e49d4f37SHong Zhang   Input Parameter:
412*e49d4f37SHong Zhang +  ts - timestepping context
413*e49d4f37SHong Zhang -  bsymptype - type of the basic symplectic scheme
414*e49d4f37SHong Zhang 
415*e49d4f37SHong Zhang   Level: intermediate
416*e49d4f37SHong Zhang @*/
417*e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticGetType(TS ts,TSBasicSymplecticType *bsymptype)
418*e49d4f37SHong Zhang {
419*e49d4f37SHong Zhang   PetscErrorCode ierr;
420*e49d4f37SHong Zhang 
421*e49d4f37SHong Zhang   PetscFunctionBegin;
422*e49d4f37SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
423*e49d4f37SHong Zhang   ierr = PetscUseMethod(ts,"TSBasicSymplecticGetType_C",(TS,TSBasicSymplecticType*),(ts,bsymptype));CHKERRQ(ierr);
424*e49d4f37SHong Zhang   PetscFunctionReturn(0);
425*e49d4f37SHong Zhang }
426*e49d4f37SHong Zhang 
427*e49d4f37SHong Zhang static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts,TSBasicSymplecticType bsymptype)
428*e49d4f37SHong Zhang {
429*e49d4f37SHong Zhang   TS_BasicSymplectic        *bsymp = (TS_BasicSymplectic*)ts->data;
430*e49d4f37SHong Zhang   BasicSymplecticSchemeLink link;
431*e49d4f37SHong Zhang   PetscBool                 match;
432*e49d4f37SHong Zhang   PetscErrorCode            ierr;
433*e49d4f37SHong Zhang 
434*e49d4f37SHong Zhang   PetscFunctionBegin;
435*e49d4f37SHong Zhang   if (bsymp->scheme) {
436*e49d4f37SHong Zhang     ierr = PetscStrcmp(bsymp->scheme->name,bsymptype,&match);CHKERRQ(ierr);
437*e49d4f37SHong Zhang     if (match) PetscFunctionReturn(0);
438*e49d4f37SHong Zhang   }
439*e49d4f37SHong Zhang   for (link = BasicSymplecticSchemeList; link; link=link->next) {
440*e49d4f37SHong Zhang     ierr = PetscStrcmp(link->sch.name,bsymptype,&match);CHKERRQ(ierr);
441*e49d4f37SHong Zhang     if (match) {
442*e49d4f37SHong Zhang       bsymp->scheme = &link->sch;
443*e49d4f37SHong Zhang       PetscFunctionReturn(0);
444*e49d4f37SHong Zhang     }
445*e49d4f37SHong Zhang   }
446*e49d4f37SHong Zhang   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",bsymptype);
447*e49d4f37SHong Zhang   PetscFunctionReturn(0);
448*e49d4f37SHong Zhang }
449*e49d4f37SHong Zhang 
450*e49d4f37SHong Zhang static PetscErrorCode  TSBasicSymplecticGetType_BasicSymplectic(TS ts,TSBasicSymplecticType *bsymptype)
451*e49d4f37SHong Zhang {
452*e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data;
453*e49d4f37SHong Zhang 
454*e49d4f37SHong Zhang   PetscFunctionBegin;
455*e49d4f37SHong Zhang   *bsymptype = bsymp->scheme->name;
456*e49d4f37SHong Zhang   PetscFunctionReturn(0);
457*e49d4f37SHong Zhang }
458*e49d4f37SHong Zhang 
459*e49d4f37SHong Zhang /*MC
460*e49d4f37SHong Zhang   TSBasicSymplectic - ODE solver using basic symplectic integration schemes
461*e49d4f37SHong Zhang 
462*e49d4f37SHong Zhang   These methods are intened for separable Hamiltonian systems
463*e49d4f37SHong Zhang 
464*e49d4f37SHong Zhang $  qdot = dH(q,p,t)/dp
465*e49d4f37SHong Zhang $  pdot = -dH(q,p,t)/dq
466*e49d4f37SHong Zhang 
467*e49d4f37SHong Zhang   where the Hamiltonian can be split into the sum of kinetic energy and potential energy
468*e49d4f37SHong Zhang 
469*e49d4f37SHong Zhang $  H(q,p,t) = T(p,t) + V(q,t).
470*e49d4f37SHong Zhang 
471*e49d4f37SHong Zhang   As a result, the system can be genearlly represented by
472*e49d4f37SHong Zhang 
473*e49d4f37SHong Zhang $  qdot = f(p,t) = dT(p,t)/dp
474*e49d4f37SHong Zhang $  pdot = g(q,t) = -dV(q,t)/dq
475*e49d4f37SHong Zhang 
476*e49d4f37SHong Zhang   and solved iteratively with
477*e49d4f37SHong Zhang 
478*e49d4f37SHong Zhang $  q_new = q_old + d_i*h*f(p_old,t_old)
479*e49d4f37SHong Zhang $  t_new = t_old + d_i*h
480*e49d4f37SHong Zhang $  p_new = p_old + c_i*h*g(p_new,t_new)
481*e49d4f37SHong Zhang $  i=0,1,...,n.
482*e49d4f37SHong Zhang 
483*e49d4f37SHong Zhang   The solution vector should contain both q and p, which correspond to (generalized) position and momentum respectively. Note that the momentum component could simply be velocity in some representations.
484*e49d4f37SHong Zhang   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 that is intended to store the user-provided RHS function.
485*e49d4f37SHong Zhang 
486*e49d4f37SHong Zhang   Reference: wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator)
487*e49d4f37SHong Zhang 
488*e49d4f37SHong Zhang   Level: beginner
489*e49d4f37SHong Zhang 
490*e49d4f37SHong Zhang .seealso:  TSCreate(), TSSetType(), TSRHSSplitSetIS(), TSRHSSplitSetRHSFunction()
491*e49d4f37SHong Zhang 
492*e49d4f37SHong Zhang M*/
493*e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts)
494*e49d4f37SHong Zhang {
495*e49d4f37SHong Zhang   TS_BasicSymplectic *bsymp;
496*e49d4f37SHong Zhang   PetscErrorCode     ierr;
497*e49d4f37SHong Zhang 
498*e49d4f37SHong Zhang   PetscFunctionBegin;
499*e49d4f37SHong Zhang   ierr = TSBasicSymplecticInitializePackage();CHKERRQ(ierr);
500*e49d4f37SHong Zhang   ierr = PetscNewLog(ts,&bsymp);CHKERRQ(ierr);
501*e49d4f37SHong Zhang   ts->data = (void*)bsymp;
502*e49d4f37SHong Zhang 
503*e49d4f37SHong Zhang   ts->ops->setup           = TSSetUp_BasicSymplectic;
504*e49d4f37SHong Zhang   ts->ops->step            = TSStep_BasicSymplectic;
505*e49d4f37SHong Zhang   ts->ops->reset           = TSReset_BasicSymplectic;
506*e49d4f37SHong Zhang   ts->ops->destroy         = TSDestroy_BasicSymplectic;
507*e49d4f37SHong Zhang   ts->ops->setfromoptions  = TSSetFromOptions_BasicSymplectic;
508*e49d4f37SHong Zhang   ts->ops->view            = TSView_BasicSymplectic;
509*e49d4f37SHong Zhang   ts->ops->interpolate     = TSInterpolate_BasicSymplectic;
510*e49d4f37SHong Zhang   ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic;
511*e49d4f37SHong Zhang 
512*e49d4f37SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticSetType_C",TSBasicSymplecticSetType_BasicSymplectic);CHKERRQ(ierr);
513*e49d4f37SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticGetType_C",TSBasicSymplecticGetType_BasicSymplectic);CHKERRQ(ierr);
514*e49d4f37SHong Zhang 
515*e49d4f37SHong Zhang   ierr = TSBasicSymplecticSetType(ts,TSBasicSymplecticDefault);CHKERRQ(ierr);
516*e49d4f37SHong Zhang   PetscFunctionReturn(0);
517*e49d4f37SHong Zhang }
518