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 CHKERRQ(TSBasicSymplecticRegister(TSBASICSYMPLECTICSIEULER,1,1,c,d)); 61 } 62 { 63 PetscReal c[2] = {0,1.0},d[2] = {0.5,0.5}; 64 CHKERRQ(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 CHKERRQ(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 CHKERRQ(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 CHKERRQ(PetscFree2(scheme->c,scheme->d)); 96 CHKERRQ(PetscFree(scheme->name)); 97 CHKERRQ(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 CHKERRQ(TSBasicSymplecticRegisterAll()); 117 CHKERRQ(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 CHKERRQ(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 PetscValidPointer(c,4); 164 PetscValidPointer(d,5); 165 166 CHKERRQ(TSBasicSymplecticInitializePackage()); 167 CHKERRQ(PetscNew(&link)); 168 scheme = &link->sch; 169 CHKERRQ(PetscStrallocpy(name,&scheme->name)); 170 scheme->order = order; 171 scheme->s = s; 172 CHKERRQ(PetscMalloc2(s,&scheme->c,s,&scheme->d)); 173 CHKERRQ(PetscArraycpy(scheme->c,c,s)); 174 CHKERRQ(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 CHKERRQ(VecGetSubVector(solution,is_q,&q)); 209 CHKERRQ(VecGetSubVector(solution,is_p,&p)); 210 CHKERRQ(VecGetSubVector(update,is_q,&q_update)); 211 CHKERRQ(VecGetSubVector(update,is_p,&p_update)); 212 213 for (iter = 0;iter<scheme->s;iter++) { 214 CHKERRQ(TSPreStage(ts,ts->ptime)); 215 /* update velocity p */ 216 if (scheme->c[iter]) { 217 CHKERRQ(TSComputeRHSFunction(subts_p,ts->ptime,q,p_update)); 218 CHKERRQ(VecAXPY(p,scheme->c[iter]*ts->time_step,p_update)); 219 } 220 /* update position q */ 221 if (scheme->d[iter]) { 222 CHKERRQ(TSComputeRHSFunction(subts_q,ts->ptime,p,q_update)); 223 CHKERRQ(VecAXPY(q,scheme->d[iter]*ts->time_step,q_update)); 224 ts->ptime = ts->ptime+scheme->d[iter]*ts->time_step; 225 } 226 CHKERRQ(TSPostStage(ts,ts->ptime,0,&solution)); 227 CHKERRQ(TSAdaptCheckStage(ts->adapt,ts,ts->ptime,solution,&stageok)); 228 if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 229 CHKERRQ(TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok)); 230 if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 231 } 232 233 ts->time_step = next_time_step; 234 CHKERRQ(VecRestoreSubVector(solution,is_q,&q)); 235 CHKERRQ(VecRestoreSubVector(solution,is_p,&p)); 236 CHKERRQ(VecRestoreSubVector(update,is_q,&q_update)); 237 CHKERRQ(VecRestoreSubVector(update,is_p,&p_update)); 238 PetscFunctionReturn(0); 239 } 240 241 static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine,DM coarse,void *ctx) 242 { 243 PetscFunctionBegin; 244 PetscFunctionReturn(0); 245 } 246 247 static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 248 { 249 PetscFunctionBegin; 250 PetscFunctionReturn(0); 251 } 252 253 static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm,DM subdm,void *ctx) 254 { 255 PetscFunctionBegin; 256 PetscFunctionReturn(0); 257 } 258 259 static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 260 { 261 PetscFunctionBegin; 262 PetscFunctionReturn(0); 263 } 264 265 static PetscErrorCode TSSetUp_BasicSymplectic(TS ts) 266 { 267 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 268 DM dm; 269 270 PetscFunctionBegin; 271 CHKERRQ(TSRHSSplitGetIS(ts,"position",&bsymp->is_q)); 272 CHKERRQ(TSRHSSplitGetIS(ts,"momentum",&bsymp->is_p)); 273 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"); 274 CHKERRQ(TSRHSSplitGetSubTS(ts,"position",&bsymp->subts_q)); 275 CHKERRQ(TSRHSSplitGetSubTS(ts,"momentum",&bsymp->subts_p)); 276 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"); 277 278 CHKERRQ(VecDuplicate(ts->vec_sol,&bsymp->update)); 279 280 CHKERRQ(TSGetAdapt(ts,&ts->adapt)); 281 CHKERRQ(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */ 282 CHKERRQ(TSGetDM(ts,&dm)); 283 if (dm) { 284 CHKERRQ(DMCoarsenHookAdd(dm,DMCoarsenHook_BasicSymplectic,DMRestrictHook_BasicSymplectic,ts)); 285 CHKERRQ(DMSubDomainHookAdd(dm,DMSubDomainHook_BasicSymplectic,DMSubDomainRestrictHook_BasicSymplectic,ts)); 286 } 287 PetscFunctionReturn(0); 288 } 289 290 static PetscErrorCode TSReset_BasicSymplectic(TS ts) 291 { 292 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 293 294 PetscFunctionBegin; 295 CHKERRQ(VecDestroy(&bsymp->update)); 296 PetscFunctionReturn(0); 297 } 298 299 static PetscErrorCode TSDestroy_BasicSymplectic(TS ts) 300 { 301 PetscFunctionBegin; 302 CHKERRQ(TSReset_BasicSymplectic(ts)); 303 CHKERRQ(PetscFree(ts->data)); 304 PetscFunctionReturn(0); 305 } 306 307 static PetscErrorCode TSSetFromOptions_BasicSymplectic(PetscOptionItems *PetscOptionsObject,TS ts) 308 { 309 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 310 311 PetscFunctionBegin; 312 CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Basic symplectic integrator options")); 313 { 314 BasicSymplecticSchemeLink link; 315 PetscInt count,choice; 316 PetscBool flg; 317 const char **namelist; 318 319 for (link=BasicSymplecticSchemeList,count=0; link; link=link->next,count++) ; 320 CHKERRQ(PetscMalloc1(count,(char***)&namelist)); 321 for (link=BasicSymplecticSchemeList,count=0; link; link=link->next,count++) namelist[count] = link->sch.name; 322 CHKERRQ(PetscOptionsEList("-ts_basicsymplectic_type","Family of basic symplectic integration method","TSBasicSymplecticSetType",(const char*const*)namelist,count,bsymp->scheme->name,&choice,&flg)); 323 if (flg) CHKERRQ(TSBasicSymplecticSetType(ts,namelist[choice])); 324 CHKERRQ(PetscFree(namelist)); 325 } 326 CHKERRQ(PetscOptionsTail()); 327 PetscFunctionReturn(0); 328 } 329 330 static PetscErrorCode TSView_BasicSymplectic(TS ts,PetscViewer viewer) 331 { 332 PetscFunctionBegin; 333 PetscFunctionReturn(0); 334 } 335 336 static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts,PetscReal t,Vec X) 337 { 338 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 339 Vec update = bsymp->update; 340 PetscReal alpha = (ts->ptime - t)/ts->time_step; 341 342 PetscFunctionBegin; 343 CHKERRQ(VecWAXPY(X,-ts->time_step,update,ts->vec_sol)); 344 CHKERRQ(VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol)); 345 PetscFunctionReturn(0); 346 } 347 348 static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 349 { 350 PetscFunctionBegin; 351 *yr = 1.0 + xr; 352 *yi = xi; 353 PetscFunctionReturn(0); 354 } 355 356 /*@C 357 TSBasicSymplecticSetType - Set the type of the basic symplectic method 358 359 Logically Collective on TS 360 361 Input Parameters: 362 + ts - timestepping context 363 - bsymptype - type of the symplectic scheme 364 365 Options Database: 366 . -ts_basicsymplectic_type <scheme> - select the scheme 367 368 Notes: 369 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 370 that is intended to store the user-provided RHS function. 371 372 Level: intermediate 373 @*/ 374 PetscErrorCode TSBasicSymplecticSetType(TS ts,TSBasicSymplecticType bsymptype) 375 { 376 PetscFunctionBegin; 377 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 378 CHKERRQ(PetscTryMethod(ts,"TSBasicSymplecticSetType_C",(TS,TSBasicSymplecticType),(ts,bsymptype))); 379 PetscFunctionReturn(0); 380 } 381 382 /*@C 383 TSBasicSymplecticGetType - Get the type of the basic symplectic method 384 385 Logically Collective on TS 386 387 Input Parameters: 388 + ts - timestepping context 389 - bsymptype - type of the basic symplectic scheme 390 391 Level: intermediate 392 @*/ 393 PetscErrorCode TSBasicSymplecticGetType(TS ts,TSBasicSymplecticType *bsymptype) 394 { 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 397 CHKERRQ(PetscUseMethod(ts,"TSBasicSymplecticGetType_C",(TS,TSBasicSymplecticType*),(ts,bsymptype))); 398 PetscFunctionReturn(0); 399 } 400 401 static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts,TSBasicSymplecticType bsymptype) 402 { 403 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 404 BasicSymplecticSchemeLink link; 405 PetscBool match; 406 407 PetscFunctionBegin; 408 if (bsymp->scheme) { 409 CHKERRQ(PetscStrcmp(bsymp->scheme->name,bsymptype,&match)); 410 if (match) PetscFunctionReturn(0); 411 } 412 for (link = BasicSymplecticSchemeList; link; link=link->next) { 413 CHKERRQ(PetscStrcmp(link->sch.name,bsymptype,&match)); 414 if (match) { 415 bsymp->scheme = &link->sch; 416 PetscFunctionReturn(0); 417 } 418 } 419 SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",bsymptype); 420 } 421 422 static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts,TSBasicSymplecticType *bsymptype) 423 { 424 TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 425 426 PetscFunctionBegin; 427 *bsymptype = bsymp->scheme->name; 428 PetscFunctionReturn(0); 429 } 430 431 /*MC 432 TSBasicSymplectic - ODE solver using basic symplectic integration schemes 433 434 These methods are intened for separable Hamiltonian systems 435 436 $ qdot = dH(q,p,t)/dp 437 $ pdot = -dH(q,p,t)/dq 438 439 where the Hamiltonian can be split into the sum of kinetic energy and potential energy 440 441 $ H(q,p,t) = T(p,t) + V(q,t). 442 443 As a result, the system can be genearlly represented by 444 445 $ qdot = f(p,t) = dT(p,t)/dp 446 $ pdot = g(q,t) = -dV(q,t)/dq 447 448 and solved iteratively with 449 450 $ q_new = q_old + d_i*h*f(p_old,t_old) 451 $ t_new = t_old + d_i*h 452 $ p_new = p_old + c_i*h*g(p_new,t_new) 453 $ i=0,1,...,n. 454 455 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. 456 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. 457 458 Reference: wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator) 459 460 Level: beginner 461 462 .seealso: TSCreate(), TSSetType(), TSRHSSplitSetIS(), TSRHSSplitSetRHSFunction() 463 464 M*/ 465 PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts) 466 { 467 TS_BasicSymplectic *bsymp; 468 469 PetscFunctionBegin; 470 CHKERRQ(TSBasicSymplecticInitializePackage()); 471 CHKERRQ(PetscNewLog(ts,&bsymp)); 472 ts->data = (void*)bsymp; 473 474 ts->ops->setup = TSSetUp_BasicSymplectic; 475 ts->ops->step = TSStep_BasicSymplectic; 476 ts->ops->reset = TSReset_BasicSymplectic; 477 ts->ops->destroy = TSDestroy_BasicSymplectic; 478 ts->ops->setfromoptions = TSSetFromOptions_BasicSymplectic; 479 ts->ops->view = TSView_BasicSymplectic; 480 ts->ops->interpolate = TSInterpolate_BasicSymplectic; 481 ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic; 482 483 CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticSetType_C",TSBasicSymplecticSetType_BasicSymplectic)); 484 CHKERRQ(PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticGetType_C",TSBasicSymplecticGetType_BasicSymplectic)); 485 486 CHKERRQ(TSBasicSymplecticSetType(ts,TSBasicSymplecticDefault)); 487 PetscFunctionReturn(0); 488 } 489