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