1e49d4f37SHong Zhang /* 2e49d4f37SHong Zhang Code for Timestepping with basic symplectic integrators for separable Hamiltonian systems 3e49d4f37SHong Zhang */ 4e49d4f37SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5e49d4f37SHong Zhang #include <petscdm.h> 6e49d4f37SHong Zhang 7e49d4f37SHong Zhang static TSBasicSymplecticType TSBasicSymplecticDefault = TSBASICSYMPLECTICSIEULER; 8e49d4f37SHong Zhang static PetscBool TSBasicSymplecticRegisterAllCalled; 9e49d4f37SHong Zhang static PetscBool TSBasicSymplecticPackageInitialized; 10e49d4f37SHong Zhang 11e49d4f37SHong Zhang typedef struct _BasicSymplecticScheme *BasicSymplecticScheme; 12e49d4f37SHong Zhang typedef struct _BasicSymplecticSchemeLink *BasicSymplecticSchemeLink; 13e49d4f37SHong Zhang 14e49d4f37SHong Zhang struct _BasicSymplecticScheme { 15e49d4f37SHong Zhang char *name; 16e49d4f37SHong Zhang PetscInt order; 17e49d4f37SHong Zhang PetscInt s; /* number of stages */ 18e49d4f37SHong Zhang PetscReal *c, *d; 19e49d4f37SHong Zhang }; 20e49d4f37SHong Zhang struct _BasicSymplecticSchemeLink { 21e49d4f37SHong Zhang struct _BasicSymplecticScheme sch; 22e49d4f37SHong Zhang BasicSymplecticSchemeLink next; 23e49d4f37SHong Zhang }; 24e49d4f37SHong Zhang static BasicSymplecticSchemeLink BasicSymplecticSchemeList; 25e49d4f37SHong Zhang typedef struct { 26e49d4f37SHong Zhang TS subts_p, subts_q; /* sub TS contexts that holds the RHSFunction pointers */ 27e49d4f37SHong Zhang IS is_p, is_q; /* IS sets for position and momentum respectively */ 28e49d4f37SHong Zhang Vec update; /* a nest work vector for generalized coordinates */ 29e49d4f37SHong Zhang BasicSymplecticScheme scheme; 30e49d4f37SHong Zhang } TS_BasicSymplectic; 31e49d4f37SHong Zhang 32e49d4f37SHong Zhang /*MC 33e49d4f37SHong Zhang TSBASICSYMPLECTICSIEULER - first order semi-implicit Euler method 34e49d4f37SHong Zhang Level: intermediate 35db781477SPatrick Sanan .seealso: `TSBASICSYMPLECTIC` 36e49d4f37SHong Zhang M*/ 37e49d4f37SHong Zhang 38e49d4f37SHong Zhang /*MC 39e49d4f37SHong Zhang TSBASICSYMPLECTICVELVERLET - second order Velocity Verlet method (leapfrog method with starting process and determing velocity and position at the same time) 40e49d4f37SHong Zhang Level: intermediate 41db781477SPatrick Sanan .seealso: `TSBASICSYMPLECTIC` 42e49d4f37SHong Zhang M*/ 43e49d4f37SHong Zhang 44e49d4f37SHong Zhang /*@C 45e49d4f37SHong Zhang TSBasicSymplecticRegisterAll - Registers all of the basic symplectic integration methods in TSBasicSymplectic 46e49d4f37SHong Zhang 47e49d4f37SHong Zhang Not Collective, but should be called by all processes which will need the schemes to be registered 48e49d4f37SHong Zhang 49e49d4f37SHong Zhang Level: advanced 50e49d4f37SHong Zhang 51db781477SPatrick Sanan .seealso: `TSBasicSymplecticRegisterDestroy()` 52e49d4f37SHong Zhang @*/ 53*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticRegisterAll(void) 54*d71ae5a4SJacob Faibussowitsch { 55e49d4f37SHong Zhang PetscFunctionBegin; 56e49d4f37SHong Zhang if (TSBasicSymplecticRegisterAllCalled) PetscFunctionReturn(0); 57e49d4f37SHong Zhang TSBasicSymplecticRegisterAllCalled = PETSC_TRUE; 58e49d4f37SHong Zhang { 59e49d4f37SHong Zhang PetscReal c[1] = {1.0}, d[1] = {1.0}; 609566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICSIEULER, 1, 1, c, d)); 61e49d4f37SHong Zhang } 62e49d4f37SHong Zhang { 63e49d4f37SHong Zhang PetscReal c[2] = {0, 1.0}, d[2] = {0.5, 0.5}; 649566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICVELVERLET, 2, 2, c, d)); 65e49d4f37SHong Zhang } 66e49d4f37SHong Zhang { 67e49d4f37SHong Zhang 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}; 689566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC3, 3, 3, c, d)); 69e49d4f37SHong Zhang } 70e49d4f37SHong Zhang { 71e49d4f37SHong Zhang #define CUBEROOTOFTWO 1.2599210498948731647672106 72e49d4f37SHong Zhang 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}; 739566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC4, 4, 4, c, d)); 74e49d4f37SHong Zhang } 75e49d4f37SHong Zhang PetscFunctionReturn(0); 76e49d4f37SHong Zhang } 77e49d4f37SHong Zhang 78e49d4f37SHong Zhang /*@C 79e49d4f37SHong Zhang TSBasicSymplecticRegisterDestroy - Frees the list of schemes that were registered by TSBasicSymplecticRegister(). 80e49d4f37SHong Zhang 81e49d4f37SHong Zhang Not Collective 82e49d4f37SHong Zhang 83e49d4f37SHong Zhang Level: advanced 84e49d4f37SHong Zhang 85db781477SPatrick Sanan .seealso: `TSBasicSymplecticRegister()`, `TSBasicSymplecticRegisterAll()` 86e49d4f37SHong Zhang @*/ 87*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticRegisterDestroy(void) 88*d71ae5a4SJacob Faibussowitsch { 89e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 90e49d4f37SHong Zhang 91e49d4f37SHong Zhang PetscFunctionBegin; 92e49d4f37SHong Zhang while ((link = BasicSymplecticSchemeList)) { 93e49d4f37SHong Zhang BasicSymplecticScheme scheme = &link->sch; 94e49d4f37SHong Zhang BasicSymplecticSchemeList = link->next; 959566063dSJacob Faibussowitsch PetscCall(PetscFree2(scheme->c, scheme->d)); 969566063dSJacob Faibussowitsch PetscCall(PetscFree(scheme->name)); 979566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 98e49d4f37SHong Zhang } 99e49d4f37SHong Zhang TSBasicSymplecticRegisterAllCalled = PETSC_FALSE; 100e49d4f37SHong Zhang PetscFunctionReturn(0); 101e49d4f37SHong Zhang } 102e49d4f37SHong Zhang 103e49d4f37SHong Zhang /*@C 104e49d4f37SHong Zhang TSBasicSymplecticInitializePackage - This function initializes everything in the TSBasicSymplectic package. It is called 1058a690491SBarry Smith from TSInitializePackage(). 106e49d4f37SHong Zhang 107e49d4f37SHong Zhang Level: developer 108e49d4f37SHong Zhang 109db781477SPatrick Sanan .seealso: `PetscInitialize()` 110e49d4f37SHong Zhang @*/ 111*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticInitializePackage(void) 112*d71ae5a4SJacob Faibussowitsch { 113e49d4f37SHong Zhang PetscFunctionBegin; 114e49d4f37SHong Zhang if (TSBasicSymplecticPackageInitialized) PetscFunctionReturn(0); 115e49d4f37SHong Zhang TSBasicSymplecticPackageInitialized = PETSC_TRUE; 1169566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegisterAll()); 1179566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSBasicSymplecticFinalizePackage)); 118e49d4f37SHong Zhang PetscFunctionReturn(0); 119e49d4f37SHong Zhang } 120e49d4f37SHong Zhang 121e49d4f37SHong Zhang /*@C 122e49d4f37SHong Zhang TSBasicSymplecticFinalizePackage - This function destroys everything in the TSBasicSymplectic package. It is 123e49d4f37SHong Zhang called from PetscFinalize(). 124e49d4f37SHong Zhang 125e49d4f37SHong Zhang Level: developer 126e49d4f37SHong Zhang 127db781477SPatrick Sanan .seealso: `PetscFinalize()` 128e49d4f37SHong Zhang @*/ 129*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticFinalizePackage(void) 130*d71ae5a4SJacob Faibussowitsch { 131e49d4f37SHong Zhang PetscFunctionBegin; 132e49d4f37SHong Zhang TSBasicSymplecticPackageInitialized = PETSC_FALSE; 1339566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegisterDestroy()); 134e49d4f37SHong Zhang PetscFunctionReturn(0); 135e49d4f37SHong Zhang } 136e49d4f37SHong Zhang 137e49d4f37SHong Zhang /*@C 138e49d4f37SHong Zhang TSBasicSymplecticRegister - register a basic symplectic integration scheme by providing the coefficients. 139e49d4f37SHong Zhang 140e49d4f37SHong Zhang Not Collective, but the same schemes should be registered on all processes on which they will be used 141e49d4f37SHong Zhang 142e49d4f37SHong Zhang Input Parameters: 143e49d4f37SHong Zhang + name - identifier for method 144e49d4f37SHong Zhang . order - approximation order of method 145e49d4f37SHong Zhang . s - number of stages, this is the dimension of the matrices below 146e49d4f37SHong Zhang . c - coefficients for updating generalized position (dimension s) 147e49d4f37SHong Zhang - d - coefficients for updating generalized momentum (dimension s) 148e49d4f37SHong Zhang 149e49d4f37SHong Zhang Notes: 150e49d4f37SHong Zhang Several symplectic methods are provided, this function is only needed to create new methods. 151e49d4f37SHong Zhang 152e49d4f37SHong Zhang Level: advanced 153e49d4f37SHong Zhang 154db781477SPatrick Sanan .seealso: `TSBasicSymplectic` 155e49d4f37SHong Zhang @*/ 156*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticRegister(TSRosWType name, PetscInt order, PetscInt s, PetscReal c[], PetscReal d[]) 157*d71ae5a4SJacob Faibussowitsch { 158e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 159e49d4f37SHong Zhang BasicSymplecticScheme scheme; 160e49d4f37SHong Zhang 161e49d4f37SHong Zhang PetscFunctionBegin; 162e49d4f37SHong Zhang PetscValidCharPointer(name, 1); 163dadcf809SJacob Faibussowitsch PetscValidRealPointer(c, 4); 164dadcf809SJacob Faibussowitsch PetscValidRealPointer(d, 5); 165e49d4f37SHong Zhang 1669566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticInitializePackage()); 1679566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 168e49d4f37SHong Zhang scheme = &link->sch; 1699566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &scheme->name)); 170e49d4f37SHong Zhang scheme->order = order; 171e49d4f37SHong Zhang scheme->s = s; 1729566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &scheme->c, s, &scheme->d)); 1739566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(scheme->c, c, s)); 1749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(scheme->d, d, s)); 175e49d4f37SHong Zhang link->next = BasicSymplecticSchemeList; 176e49d4f37SHong Zhang BasicSymplecticSchemeList = link; 177e49d4f37SHong Zhang PetscFunctionReturn(0); 178e49d4f37SHong Zhang } 179e49d4f37SHong Zhang 180e49d4f37SHong Zhang /* 181e49d4f37SHong Zhang The simplified form of the equations are: 182e49d4f37SHong Zhang 183e49d4f37SHong Zhang $ p_{i+1} = p_i + c_i*g(q_i)*h 184e49d4f37SHong Zhang $ q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h 185e49d4f37SHong Zhang 186e49d4f37SHong Zhang Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p. 187e49d4f37SHong Zhang 188e49d4f37SHong Zhang To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps: 189e49d4f37SHong Zhang 190e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration multiplied by c_1 191e49d4f37SHong Zhang - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1 192e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2 193e49d4f37SHong Zhang - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2 194e49d4f37SHong Zhang 195e49d4f37SHong Zhang */ 196*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_BasicSymplectic(TS ts) 197*d71ae5a4SJacob Faibussowitsch { 198e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 199e49d4f37SHong Zhang BasicSymplecticScheme scheme = bsymp->scheme; 200e49d4f37SHong Zhang Vec solution = ts->vec_sol, update = bsymp->update, q, p, q_update, p_update; 201e49d4f37SHong Zhang IS is_q = bsymp->is_q, is_p = bsymp->is_p; 202e49d4f37SHong Zhang TS subts_q = bsymp->subts_q, subts_p = bsymp->subts_p; 203e49d4f37SHong Zhang PetscBool stageok; 204e49d4f37SHong Zhang PetscReal next_time_step = ts->time_step; 205e49d4f37SHong Zhang PetscInt iter; 206e49d4f37SHong Zhang 207e49d4f37SHong Zhang PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(solution, is_q, &q)); 2099566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(solution, is_p, &p)); 2109566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(update, is_q, &q_update)); 2119566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(update, is_p, &p_update)); 212e49d4f37SHong Zhang 213e49d4f37SHong Zhang for (iter = 0; iter < scheme->s; iter++) { 2149566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 215e49d4f37SHong Zhang /* update velocity p */ 216e49d4f37SHong Zhang if (scheme->c[iter]) { 2179566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(subts_p, ts->ptime, q, p_update)); 2189566063dSJacob Faibussowitsch PetscCall(VecAXPY(p, scheme->c[iter] * ts->time_step, p_update)); 219e49d4f37SHong Zhang } 220e49d4f37SHong Zhang /* update position q */ 221e49d4f37SHong Zhang if (scheme->d[iter]) { 2229566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(subts_q, ts->ptime, p, q_update)); 2239566063dSJacob Faibussowitsch PetscCall(VecAXPY(q, scheme->d[iter] * ts->time_step, q_update)); 224e49d4f37SHong Zhang ts->ptime = ts->ptime + scheme->d[iter] * ts->time_step; 225e49d4f37SHong Zhang } 2269566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &solution)); 2279566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok)); 2289371c9d4SSatish Balay if (!stageok) { 2299371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 2309371c9d4SSatish Balay PetscFunctionReturn(0); 2319371c9d4SSatish Balay } 2329566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok)); 2339371c9d4SSatish Balay if (!stageok) { 2349371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 2359371c9d4SSatish Balay PetscFunctionReturn(0); 2369371c9d4SSatish Balay } 237e49d4f37SHong Zhang } 238e49d4f37SHong Zhang 239e49d4f37SHong Zhang ts->time_step = next_time_step; 2409566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(solution, is_q, &q)); 2419566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(solution, is_p, &p)); 2429566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(update, is_q, &q_update)); 2439566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(update, is_p, &p_update)); 244e49d4f37SHong Zhang PetscFunctionReturn(0); 245e49d4f37SHong Zhang } 246e49d4f37SHong Zhang 247*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, void *ctx) 248*d71ae5a4SJacob Faibussowitsch { 249e49d4f37SHong Zhang PetscFunctionBegin; 250e49d4f37SHong Zhang PetscFunctionReturn(0); 251e49d4f37SHong Zhang } 252e49d4f37SHong Zhang 253*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 254*d71ae5a4SJacob Faibussowitsch { 255e49d4f37SHong Zhang PetscFunctionBegin; 256e49d4f37SHong Zhang PetscFunctionReturn(0); 257e49d4f37SHong Zhang } 258e49d4f37SHong Zhang 259*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, void *ctx) 260*d71ae5a4SJacob Faibussowitsch { 261e49d4f37SHong Zhang PetscFunctionBegin; 262e49d4f37SHong Zhang PetscFunctionReturn(0); 263e49d4f37SHong Zhang } 264e49d4f37SHong Zhang 265*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 266*d71ae5a4SJacob Faibussowitsch { 267e49d4f37SHong Zhang PetscFunctionBegin; 268e49d4f37SHong Zhang PetscFunctionReturn(0); 269e49d4f37SHong Zhang } 270e49d4f37SHong Zhang 271*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BasicSymplectic(TS ts) 272*d71ae5a4SJacob Faibussowitsch { 273e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 274e49d4f37SHong Zhang DM dm; 275e49d4f37SHong Zhang 276e49d4f37SHong Zhang PetscFunctionBegin; 2779566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q)); 2789566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p)); 2793c633725SBarry Smith 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"); 2809566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q)); 2819566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p)); 2823c633725SBarry Smith 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"); 283e49d4f37SHong Zhang 2849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update)); 285e49d4f37SHong Zhang 2869566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 2879566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */ 2889566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 289e49d4f37SHong Zhang if (dm) { 2909566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts)); 2919566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts)); 292e49d4f37SHong Zhang } 293e49d4f37SHong Zhang PetscFunctionReturn(0); 294e49d4f37SHong Zhang } 295e49d4f37SHong Zhang 296*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BasicSymplectic(TS ts) 297*d71ae5a4SJacob Faibussowitsch { 298e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 299e49d4f37SHong Zhang 300e49d4f37SHong Zhang PetscFunctionBegin; 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bsymp->update)); 302e49d4f37SHong Zhang PetscFunctionReturn(0); 303e49d4f37SHong Zhang } 304e49d4f37SHong Zhang 305*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BasicSymplectic(TS ts) 306*d71ae5a4SJacob Faibussowitsch { 307e49d4f37SHong Zhang PetscFunctionBegin; 3089566063dSJacob Faibussowitsch PetscCall(TSReset_BasicSymplectic(ts)); 3092e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL)); 3102e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL)); 3119566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 312e49d4f37SHong Zhang PetscFunctionReturn(0); 313e49d4f37SHong Zhang } 314e49d4f37SHong Zhang 315*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems *PetscOptionsObject) 316*d71ae5a4SJacob Faibussowitsch { 317e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 318e49d4f37SHong Zhang 319e49d4f37SHong Zhang PetscFunctionBegin; 320d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options"); 321e49d4f37SHong Zhang { 322e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 323e49d4f37SHong Zhang PetscInt count, choice; 324e49d4f37SHong Zhang PetscBool flg; 325e49d4f37SHong Zhang const char **namelist; 326e49d4f37SHong Zhang 3279371c9d4SSatish Balay for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) 3289371c9d4SSatish Balay ; 3299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 330e49d4f37SHong Zhang for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) namelist[count] = link->sch.name; 3319566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_basicsymplectic_type", "Family of basic symplectic integration method", "TSBasicSymplecticSetType", (const char *const *)namelist, count, bsymp->scheme->name, &choice, &flg)); 3329566063dSJacob Faibussowitsch if (flg) PetscCall(TSBasicSymplecticSetType(ts, namelist[choice])); 3339566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 334e49d4f37SHong Zhang } 335d0609cedSBarry Smith PetscOptionsHeadEnd(); 336e49d4f37SHong Zhang PetscFunctionReturn(0); 337e49d4f37SHong Zhang } 338e49d4f37SHong Zhang 339*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BasicSymplectic(TS ts, PetscViewer viewer) 340*d71ae5a4SJacob Faibussowitsch { 341e49d4f37SHong Zhang PetscFunctionBegin; 342e49d4f37SHong Zhang PetscFunctionReturn(0); 343e49d4f37SHong Zhang } 344e49d4f37SHong Zhang 345*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X) 346*d71ae5a4SJacob Faibussowitsch { 347e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 348e49d4f37SHong Zhang Vec update = bsymp->update; 349e49d4f37SHong Zhang PetscReal alpha = (ts->ptime - t) / ts->time_step; 350e49d4f37SHong Zhang 351e49d4f37SHong Zhang PetscFunctionBegin; 3529566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol)); 3539566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol)); 354e49d4f37SHong Zhang PetscFunctionReturn(0); 355e49d4f37SHong Zhang } 356e49d4f37SHong Zhang 357*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) 358*d71ae5a4SJacob Faibussowitsch { 359e49d4f37SHong Zhang PetscFunctionBegin; 360e49d4f37SHong Zhang *yr = 1.0 + xr; 361e49d4f37SHong Zhang *yi = xi; 362e49d4f37SHong Zhang PetscFunctionReturn(0); 363e49d4f37SHong Zhang } 364e49d4f37SHong Zhang 365e49d4f37SHong Zhang /*@C 366e49d4f37SHong Zhang TSBasicSymplecticSetType - Set the type of the basic symplectic method 367e49d4f37SHong Zhang 368e49d4f37SHong Zhang Logically Collective on TS 369e49d4f37SHong Zhang 370d8d19677SJose E. Roman Input Parameters: 371e49d4f37SHong Zhang + ts - timestepping context 372e49d4f37SHong Zhang - bsymptype - type of the symplectic scheme 373e49d4f37SHong Zhang 374e49d4f37SHong Zhang Options Database: 375147403d9SBarry Smith . -ts_basicsymplectic_type <scheme> - select the scheme 376e49d4f37SHong Zhang 377e49d4f37SHong Zhang Notes: 378147403d9SBarry Smith 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 379147403d9SBarry Smith that is intended to store the user-provided RHS function. 380e49d4f37SHong Zhang 381e49d4f37SHong Zhang Level: intermediate 382e49d4f37SHong Zhang @*/ 383*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype) 384*d71ae5a4SJacob Faibussowitsch { 385e49d4f37SHong Zhang PetscFunctionBegin; 386e49d4f37SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 387cac4c232SBarry Smith PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype)); 388e49d4f37SHong Zhang PetscFunctionReturn(0); 389e49d4f37SHong Zhang } 390e49d4f37SHong Zhang 391e49d4f37SHong Zhang /*@C 392e49d4f37SHong Zhang TSBasicSymplecticGetType - Get the type of the basic symplectic method 393e49d4f37SHong Zhang 394e49d4f37SHong Zhang Logically Collective on TS 395e49d4f37SHong Zhang 396d8d19677SJose E. Roman Input Parameters: 397e49d4f37SHong Zhang + ts - timestepping context 398e49d4f37SHong Zhang - bsymptype - type of the basic symplectic scheme 399e49d4f37SHong Zhang 400e49d4f37SHong Zhang Level: intermediate 401e49d4f37SHong Zhang @*/ 402*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype) 403*d71ae5a4SJacob Faibussowitsch { 404e49d4f37SHong Zhang PetscFunctionBegin; 405e49d4f37SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 406cac4c232SBarry Smith PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype)); 407e49d4f37SHong Zhang PetscFunctionReturn(0); 408e49d4f37SHong Zhang } 409e49d4f37SHong Zhang 410*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype) 411*d71ae5a4SJacob Faibussowitsch { 412e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 413e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 414e49d4f37SHong Zhang PetscBool match; 415e49d4f37SHong Zhang 416e49d4f37SHong Zhang PetscFunctionBegin; 417e49d4f37SHong Zhang if (bsymp->scheme) { 4189566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match)); 419e49d4f37SHong Zhang if (match) PetscFunctionReturn(0); 420e49d4f37SHong Zhang } 421e49d4f37SHong Zhang for (link = BasicSymplecticSchemeList; link; link = link->next) { 4229566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match)); 423e49d4f37SHong Zhang if (match) { 424e49d4f37SHong Zhang bsymp->scheme = &link->sch; 425e49d4f37SHong Zhang PetscFunctionReturn(0); 426e49d4f37SHong Zhang } 427e49d4f37SHong Zhang } 42898921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", bsymptype); 429e49d4f37SHong Zhang } 430e49d4f37SHong Zhang 431*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts, TSBasicSymplecticType *bsymptype) 432*d71ae5a4SJacob Faibussowitsch { 433e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 434e49d4f37SHong Zhang 435e49d4f37SHong Zhang PetscFunctionBegin; 436e49d4f37SHong Zhang *bsymptype = bsymp->scheme->name; 437e49d4f37SHong Zhang PetscFunctionReturn(0); 438e49d4f37SHong Zhang } 439e49d4f37SHong Zhang 440e49d4f37SHong Zhang /*MC 441e49d4f37SHong Zhang TSBasicSymplectic - ODE solver using basic symplectic integration schemes 442e49d4f37SHong Zhang 443e49d4f37SHong Zhang These methods are intened for separable Hamiltonian systems 444e49d4f37SHong Zhang 445e49d4f37SHong Zhang $ qdot = dH(q,p,t)/dp 446e49d4f37SHong Zhang $ pdot = -dH(q,p,t)/dq 447e49d4f37SHong Zhang 448e49d4f37SHong Zhang where the Hamiltonian can be split into the sum of kinetic energy and potential energy 449e49d4f37SHong Zhang 450e49d4f37SHong Zhang $ H(q,p,t) = T(p,t) + V(q,t). 451e49d4f37SHong Zhang 452e49d4f37SHong Zhang As a result, the system can be genearlly represented by 453e49d4f37SHong Zhang 454e49d4f37SHong Zhang $ qdot = f(p,t) = dT(p,t)/dp 455e49d4f37SHong Zhang $ pdot = g(q,t) = -dV(q,t)/dq 456e49d4f37SHong Zhang 457e49d4f37SHong Zhang and solved iteratively with 458e49d4f37SHong Zhang 459e49d4f37SHong Zhang $ q_new = q_old + d_i*h*f(p_old,t_old) 460e49d4f37SHong Zhang $ t_new = t_old + d_i*h 461e49d4f37SHong Zhang $ p_new = p_old + c_i*h*g(p_new,t_new) 462e49d4f37SHong Zhang $ i=0,1,...,n. 463e49d4f37SHong Zhang 464e49d4f37SHong Zhang 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. 465e49d4f37SHong Zhang 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. 466e49d4f37SHong Zhang 467e49d4f37SHong Zhang Reference: wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator) 468e49d4f37SHong Zhang 469e49d4f37SHong Zhang Level: beginner 470e49d4f37SHong Zhang 471db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()` 472e49d4f37SHong Zhang 473e49d4f37SHong Zhang M*/ 474*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts) 475*d71ae5a4SJacob Faibussowitsch { 476e49d4f37SHong Zhang TS_BasicSymplectic *bsymp; 477e49d4f37SHong Zhang 478e49d4f37SHong Zhang PetscFunctionBegin; 4799566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticInitializePackage()); 4804dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bsymp)); 481e49d4f37SHong Zhang ts->data = (void *)bsymp; 482e49d4f37SHong Zhang 483e49d4f37SHong Zhang ts->ops->setup = TSSetUp_BasicSymplectic; 484e49d4f37SHong Zhang ts->ops->step = TSStep_BasicSymplectic; 485e49d4f37SHong Zhang ts->ops->reset = TSReset_BasicSymplectic; 486e49d4f37SHong Zhang ts->ops->destroy = TSDestroy_BasicSymplectic; 487e49d4f37SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_BasicSymplectic; 488e49d4f37SHong Zhang ts->ops->view = TSView_BasicSymplectic; 489e49d4f37SHong Zhang ts->ops->interpolate = TSInterpolate_BasicSymplectic; 490e49d4f37SHong Zhang ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic; 491e49d4f37SHong Zhang 4929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", TSBasicSymplecticSetType_BasicSymplectic)); 4939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", TSBasicSymplecticGetType_BasicSymplectic)); 494e49d4f37SHong Zhang 4959566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticSetType(ts, TSBasicSymplecticDefault)); 496e49d4f37SHong Zhang PetscFunctionReturn(0); 497e49d4f37SHong Zhang } 498