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