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