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 @*/ 539371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticRegisterAll(void) { 54e49d4f37SHong Zhang PetscFunctionBegin; 55e49d4f37SHong Zhang if (TSBasicSymplecticRegisterAllCalled) PetscFunctionReturn(0); 56e49d4f37SHong Zhang TSBasicSymplecticRegisterAllCalled = PETSC_TRUE; 57e49d4f37SHong Zhang { 58e49d4f37SHong Zhang PetscReal c[1] = {1.0}, d[1] = {1.0}; 599566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICSIEULER, 1, 1, c, d)); 60e49d4f37SHong Zhang } 61e49d4f37SHong Zhang { 62e49d4f37SHong Zhang PetscReal c[2] = {0, 1.0}, d[2] = {0.5, 0.5}; 639566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTICVELVERLET, 2, 2, c, d)); 64e49d4f37SHong Zhang } 65e49d4f37SHong Zhang { 66e49d4f37SHong 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}; 679566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC3, 3, 3, c, d)); 68e49d4f37SHong Zhang } 69e49d4f37SHong Zhang { 70e49d4f37SHong Zhang #define CUBEROOTOFTWO 1.2599210498948731647672106 71e49d4f37SHong 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}; 729566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegister(TSBASICSYMPLECTIC4, 4, 4, c, d)); 73e49d4f37SHong Zhang } 74e49d4f37SHong Zhang PetscFunctionReturn(0); 75e49d4f37SHong Zhang } 76e49d4f37SHong Zhang 77e49d4f37SHong Zhang /*@C 78e49d4f37SHong Zhang TSBasicSymplecticRegisterDestroy - Frees the list of schemes that were registered by TSBasicSymplecticRegister(). 79e49d4f37SHong Zhang 80e49d4f37SHong Zhang Not Collective 81e49d4f37SHong Zhang 82e49d4f37SHong Zhang Level: advanced 83e49d4f37SHong Zhang 84db781477SPatrick Sanan .seealso: `TSBasicSymplecticRegister()`, `TSBasicSymplecticRegisterAll()` 85e49d4f37SHong Zhang @*/ 869371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticRegisterDestroy(void) { 87e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 88e49d4f37SHong Zhang 89e49d4f37SHong Zhang PetscFunctionBegin; 90e49d4f37SHong Zhang while ((link = BasicSymplecticSchemeList)) { 91e49d4f37SHong Zhang BasicSymplecticScheme scheme = &link->sch; 92e49d4f37SHong Zhang BasicSymplecticSchemeList = link->next; 939566063dSJacob Faibussowitsch PetscCall(PetscFree2(scheme->c, scheme->d)); 949566063dSJacob Faibussowitsch PetscCall(PetscFree(scheme->name)); 959566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 96e49d4f37SHong Zhang } 97e49d4f37SHong Zhang TSBasicSymplecticRegisterAllCalled = PETSC_FALSE; 98e49d4f37SHong Zhang PetscFunctionReturn(0); 99e49d4f37SHong Zhang } 100e49d4f37SHong Zhang 101e49d4f37SHong Zhang /*@C 102e49d4f37SHong Zhang TSBasicSymplecticInitializePackage - This function initializes everything in the TSBasicSymplectic package. It is called 1038a690491SBarry Smith from TSInitializePackage(). 104e49d4f37SHong Zhang 105e49d4f37SHong Zhang Level: developer 106e49d4f37SHong Zhang 107db781477SPatrick Sanan .seealso: `PetscInitialize()` 108e49d4f37SHong Zhang @*/ 1099371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticInitializePackage(void) { 110e49d4f37SHong Zhang PetscFunctionBegin; 111e49d4f37SHong Zhang if (TSBasicSymplecticPackageInitialized) PetscFunctionReturn(0); 112e49d4f37SHong Zhang TSBasicSymplecticPackageInitialized = PETSC_TRUE; 1139566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegisterAll()); 1149566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSBasicSymplecticFinalizePackage)); 115e49d4f37SHong Zhang PetscFunctionReturn(0); 116e49d4f37SHong Zhang } 117e49d4f37SHong Zhang 118e49d4f37SHong Zhang /*@C 119e49d4f37SHong Zhang TSBasicSymplecticFinalizePackage - This function destroys everything in the TSBasicSymplectic package. It is 120e49d4f37SHong Zhang called from PetscFinalize(). 121e49d4f37SHong Zhang 122e49d4f37SHong Zhang Level: developer 123e49d4f37SHong Zhang 124db781477SPatrick Sanan .seealso: `PetscFinalize()` 125e49d4f37SHong Zhang @*/ 1269371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticFinalizePackage(void) { 127e49d4f37SHong Zhang PetscFunctionBegin; 128e49d4f37SHong Zhang TSBasicSymplecticPackageInitialized = PETSC_FALSE; 1299566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticRegisterDestroy()); 130e49d4f37SHong Zhang PetscFunctionReturn(0); 131e49d4f37SHong Zhang } 132e49d4f37SHong Zhang 133e49d4f37SHong Zhang /*@C 134e49d4f37SHong Zhang TSBasicSymplecticRegister - register a basic symplectic integration scheme by providing the coefficients. 135e49d4f37SHong Zhang 136e49d4f37SHong Zhang Not Collective, but the same schemes should be registered on all processes on which they will be used 137e49d4f37SHong Zhang 138e49d4f37SHong Zhang Input Parameters: 139e49d4f37SHong Zhang + name - identifier for method 140e49d4f37SHong Zhang . order - approximation order of method 141e49d4f37SHong Zhang . s - number of stages, this is the dimension of the matrices below 142e49d4f37SHong Zhang . c - coefficients for updating generalized position (dimension s) 143e49d4f37SHong Zhang - d - coefficients for updating generalized momentum (dimension s) 144e49d4f37SHong Zhang 145e49d4f37SHong Zhang Notes: 146e49d4f37SHong Zhang Several symplectic methods are provided, this function is only needed to create new methods. 147e49d4f37SHong Zhang 148e49d4f37SHong Zhang Level: advanced 149e49d4f37SHong Zhang 150db781477SPatrick Sanan .seealso: `TSBasicSymplectic` 151e49d4f37SHong Zhang @*/ 1529371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticRegister(TSRosWType name, PetscInt order, PetscInt s, PetscReal c[], PetscReal d[]) { 153e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 154e49d4f37SHong Zhang BasicSymplecticScheme scheme; 155e49d4f37SHong Zhang 156e49d4f37SHong Zhang PetscFunctionBegin; 157e49d4f37SHong Zhang PetscValidCharPointer(name, 1); 158dadcf809SJacob Faibussowitsch PetscValidRealPointer(c, 4); 159dadcf809SJacob Faibussowitsch PetscValidRealPointer(d, 5); 160e49d4f37SHong Zhang 1619566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticInitializePackage()); 1629566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 163e49d4f37SHong Zhang scheme = &link->sch; 1649566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &scheme->name)); 165e49d4f37SHong Zhang scheme->order = order; 166e49d4f37SHong Zhang scheme->s = s; 1679566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &scheme->c, s, &scheme->d)); 1689566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(scheme->c, c, s)); 1699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(scheme->d, d, s)); 170e49d4f37SHong Zhang link->next = BasicSymplecticSchemeList; 171e49d4f37SHong Zhang BasicSymplecticSchemeList = link; 172e49d4f37SHong Zhang PetscFunctionReturn(0); 173e49d4f37SHong Zhang } 174e49d4f37SHong Zhang 175e49d4f37SHong Zhang /* 176e49d4f37SHong Zhang The simplified form of the equations are: 177e49d4f37SHong Zhang 178e49d4f37SHong Zhang $ p_{i+1} = p_i + c_i*g(q_i)*h 179e49d4f37SHong Zhang $ q_{i+1} = q_i + d_i*f(p_{i+1},t_{i+1})*h 180e49d4f37SHong Zhang 181e49d4f37SHong Zhang Several symplectic integrators are given below. An illustrative way to use them is to consider a particle with position q and velocity p. 182e49d4f37SHong Zhang 183e49d4f37SHong Zhang To apply a timestep with values c_{1,2},d_{1,2} to the particle, carry out the following steps: 184e49d4f37SHong Zhang 185e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration multiplied by c_1 186e49d4f37SHong Zhang - Update the position of the particle by adding to it its (updated) velocity multiplied by d_1 187e49d4f37SHong Zhang - Update the velocity of the particle by adding to it its acceleration (at the updated position) multiplied by c_2 188e49d4f37SHong Zhang - Update the position of the particle by adding to it its (double-updated) velocity multiplied by d_2 189e49d4f37SHong Zhang 190e49d4f37SHong Zhang */ 1919371c9d4SSatish Balay static PetscErrorCode TSStep_BasicSymplectic(TS ts) { 192e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 193e49d4f37SHong Zhang BasicSymplecticScheme scheme = bsymp->scheme; 194e49d4f37SHong Zhang Vec solution = ts->vec_sol, update = bsymp->update, q, p, q_update, p_update; 195e49d4f37SHong Zhang IS is_q = bsymp->is_q, is_p = bsymp->is_p; 196e49d4f37SHong Zhang TS subts_q = bsymp->subts_q, subts_p = bsymp->subts_p; 197e49d4f37SHong Zhang PetscBool stageok; 198e49d4f37SHong Zhang PetscReal next_time_step = ts->time_step; 199e49d4f37SHong Zhang PetscInt iter; 200e49d4f37SHong Zhang 201e49d4f37SHong Zhang PetscFunctionBegin; 2029566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(solution, is_q, &q)); 2039566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(solution, is_p, &p)); 2049566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(update, is_q, &q_update)); 2059566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(update, is_p, &p_update)); 206e49d4f37SHong Zhang 207e49d4f37SHong Zhang for (iter = 0; iter < scheme->s; iter++) { 2089566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ts->ptime)); 209e49d4f37SHong Zhang /* update velocity p */ 210e49d4f37SHong Zhang if (scheme->c[iter]) { 2119566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(subts_p, ts->ptime, q, p_update)); 2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(p, scheme->c[iter] * ts->time_step, p_update)); 213e49d4f37SHong Zhang } 214e49d4f37SHong Zhang /* update position q */ 215e49d4f37SHong Zhang if (scheme->d[iter]) { 2169566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(subts_q, ts->ptime, p, q_update)); 2179566063dSJacob Faibussowitsch PetscCall(VecAXPY(q, scheme->d[iter] * ts->time_step, q_update)); 218e49d4f37SHong Zhang ts->ptime = ts->ptime + scheme->d[iter] * ts->time_step; 219e49d4f37SHong Zhang } 2209566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ts->ptime, 0, &solution)); 2219566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime, solution, &stageok)); 2229371c9d4SSatish Balay if (!stageok) { 2239371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 2249371c9d4SSatish Balay PetscFunctionReturn(0); 2259371c9d4SSatish Balay } 2269566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts, ts->ptime + ts->time_step, update, &stageok)); 2279371c9d4SSatish Balay if (!stageok) { 2289371c9d4SSatish Balay ts->reason = TS_DIVERGED_STEP_REJECTED; 2299371c9d4SSatish Balay PetscFunctionReturn(0); 2309371c9d4SSatish Balay } 231e49d4f37SHong Zhang } 232e49d4f37SHong Zhang 233e49d4f37SHong Zhang ts->time_step = next_time_step; 2349566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(solution, is_q, &q)); 2359566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(solution, is_p, &p)); 2369566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(update, is_q, &q_update)); 2379566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(update, is_p, &p_update)); 238e49d4f37SHong Zhang PetscFunctionReturn(0); 239e49d4f37SHong Zhang } 240e49d4f37SHong Zhang 2419371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine, DM coarse, void *ctx) { 242e49d4f37SHong Zhang PetscFunctionBegin; 243e49d4f37SHong Zhang PetscFunctionReturn(0); 244e49d4f37SHong Zhang } 245e49d4f37SHong Zhang 2469371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 247e49d4f37SHong Zhang PetscFunctionBegin; 248e49d4f37SHong Zhang PetscFunctionReturn(0); 249e49d4f37SHong Zhang } 250e49d4f37SHong Zhang 2519371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm, DM subdm, void *ctx) { 252e49d4f37SHong Zhang PetscFunctionBegin; 253e49d4f37SHong Zhang PetscFunctionReturn(0); 254e49d4f37SHong Zhang } 255e49d4f37SHong Zhang 2569371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 257e49d4f37SHong Zhang PetscFunctionBegin; 258e49d4f37SHong Zhang PetscFunctionReturn(0); 259e49d4f37SHong Zhang } 260e49d4f37SHong Zhang 2619371c9d4SSatish Balay static PetscErrorCode TSSetUp_BasicSymplectic(TS ts) { 262e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 263e49d4f37SHong Zhang DM dm; 264e49d4f37SHong Zhang 265e49d4f37SHong Zhang PetscFunctionBegin; 2669566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "position", &bsymp->is_q)); 2679566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "momentum", &bsymp->is_p)); 2683c633725SBarry 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"); 2699566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "position", &bsymp->subts_q)); 2709566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "momentum", &bsymp->subts_p)); 2713c633725SBarry 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"); 272e49d4f37SHong Zhang 2739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &bsymp->update)); 274e49d4f37SHong Zhang 2759566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 2769566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */ 2779566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 278e49d4f37SHong Zhang if (dm) { 2799566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_BasicSymplectic, DMRestrictHook_BasicSymplectic, ts)); 2809566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_BasicSymplectic, DMSubDomainRestrictHook_BasicSymplectic, ts)); 281e49d4f37SHong Zhang } 282e49d4f37SHong Zhang PetscFunctionReturn(0); 283e49d4f37SHong Zhang } 284e49d4f37SHong Zhang 2859371c9d4SSatish Balay static PetscErrorCode TSReset_BasicSymplectic(TS ts) { 286e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 287e49d4f37SHong Zhang 288e49d4f37SHong Zhang PetscFunctionBegin; 2899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bsymp->update)); 290e49d4f37SHong Zhang PetscFunctionReturn(0); 291e49d4f37SHong Zhang } 292e49d4f37SHong Zhang 2939371c9d4SSatish Balay static PetscErrorCode TSDestroy_BasicSymplectic(TS ts) { 294e49d4f37SHong Zhang PetscFunctionBegin; 2959566063dSJacob Faibussowitsch PetscCall(TSReset_BasicSymplectic(ts)); 2962e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", NULL)); 2972e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", NULL)); 2989566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 299e49d4f37SHong Zhang PetscFunctionReturn(0); 300e49d4f37SHong Zhang } 301e49d4f37SHong Zhang 3029371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_BasicSymplectic(TS ts, PetscOptionItems *PetscOptionsObject) { 303e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 304e49d4f37SHong Zhang 305e49d4f37SHong Zhang PetscFunctionBegin; 306d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Basic symplectic integrator options"); 307e49d4f37SHong Zhang { 308e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 309e49d4f37SHong Zhang PetscInt count, choice; 310e49d4f37SHong Zhang PetscBool flg; 311e49d4f37SHong Zhang const char **namelist; 312e49d4f37SHong Zhang 3139371c9d4SSatish Balay for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) 3149371c9d4SSatish Balay ; 3159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 316e49d4f37SHong Zhang for (link = BasicSymplecticSchemeList, count = 0; link; link = link->next, count++) namelist[count] = link->sch.name; 3179566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_basicsymplectic_type", "Family of basic symplectic integration method", "TSBasicSymplecticSetType", (const char *const *)namelist, count, bsymp->scheme->name, &choice, &flg)); 3189566063dSJacob Faibussowitsch if (flg) PetscCall(TSBasicSymplecticSetType(ts, namelist[choice])); 3199566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 320e49d4f37SHong Zhang } 321d0609cedSBarry Smith PetscOptionsHeadEnd(); 322e49d4f37SHong Zhang PetscFunctionReturn(0); 323e49d4f37SHong Zhang } 324e49d4f37SHong Zhang 3259371c9d4SSatish Balay static PetscErrorCode TSView_BasicSymplectic(TS ts, PetscViewer viewer) { 326e49d4f37SHong Zhang PetscFunctionBegin; 327e49d4f37SHong Zhang PetscFunctionReturn(0); 328e49d4f37SHong Zhang } 329e49d4f37SHong Zhang 3309371c9d4SSatish Balay static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts, PetscReal t, Vec X) { 331e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 332e49d4f37SHong Zhang Vec update = bsymp->update; 333e49d4f37SHong Zhang PetscReal alpha = (ts->ptime - t) / ts->time_step; 334e49d4f37SHong Zhang 335e49d4f37SHong Zhang PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X, -ts->time_step, update, ts->vec_sol)); 3379566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X, 1.0 - alpha, alpha, ts->vec_sol)); 338e49d4f37SHong Zhang PetscFunctionReturn(0); 339e49d4f37SHong Zhang } 340e49d4f37SHong Zhang 3419371c9d4SSatish Balay static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) { 342e49d4f37SHong Zhang PetscFunctionBegin; 343e49d4f37SHong Zhang *yr = 1.0 + xr; 344e49d4f37SHong Zhang *yi = xi; 345e49d4f37SHong Zhang PetscFunctionReturn(0); 346e49d4f37SHong Zhang } 347e49d4f37SHong Zhang 348e49d4f37SHong Zhang /*@C 349e49d4f37SHong Zhang TSBasicSymplecticSetType - Set the type of the basic symplectic method 350e49d4f37SHong Zhang 351e49d4f37SHong Zhang Logically Collective on TS 352e49d4f37SHong Zhang 353d8d19677SJose E. Roman Input Parameters: 354e49d4f37SHong Zhang + ts - timestepping context 355e49d4f37SHong Zhang - bsymptype - type of the symplectic scheme 356e49d4f37SHong Zhang 357e49d4f37SHong Zhang Options Database: 358147403d9SBarry Smith . -ts_basicsymplectic_type <scheme> - select the scheme 359e49d4f37SHong Zhang 360e49d4f37SHong Zhang Notes: 361147403d9SBarry 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 362147403d9SBarry Smith that is intended to store the user-provided RHS function. 363e49d4f37SHong Zhang 364e49d4f37SHong Zhang Level: intermediate 365e49d4f37SHong Zhang @*/ 3669371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticSetType(TS ts, TSBasicSymplecticType bsymptype) { 367e49d4f37SHong Zhang PetscFunctionBegin; 368e49d4f37SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 369cac4c232SBarry Smith PetscTryMethod(ts, "TSBasicSymplecticSetType_C", (TS, TSBasicSymplecticType), (ts, bsymptype)); 370e49d4f37SHong Zhang PetscFunctionReturn(0); 371e49d4f37SHong Zhang } 372e49d4f37SHong Zhang 373e49d4f37SHong Zhang /*@C 374e49d4f37SHong Zhang TSBasicSymplecticGetType - Get the type of the basic symplectic method 375e49d4f37SHong Zhang 376e49d4f37SHong Zhang Logically Collective on TS 377e49d4f37SHong Zhang 378d8d19677SJose E. Roman Input Parameters: 379e49d4f37SHong Zhang + ts - timestepping context 380e49d4f37SHong Zhang - bsymptype - type of the basic symplectic scheme 381e49d4f37SHong Zhang 382e49d4f37SHong Zhang Level: intermediate 383e49d4f37SHong Zhang @*/ 3849371c9d4SSatish Balay PetscErrorCode TSBasicSymplecticGetType(TS ts, TSBasicSymplecticType *bsymptype) { 385e49d4f37SHong Zhang PetscFunctionBegin; 386e49d4f37SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 387cac4c232SBarry Smith PetscUseMethod(ts, "TSBasicSymplecticGetType_C", (TS, TSBasicSymplecticType *), (ts, bsymptype)); 388e49d4f37SHong Zhang PetscFunctionReturn(0); 389e49d4f37SHong Zhang } 390e49d4f37SHong Zhang 3919371c9d4SSatish Balay static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts, TSBasicSymplecticType bsymptype) { 392e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 393e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 394e49d4f37SHong Zhang PetscBool match; 395e49d4f37SHong Zhang 396e49d4f37SHong Zhang PetscFunctionBegin; 397e49d4f37SHong Zhang if (bsymp->scheme) { 3989566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(bsymp->scheme->name, bsymptype, &match)); 399e49d4f37SHong Zhang if (match) PetscFunctionReturn(0); 400e49d4f37SHong Zhang } 401e49d4f37SHong Zhang for (link = BasicSymplecticSchemeList; link; link = link->next) { 4029566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->sch.name, bsymptype, &match)); 403e49d4f37SHong Zhang if (match) { 404e49d4f37SHong Zhang bsymp->scheme = &link->sch; 405e49d4f37SHong Zhang PetscFunctionReturn(0); 406e49d4f37SHong Zhang } 407e49d4f37SHong Zhang } 40898921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", bsymptype); 409e49d4f37SHong Zhang } 410e49d4f37SHong Zhang 4119371c9d4SSatish Balay static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts, TSBasicSymplecticType *bsymptype) { 412e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic *)ts->data; 413e49d4f37SHong Zhang 414e49d4f37SHong Zhang PetscFunctionBegin; 415e49d4f37SHong Zhang *bsymptype = bsymp->scheme->name; 416e49d4f37SHong Zhang PetscFunctionReturn(0); 417e49d4f37SHong Zhang } 418e49d4f37SHong Zhang 419e49d4f37SHong Zhang /*MC 420e49d4f37SHong Zhang TSBasicSymplectic - ODE solver using basic symplectic integration schemes 421e49d4f37SHong Zhang 422e49d4f37SHong Zhang These methods are intened for separable Hamiltonian systems 423e49d4f37SHong Zhang 424e49d4f37SHong Zhang $ qdot = dH(q,p,t)/dp 425e49d4f37SHong Zhang $ pdot = -dH(q,p,t)/dq 426e49d4f37SHong Zhang 427e49d4f37SHong Zhang where the Hamiltonian can be split into the sum of kinetic energy and potential energy 428e49d4f37SHong Zhang 429e49d4f37SHong Zhang $ H(q,p,t) = T(p,t) + V(q,t). 430e49d4f37SHong Zhang 431e49d4f37SHong Zhang As a result, the system can be genearlly represented by 432e49d4f37SHong Zhang 433e49d4f37SHong Zhang $ qdot = f(p,t) = dT(p,t)/dp 434e49d4f37SHong Zhang $ pdot = g(q,t) = -dV(q,t)/dq 435e49d4f37SHong Zhang 436e49d4f37SHong Zhang and solved iteratively with 437e49d4f37SHong Zhang 438e49d4f37SHong Zhang $ q_new = q_old + d_i*h*f(p_old,t_old) 439e49d4f37SHong Zhang $ t_new = t_old + d_i*h 440e49d4f37SHong Zhang $ p_new = p_old + c_i*h*g(p_new,t_new) 441e49d4f37SHong Zhang $ i=0,1,...,n. 442e49d4f37SHong Zhang 443e49d4f37SHong 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. 444e49d4f37SHong 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. 445e49d4f37SHong Zhang 446e49d4f37SHong Zhang Reference: wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator) 447e49d4f37SHong Zhang 448e49d4f37SHong Zhang Level: beginner 449e49d4f37SHong Zhang 450db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()` 451e49d4f37SHong Zhang 452e49d4f37SHong Zhang M*/ 4539371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts) { 454e49d4f37SHong Zhang TS_BasicSymplectic *bsymp; 455e49d4f37SHong Zhang 456e49d4f37SHong Zhang PetscFunctionBegin; 4579566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticInitializePackage()); 458*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bsymp)); 459e49d4f37SHong Zhang ts->data = (void *)bsymp; 460e49d4f37SHong Zhang 461e49d4f37SHong Zhang ts->ops->setup = TSSetUp_BasicSymplectic; 462e49d4f37SHong Zhang ts->ops->step = TSStep_BasicSymplectic; 463e49d4f37SHong Zhang ts->ops->reset = TSReset_BasicSymplectic; 464e49d4f37SHong Zhang ts->ops->destroy = TSDestroy_BasicSymplectic; 465e49d4f37SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_BasicSymplectic; 466e49d4f37SHong Zhang ts->ops->view = TSView_BasicSymplectic; 467e49d4f37SHong Zhang ts->ops->interpolate = TSInterpolate_BasicSymplectic; 468e49d4f37SHong Zhang ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic; 469e49d4f37SHong Zhang 4709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticSetType_C", TSBasicSymplecticSetType_BasicSymplectic)); 4719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBasicSymplecticGetType_C", TSBasicSymplecticGetType_BasicSymplectic)); 472e49d4f37SHong Zhang 4739566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticSetType(ts, TSBasicSymplecticDefault)); 474e49d4f37SHong Zhang PetscFunctionReturn(0); 475e49d4f37SHong Zhang } 476