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