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 35 Level: intermediate 36 37 .seealso: [](ch_ts), `TSBASICSYMPLECTIC` 38 M*/ 39 40 /*MC 41 TSBASICSYMPLECTICVELVERLET - second order Velocity Verlet method (leapfrog method with starting process and determining velocity and position at the same time) 42 43 Level: intermediate 44 45 .seealso: [](ch_ts), `TSBASICSYMPLECTIC` 46 M*/ 47 48 /*@C 49 TSBasicSymplecticRegisterAll - Registers all of the basic symplectic integration methods in `TSBASICSYMPLECTIC` 50 51 Not Collective, but should be called by all processes which will need the schemes to be registered 52 53 Level: advanced 54 55 .seealso: [](ch_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 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 PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC4, 4, 4, c, d)); 78 } 79 PetscFunctionReturn(PETSC_SUCCESS); 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 .seealso: [](ch_ts), `TSBasicSymplecticRegister()`, `TSBasicSymplecticRegisterAll()`, `TSBASICSYMPLECTIC` 90 @*/ 91 PetscErrorCode TSBasicSymplecticRegisterDestroy(void) 92 { 93 BasicSymplecticSchemeLink link; 94 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 } 106 107 /*@C 108 TSBasicSymplecticInitializePackage - This function initializes everything in the `TSBASICSYMPLECTIC` package. It is called 109 from `TSInitializePackage()`. 110 111 Level: developer 112 113 .seealso: [](ch_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 } 124 125 /*@C 126 TSBasicSymplecticFinalizePackage - This function destroys everything in the `TSBASICSYMPLECTIC` package. It is 127 called from `PetscFinalize()`. 128 129 Level: developer 130 131 .seealso: [](ch_ts), `PetscFinalize()`, `TSBASICSYMPLECTIC` 132 @*/ 133 PetscErrorCode TSBasicSymplecticFinalizePackage(void) 134 { 135 PetscFunctionBegin; 136 TSBasicSymplecticPackageInitialized = PETSC_FALSE; 137 PetscCall(TSBasicSymplecticRegisterDestroy()); 138 PetscFunctionReturn(PETSC_SUCCESS); 139 } 140 141 /*@C 142 TSBasicSymplecticRegister - register a basic symplectic integration scheme by providing the coefficients. 143 144 Not Collective, but the same schemes should be registered on all processes on which they will be used 145 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) 152 153 Level: advanced 154 155 Notes: 156 Several symplectic methods are provided, this function is only needed to create new methods. 157 158 .seealso: [](ch_ts), `TSBASICSYMPLECTIC` 159 @*/ 160 PetscErrorCode TSBasicSymplecticRegister(TSRosWType name, PetscInt order, PetscInt s, PetscReal c[], PetscReal d[]) 161 { 162 BasicSymplecticSchemeLink link; 163 BasicSymplecticScheme scheme; 164 165 PetscFunctionBegin; 166 PetscValidCharPointer(name, 1); 167 PetscValidRealPointer(c, 4); 168 PetscValidRealPointer(d, 5); 169 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 } 183 184 /* 185 The simplified form of the equations are: 186 187 .vb 188 p_{i+1} = p_i + c_i*g(q_i)*h 189 q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h 190 .ve 191 192 Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p. 193 194 To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps: 195 .vb 196 - Update the velocity of the particle by adding to it its acceleration multiplied by c_1 197 - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1 198 - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2 199 - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2 200 .ve 201 202 */ 203 static PetscErrorCode TSStep_BasicSymplectic(TS ts) 204 { 205 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 206 BasicSymplecticScheme scheme = bsymp->scheme; 207 Vec solution = ts->vec_sol, update = bsymp->update, q, p, q_update, p_update; 208 IS is_q = bsymp->is_q, is_p = bsymp->is_p; 209 TS subts_q = bsymp->subts_q, subts_p = bsymp->subts_p; 210 PetscBool stageok; 211 PetscReal next_time_step = ts->time_step; 212 PetscInt iter; 213 214 PetscFunctionBegin; 215 PetscCall(VecGetSubVector(solution, is_q, &q)); 216 PetscCall(VecGetSubVector(solution, is_p, &p)); 217 PetscCall(VecGetSubVector(update, is_q, &q_update)); 218 PetscCall(VecGetSubVector(update, is_p, &p_update)); 219 220 for (iter = 0; iter < scheme->s; iter++) { 221 PetscCall(TSPreStage(ts, ts->ptime)); 222 /* update velocity p */ 223 if (scheme->c[iter]) { 224 PetscCall(TSComputeRHSFunction(subts_p, ts->ptime, q, p_update)); 225 PetscCall(VecAXPY(p, scheme->c[iter] * ts->time_step, p_update)); 226 } 227 /* update position q */ 228 if (scheme->d[iter]) { 229 PetscCall(TSComputeRHSFunction(subts_q, ts->ptime, p, q_update)); 230 PetscCall(VecAXPY(q, scheme->d[iter] * ts->time_step, q_update)); 231 ts->ptime = ts->ptime + scheme->d[iter] * ts->time_step; 232 } 233 PetscCall(TSPostStage(ts, ts->ptime, 0, &solution)); 234 PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok)); 235 if (!stageok) { 236 ts->reason = TS_DIVERGED_STEP_REJECTED; 237 PetscFunctionReturn(PETSC_SUCCESS); 238 } 239 PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok)); 240 if (!stageok) { 241 ts->reason = TS_DIVERGED_STEP_REJECTED; 242 PetscFunctionReturn(PETSC_SUCCESS); 243 } 244 } 245 246 ts->time_step = next_time_step; 247 PetscCall(VecRestoreSubVector(solution, is_q, &q)); 248 PetscCall(VecRestoreSubVector(solution, is_p, &p)); 249 PetscCall(VecRestoreSubVector(update, is_q, &q_update)); 250 PetscCall(VecRestoreSubVector(update, is_p, &p_update)); 251 PetscFunctionReturn(PETSC_SUCCESS); 252 } 253 254 static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, void *ctx) 255 { 256 PetscFunctionBegin; 257 PetscFunctionReturn(PETSC_SUCCESS); 258 } 259 260 static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 261 { 262 PetscFunctionBegin; 263 PetscFunctionReturn(PETSC_SUCCESS); 264 } 265 266 static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, void *ctx) 267 { 268 PetscFunctionBegin; 269 PetscFunctionReturn(PETSC_SUCCESS); 270 } 271 272 static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 273 { 274 PetscFunctionBegin; 275 PetscFunctionReturn(PETSC_SUCCESS); 276 } 277 278 static PetscErrorCode TSSetUp_BasicSymplectic(TS ts) 279 { 280 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 281 DM dm; 282 283 PetscFunctionBegin; 284 PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q)); 285 PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p)); 286 PetscCheck(bsymp->is_q && bsymp->is_p, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up RHSSplits with TSRHSSplitSetIS() using split names position and momentum respectively in order to use -ts_type basicsymplectic"); 287 PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q)); 288 PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p)); 289 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"); 290 291 PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update)); 292 293 PetscCall(TSGetAdapt(ts, &ts->adapt)); 294 PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */ 295 PetscCall(TSGetDM(ts, &dm)); 296 if (dm) { 297 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts)); 298 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts)); 299 } 300 PetscFunctionReturn(PETSC_SUCCESS); 301 } 302 303 static PetscErrorCode TSReset_BasicSymplectic(TS ts) 304 { 305 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 306 307 PetscFunctionBegin; 308 PetscCall(VecDestroy(&bsymp->update)); 309 PetscFunctionReturn(PETSC_SUCCESS); 310 } 311 312 static PetscErrorCode TSDestroy_BasicSymplectic(TS ts) 313 { 314 PetscFunctionBegin; 315 PetscCall(TSReset_BasicSymplectic(ts)); 316 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL)); 317 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL)); 318 PetscCall(PetscFree(ts->data)); 319 PetscFunctionReturn(PETSC_SUCCESS); 320 } 321 322 static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems *PetscOptionsObject) 323 { 324 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 325 326 PetscFunctionBegin; 327 PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options"); 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 ; 336 PetscCall(PetscMalloc1(count, (char ***)&namelist)); 337 for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) namelist[count] = link->sch.name; 338 PetscCall(PetscOptionsEList("-ts_basicsymplectic_type", "Family of basic symplectic integration method", "TSBasicSymplecticSetType", (const char *const *)namelist, count, bsymp->scheme->name, &choice, &flg)); 339 if (flg) PetscCall(TSBasicSymplecticSetType(ts, namelist[choice])); 340 PetscCall(PetscFree(namelist)); 341 } 342 PetscOptionsHeadEnd(); 343 PetscFunctionReturn(PETSC_SUCCESS); 344 } 345 346 static PetscErrorCode TSView_BasicSymplectic(TS ts, PetscViewer viewer) 347 { 348 PetscFunctionBegin; 349 PetscFunctionReturn(PETSC_SUCCESS); 350 } 351 352 static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X) 353 { 354 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 355 Vec update = bsymp->update; 356 PetscReal alpha = (ts->ptime - t) / ts->time_step; 357 358 PetscFunctionBegin; 359 PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol)); 360 PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol)); 361 PetscFunctionReturn(PETSC_SUCCESS); 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(PETSC_SUCCESS); 370 } 371 372 /*@C 373 TSBasicSymplecticSetType - Set the type of the basic symplectic method 374 375 Logically Collective 376 377 Input Parameters: 378 + ts - timestepping context 379 - bsymptype - type of the symplectic scheme 380 381 Options Database Key: 382 . -ts_basicsymplectic_type <scheme> - select the scheme 383 384 Level: intermediate 385 386 Note: 387 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` 388 that is intended to store the user-provided RHS function. 389 390 .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()` 391 @*/ 392 PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype) 393 { 394 PetscFunctionBegin; 395 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 396 PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype)); 397 PetscFunctionReturn(PETSC_SUCCESS); 398 } 399 400 /*@C 401 TSBasicSymplecticGetType - Get the type of the basic symplectic method 402 403 Logically Collective 404 405 Input Parameters: 406 + ts - timestepping context 407 - bsymptype - type of the basic symplectic scheme 408 409 Level: intermediate 410 411 .seealso: [](ch_ts), `TSBASICSYMPLECTIC`, `TSBasicSymplecticType`, `TSBasicSymplecticSetType()` 412 @*/ 413 PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype) 414 { 415 PetscFunctionBegin; 416 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 417 PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype)); 418 PetscFunctionReturn(PETSC_SUCCESS); 419 } 420 421 static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype) 422 { 423 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 424 BasicSymplecticSchemeLink link; 425 PetscBool match; 426 427 PetscFunctionBegin; 428 if (bsymp->scheme) { 429 PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match)); 430 if (match) PetscFunctionReturn(PETSC_SUCCESS); 431 } 432 for (link = BasicSymplecticSchemeList; link; link = link->next) { 433 PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match)); 434 if (match) { 435 bsymp->scheme = &link->sch; 436 PetscFunctionReturn(PETSC_SUCCESS); 437 } 438 } 439 SETERRQ(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(PETSC_SUCCESS); 449 } 450 451 /*MC 452 TSBASICSYMPLECTIC - ODE solver using basic symplectic integration schemes 453 454 These methods are intended for separable Hamiltonian systems 455 .vb 456 qdot = dH(q,p,t)/dp 457 pdot = -dH(q,p,t)/dq 458 .ve 459 460 where the Hamiltonian can be split into the sum of kinetic energy and potential energy 461 .vb 462 H(q,p,t) = T(p,t) + V(q,t). 463 .ve 464 465 As a result, the system can be genearlly represented by 466 .vb 467 qdot = f(p,t) = dT(p,t)/dp 468 pdot = g(q,t) = -dV(q,t)/dq 469 .ve 470 471 and solved iteratively with 472 .vb 473 q_new = q_old + d_i*h*f(p_old,t_old) 474 t_new = t_old + d_i*h 475 p_new = p_old + c_i*h*g(p_new,t_new) 476 i=0,1,...,n. 477 .ve 478 479 The solution vector should contain both q and p, which correspond to (generalized) position and momentum respectively. Note that the momentum component 480 could simply be velocity in some representations. The symplectic solver always expects a two-way splitting with the split names being "position" and "momentum". 481 Each split is associated with an `IS` object and a sub-`TS` that is intended to store the user-provided RHS function. 482 483 Level: beginner 484 485 Reference: 486 . * - wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator) 487 488 .seealso: [](ch_ts), `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()`, `TSType` 489 M*/ 490 PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts) 491 { 492 TS_BasicSymplectic *bsymp; 493 494 PetscFunctionBegin; 495 PetscCall(TSBasicSymplecticInitializePackage()); 496 PetscCall(PetscNew(&bsymp)); 497 ts->data = (void *)bsymp; 498 499 ts->ops->setup = TSSetUp_BasicSymplectic; 500 ts->ops->step = TSStep_BasicSymplectic; 501 ts->ops->reset = TSReset_BasicSymplectic; 502 ts->ops->destroy = TSDestroy_BasicSymplectic; 503 ts->ops->setfromoptions = TSSetFromOptions_BasicSymplectic; 504 ts->ops->view = TSView_BasicSymplectic; 505 ts->ops->interpolate = TSInterpolate_BasicSymplectic; 506 ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic; 507 508 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", TSBasicSymplecticSetType_BasicSymplectic)); 509 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", TSBasicSymplecticGetType_BasicSymplectic)); 510 511 PetscCall(TSBasicSymplecticSetType(ts, TSBasicSymplecticDefault)); 512 PetscFunctionReturn(PETSC_SUCCESS); 513 } 514