Actual source code: basicsymplectic.c

  1: /*
  2:   Code for Timestepping with basic symplectic integrators for separable Hamiltonian systems
  3: */
  4: #include <petsc/private/tsimpl.h>
  5: #include <petscdm.h>

  7: static TSBasicSymplecticType TSBasicSymplecticDefault = TSBASICSYMPLECTICSIEULER;
  8: static PetscBool             TSBasicSymplecticRegisterAllCalled;
  9: static PetscBool             TSBasicSymplecticPackageInitialized;

 11: typedef struct _BasicSymplecticScheme     *BasicSymplecticScheme;
 12: typedef struct _BasicSymplecticSchemeLink *BasicSymplecticSchemeLink;

 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;

 32: /*MC
 33:   TSBASICSYMPLECTICSIEULER - first order semi-implicit Euler method

 35:   Level: intermediate

 37: .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`
 38: M*/

 40: /*MC
 41:   TSBASICSYMPLECTICVELVERLET - second order Velocity Verlet method (leapfrog method with starting process and determing velocity and position at the same time)

 43: Level: intermediate

 45: .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`
 46: M*/

 48: /*@C
 49:   TSBasicSymplecticRegisterAll - Registers all of the basic symplectic integration methods in `TSBASICSYMPLECTIC`

 51:   Not Collective, but should be called by all processes which will need the schemes to be registered

 53:   Level: advanced

 55: .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegisterDestroy()`
 56: @*/
 57: PetscErrorCode TSBasicSymplecticRegisterAll(void)
 58: {
 59:   PetscFunctionBegin;
 60:   if (TSBasicSymplecticRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
 61:   TSBasicSymplecticRegisterAllCalled = PETSC_TRUE;
 62:   {
 63:     PetscReal c[1] = {1.0}, d[1] = {1.0};
 64:     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICSIEULER, 1, 1, c, d));
 65:   }
 66:   {
 67:     PetscReal c[2] = {0, 1.0}, d[2] = {0.5, 0.5};
 68:     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICVELVERLET, 2, 2, c, d));
 69:   }
 70:   {
 71:     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:     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC3, 3, 3, c, d));
 73:   }
 74:   {
 75: #define CUBE../../../../..OFTWO 1.2599210498948731647672106
 76:     PetscReal c[4] = {1.0 / 2.0 / (2.0 - CUBE../../../../..OFTWO), (1.0 - CUBE../../../../..OFTWO) / 2.0 / (2.0 - CUBE../../../../..OFTWO), (1.0 - CUBE../../../../..OFTWO) / 2.0 / (2.0 - CUBE../../../../..OFTWO), 1.0 / 2.0 / (2.0 - CUBE../../../../..OFTWO)}, d[4] = {1.0 / (2.0 - CUBE../../../../..OFTWO), -CUBE../../../../..OFTWO / (2.0 - CUBE../../../../..OFTWO), 1.0 / (2.0 - CUBE../../../../..OFTWO), 0};
 77:     PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC4, 4, 4, c, d));
 78:   }
 79:   PetscFunctionReturn(PETSC_SUCCESS);
 80: }

 82: /*@C
 83:    TSBasicSymplecticRegisterDestroy - Frees the list of schemes that were registered by `TSBasicSymplecticRegister()`.

 85:    Not Collective

 87:    Level: advanced

 89: .seealso: [](chapter_ts), `TSBasicSymplecticRegister()`, `TSBasicSymplecticRegisterAll()`, `TSBASICSYMPLECTIC`
 90: @*/
 91: PetscErrorCode TSBasicSymplecticRegisterDestroy(void)
 92: {
 93:   BasicSymplecticSchemeLink link;

 95:   PetscFunctionBegin;
 96:   while ((link = BasicSymplecticSchemeList)) {
 97:     BasicSymplecticScheme scheme = &link->sch;
 98:     BasicSymplecticSchemeList    = link->next;
 99:     PetscCall(PetscFree2(scheme->c, scheme->d));
100:     PetscCall(PetscFree(scheme->name));
101:     PetscCall(PetscFree(link));
102:   }
103:   TSBasicSymplecticRegisterAllCalled = PETSC_FALSE;
104:   PetscFunctionReturn(PETSC_SUCCESS);
105: }

107: /*@C
108:   TSBasicSymplecticInitializePackage - This function initializes everything in the `TSBASICSYMPLECTIC` package. It is called
109:   from `TSInitializePackage()`.

111:   Level: developer

113: .seealso: [](chapter_ts), `PetscInitialize()`, `TSBASICSYMPLECTIC`
114: @*/
115: PetscErrorCode TSBasicSymplecticInitializePackage(void)
116: {
117:   PetscFunctionBegin;
118:   if (TSBasicSymplecticPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
119:   TSBasicSymplecticPackageInitialized = PETSC_TRUE;
120:   PetscCall(TSBasicSymplecticRegisterAll());
121:   PetscCall(PetscRegisterFinalize(TSBasicSymplecticFinalizePackage));
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }

125: /*@C
126:   TSBasicSymplecticFinalizePackage - This function destroys everything in the `TSBASICSYMPLECTIC` package. It is
127:   called from `PetscFinalize()`.

129:   Level: developer

131: .seealso: [](chapter_ts), `PetscFinalize()`, `TSBASICSYMPLECTIC`
132: @*/
133: PetscErrorCode TSBasicSymplecticFinalizePackage(void)
134: {
135:   PetscFunctionBegin;
136:   TSBasicSymplecticPackageInitialized = PETSC_FALSE;
137:   PetscCall(TSBasicSymplecticRegisterDestroy());
138:   PetscFunctionReturn(PETSC_SUCCESS);
139: }

141: /*@C
142:    TSBasicSymplecticRegister - register a basic symplectic integration scheme by providing the coefficients.

144:    Not Collective, but the same schemes should be registered on all processes on which they will be used

146:    Input Parameters:
147: +  name - identifier for method
148: .  order - approximation order of method
149: .  s - number of stages, this is the dimension of the matrices below
150: .  c - coefficients for updating generalized position (dimension s)
151: -  d - coefficients for updating generalized momentum (dimension s)

153:    Level: advanced

155:    Notes:
156:    Several symplectic methods are provided, this function is only needed to create new methods.

158: .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`
159: @*/
160: PetscErrorCode TSBasicSymplecticRegister(TSRosWType name, PetscInt order, PetscInt s, PetscReal c[], PetscReal d[])
161: {
162:   BasicSymplecticSchemeLink link;
163:   BasicSymplecticScheme     scheme;

165:   PetscFunctionBegin;

170:   PetscCall(TSBasicSymplecticInitializePackage());
171:   PetscCall(PetscNew(&link));
172:   scheme = &link->sch;
173:   PetscCall(PetscStrallocpy(name, &scheme->name));
174:   scheme->order = order;
175:   scheme->s     = s;
176:   PetscCall(PetscMalloc2(s, &scheme->c, s, &scheme->d));
177:   PetscCall(PetscArraycpy(scheme->c, c, s));
178:   PetscCall(PetscArraycpy(scheme->d, d, s));
179:   link->next                = BasicSymplecticSchemeList;
180:   BasicSymplecticSchemeList = link;
181:   PetscFunctionReturn(PETSC_SUCCESS);
182: }

184: /*
185: The simplified form of the equations are:

187: $ p_{i+1} = p_i + c_i*g(q_i)*h
188: $ q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h

190: Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p.

192: To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps:

194: - Update the velocity of the particle by adding to it its acceleration multiplied by c_1
195: - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1
196: - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2
197: - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2

199: */
200: static PetscErrorCode TSStep_BasicSymplectic(TS ts)
201: {
202:   TS_BasicSymplectic   *bsymp    = (TS_BasicSymplectic *)ts->data;
203:   BasicSymplecticScheme scheme   = bsymp->scheme;
204:   Vec                   solution = ts->vec_sol, update = bsymp->update, q, p, q_update, p_update;
205:   IS                    is_q = bsymp->is_q, is_p = bsymp->is_p;
206:   TS                    subts_q = bsymp->subts_q, subts_p = bsymp->subts_p;
207:   PetscBool             stageok;
208:   PetscReal             next_time_step = ts->time_step;
209:   PetscInt              iter;

211:   PetscFunctionBegin;
212:   PetscCall(VecGetSubVector(solution, is_q, &q));
213:   PetscCall(VecGetSubVector(solution, is_p, &p));
214:   PetscCall(VecGetSubVector(update, is_q, &q_update));
215:   PetscCall(VecGetSubVector(update, is_p, &p_update));

217:   for (iter = 0; iter < scheme->s; iter++) {
218:     PetscCall(TSPreStage(ts, ts->ptime));
219:     /* update velocity p */
220:     if (scheme->c[iter]) {
221:       PetscCall(TSComputeRHSFunction(subts_p, ts->ptime, q, p_update));
222:       PetscCall(VecAXPY(p, scheme->c[iter] * ts->time_step, p_update));
223:     }
224:     /* update position q */
225:     if (scheme->d[iter]) {
226:       PetscCall(TSComputeRHSFunction(subts_q, ts->ptime, p, q_update));
227:       PetscCall(VecAXPY(q, scheme->d[iter] * ts->time_step, q_update));
228:       ts->ptime = ts->ptime + scheme->d[iter] * ts->time_step;
229:     }
230:     PetscCall(TSPostStage(ts, ts->ptime, 0, &solution));
231:     PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok));
232:     if (!stageok) {
233:       ts->reason = TS_DIVERGED_STEP_REJECTED;
234:       PetscFunctionReturn(PETSC_SUCCESS);
235:     }
236:     PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok));
237:     if (!stageok) {
238:       ts->reason = TS_DIVERGED_STEP_REJECTED;
239:       PetscFunctionReturn(PETSC_SUCCESS);
240:     }
241:   }

243:   ts->time_step = next_time_step;
244:   PetscCall(VecRestoreSubVector(solution, is_q, &q));
245:   PetscCall(VecRestoreSubVector(solution, is_p, &p));
246:   PetscCall(VecRestoreSubVector(update, is_q, &q_update));
247:   PetscCall(VecRestoreSubVector(update, is_p, &p_update));
248:   PetscFunctionReturn(PETSC_SUCCESS);
249: }

251: static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, void *ctx)
252: {
253:   PetscFunctionBegin;
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

257: static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
258: {
259:   PetscFunctionBegin;
260:   PetscFunctionReturn(PETSC_SUCCESS);
261: }

263: static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, void *ctx)
264: {
265:   PetscFunctionBegin;
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
270: {
271:   PetscFunctionBegin;
272:   PetscFunctionReturn(PETSC_SUCCESS);
273: }

275: static PetscErrorCode TSSetUp_BasicSymplectic(TS ts)
276: {
277:   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;
278:   DM                  dm;

280:   PetscFunctionBegin;
281:   PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q));
282:   PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p));
283:   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");
284:   PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q));
285:   PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p));
286:   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");

288:   PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update));

290:   PetscCall(TSGetAdapt(ts, &ts->adapt));
291:   PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */
292:   PetscCall(TSGetDM(ts, &dm));
293:   if (dm) {
294:     PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts));
295:     PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts));
296:   }
297:   PetscFunctionReturn(PETSC_SUCCESS);
298: }

300: static PetscErrorCode TSReset_BasicSymplectic(TS ts)
301: {
302:   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;

304:   PetscFunctionBegin;
305:   PetscCall(VecDestroy(&bsymp->update));
306:   PetscFunctionReturn(PETSC_SUCCESS);
307: }

309: static PetscErrorCode TSDestroy_BasicSymplectic(TS ts)
310: {
311:   PetscFunctionBegin;
312:   PetscCall(TSReset_BasicSymplectic(ts));
313:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL));
314:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL));
315:   PetscCall(PetscFree(ts->data));
316:   PetscFunctionReturn(PETSC_SUCCESS);
317: }

319: static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems *PetscOptionsObject)
320: {
321:   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;

323:   PetscFunctionBegin;
324:   PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options");
325:   {
326:     BasicSymplecticSchemeLink link;
327:     PetscInt                  count, choice;
328:     PetscBool                 flg;
329:     const char              **namelist;

331:     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++)
332:       ;
333:     PetscCall(PetscMalloc1(count, (char ***)&namelist));
334:     for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) namelist[count] = link->sch.name;
335:     PetscCall(PetscOptionsEList("-ts_basicsymplectic_type", "Family of basic symplectic integration method", "TSBasicSymplecticSetType", (const char *const *)namelist, count, bsymp->scheme->name, &choice, &flg));
336:     if (flg) PetscCall(TSBasicSymplecticSetType(ts, namelist[choice]));
337:     PetscCall(PetscFree(namelist));
338:   }
339:   PetscOptionsHeadEnd();
340:   PetscFunctionReturn(PETSC_SUCCESS);
341: }

343: static PetscErrorCode TSView_BasicSymplectic(TS ts, PetscViewer viewer)
344: {
345:   PetscFunctionBegin;
346:   PetscFunctionReturn(PETSC_SUCCESS);
347: }

349: static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X)
350: {
351:   TS_BasicSymplectic *bsymp  = (TS_BasicSymplectic *)ts->data;
352:   Vec                 update = bsymp->update;
353:   PetscReal           alpha  = (ts->ptime - t) / ts->time_step;

355:   PetscFunctionBegin;
356:   PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol));
357:   PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol));
358:   PetscFunctionReturn(PETSC_SUCCESS);
359: }

361: static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
362: {
363:   PetscFunctionBegin;
364:   *yr = 1.0 + xr;
365:   *yi = xi;
366:   PetscFunctionReturn(PETSC_SUCCESS);
367: }

369: /*@C
370:   TSBasicSymplecticSetType - Set the type of the basic symplectic method

372:   Logically Collective

374:   Input Parameters:
375: +  ts - timestepping context
376: -  bsymptype - type of the symplectic scheme

378:   Options Database Key:
379: .  -ts_basicsymplectic_type <scheme> - select the scheme

381:   Level: intermediate

383:   Note:
384:     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`
385:     that is intended to store the user-provided RHS function.

387: .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
388: @*/
389: PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype)
390: {
391:   PetscFunctionBegin;
393:   PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype));
394:   PetscFunctionReturn(PETSC_SUCCESS);
395: }

397: /*@C
398:   TSBasicSymplecticGetType - Get the type of the basic symplectic method

400:   Logically Collective

402:   Input Parameters:
403: +  ts - timestepping context
404: -  bsymptype - type of the basic symplectic scheme

406:   Level: intermediate

408: .seealso: [](chapter_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()`
409: @*/
410: PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype)
411: {
412:   PetscFunctionBegin;
414:   PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype));
415:   PetscFunctionReturn(PETSC_SUCCESS);
416: }

418: static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype)
419: {
420:   TS_BasicSymplectic       *bsymp = (TS_BasicSymplectic *)ts->data;
421:   BasicSymplecticSchemeLink link;
422:   PetscBool                 match;

424:   PetscFunctionBegin;
425:   if (bsymp->scheme) {
426:     PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match));
427:     if (match) PetscFunctionReturn(PETSC_SUCCESS);
428:   }
429:   for (link = BasicSymplecticSchemeList; link; link = link->next) {
430:     PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match));
431:     if (match) {
432:       bsymp->scheme = &link->sch;
433:       PetscFunctionReturn(PETSC_SUCCESS);
434:     }
435:   }
436:   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", bsymptype);
437: }

439: static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts, TSBasicSymplecticType *bsymptype)
440: {
441:   TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data;

443:   PetscFunctionBegin;
444:   *bsymptype = bsymp->scheme->name;
445:   PetscFunctionReturn(PETSC_SUCCESS);
446: }

448: /*MC
449:   TSBASICSYMPLECTIC - ODE solver using basic symplectic integration schemes

451:   These methods are intended for separable Hamiltonian systems
452: .vb
453:   qdot = dH(q,p,t)/dp
454:   pdot = -dH(q,p,t)/dq
455: .ve

457:   where the Hamiltonian can be split into the sum of kinetic energy and potential energy
458: .vb
459:   H(q,p,t) = T(p,t) + V(q,t).
460: .ve

462:   As a result, the system can be genearlly represented by
463: .vb
464:   qdot = f(p,t) = dT(p,t)/dp
465:   pdot = g(q,t) = -dV(q,t)/dq
466: .ve

468:   and solved iteratively with
469: .vb
470:   q_new = q_old + d_i*h*f(p_old,t_old)
471:   t_new = t_old + d_i*h
472:   p_new = p_old + c_i*h*g(p_new,t_new)
473:   i=0,1,...,n.
474: .ve

476:   The solution vector should contain both q and p, which correspond to (generalized) position and momentum respectively. Note that the momentum component
477:   could simply be velocity in some representations. The symplectic solver always expects a two-way splitting with the split names being "position" and "momentum".
478:   Each split is associated with an `IS` object and a sub-`TS` that is intended to store the user-provided RHS function.

480:   Level: beginner

482:   Reference:
483: . * -  wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator)

485: .seealso: [](chapter_ts), `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()`, `TSType`
486: M*/
487: PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts)
488: {
489:   TS_BasicSymplectic *bsymp;

491:   PetscFunctionBegin;
492:   PetscCall(TSBasicSymplecticInitializePackage());
493:   PetscCall(PetscNew(&bsymp));
494:   ts->data = (void *)bsymp;

496:   ts->ops->setup           = TSSetUp_BasicSymplectic;
497:   ts->ops->step            = TSStep_BasicSymplectic;
498:   ts->ops->reset           = TSReset_BasicSymplectic;
499:   ts->ops->destroy         = TSDestroy_BasicSymplectic;
500:   ts->ops->setfromoptions  = TSSetFromOptions_BasicSymplectic;
501:   ts->ops->view            = TSView_BasicSymplectic;
502:   ts->ops->interpolate     = TSInterpolate_BasicSymplectic;
503:   ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic;

505:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", TSBasicSymplecticSetType_BasicSymplectic));
506:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", TSBasicSymplecticGetType_BasicSymplectic));

508:   PetscCall(TSBasicSymplecticSetType(ts, TSBasicSymplecticDefault));
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }