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 @*/ 53e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticRegisterAll(void) 54e49d4f37SHong Zhang { 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 @*/ 87e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticRegisterDestroy(void) 88e49d4f37SHong Zhang { 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 @*/ 111e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticInitializePackage(void) 112e49d4f37SHong Zhang { 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 @*/ 129e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticFinalizePackage(void) 130e49d4f37SHong Zhang { 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 @*/ 156e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticRegister(TSRosWType name,PetscInt order,PetscInt s,PetscReal c[],PetscReal d[]) 157e49d4f37SHong Zhang { 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 */ 196e49d4f37SHong Zhang static PetscErrorCode TSStep_BasicSymplectic(TS ts) 197e49d4f37SHong Zhang { 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)); 228e49d4f37SHong Zhang if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 2299566063dSJacob Faibussowitsch PetscCall(TSFunctionDomainError(ts,ts->ptime+ts->time_step,update,&stageok)); 230e49d4f37SHong Zhang if (!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 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 241e49d4f37SHong Zhang static PetscErrorCode DMCoarsenHook_BasicSymplectic(DM fine,DM coarse,void *ctx) 242e49d4f37SHong Zhang { 243e49d4f37SHong Zhang PetscFunctionBegin; 244e49d4f37SHong Zhang PetscFunctionReturn(0); 245e49d4f37SHong Zhang } 246e49d4f37SHong Zhang 247e49d4f37SHong Zhang static PetscErrorCode DMRestrictHook_BasicSymplectic(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 248e49d4f37SHong Zhang { 249e49d4f37SHong Zhang PetscFunctionBegin; 250e49d4f37SHong Zhang PetscFunctionReturn(0); 251e49d4f37SHong Zhang } 252e49d4f37SHong Zhang 253e49d4f37SHong Zhang static PetscErrorCode DMSubDomainHook_BasicSymplectic(DM dm,DM subdm,void *ctx) 254e49d4f37SHong Zhang { 255e49d4f37SHong Zhang PetscFunctionBegin; 256e49d4f37SHong Zhang PetscFunctionReturn(0); 257e49d4f37SHong Zhang } 258e49d4f37SHong Zhang 259e49d4f37SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_BasicSymplectic(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 260e49d4f37SHong Zhang { 261e49d4f37SHong Zhang PetscFunctionBegin; 262e49d4f37SHong Zhang PetscFunctionReturn(0); 263e49d4f37SHong Zhang } 264e49d4f37SHong Zhang 265e49d4f37SHong Zhang static PetscErrorCode TSSetUp_BasicSymplectic(TS ts) 266e49d4f37SHong Zhang { 267e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 268e49d4f37SHong Zhang DM dm; 269e49d4f37SHong Zhang 270e49d4f37SHong Zhang PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"position",&bsymp->is_q)); 2729566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"momentum",&bsymp->is_p)); 2733c633725SBarry 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"); 2749566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"position",&bsymp->subts_q)); 2759566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"momentum",&bsymp->subts_p)); 2763c633725SBarry 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"); 277e49d4f37SHong Zhang 2789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&bsymp->update)); 279e49d4f37SHong Zhang 2809566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&ts->adapt)); 2819566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); /* make sure to use fixed time stepping */ 2829566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 283e49d4f37SHong Zhang if (dm) { 2849566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_BasicSymplectic,DMRestrictHook_BasicSymplectic,ts)); 2859566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_BasicSymplectic,DMSubDomainRestrictHook_BasicSymplectic,ts)); 286e49d4f37SHong Zhang } 287e49d4f37SHong Zhang PetscFunctionReturn(0); 288e49d4f37SHong Zhang } 289e49d4f37SHong Zhang 290e49d4f37SHong Zhang static PetscErrorCode TSReset_BasicSymplectic(TS ts) 291e49d4f37SHong Zhang { 292e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 293e49d4f37SHong Zhang 294e49d4f37SHong Zhang PetscFunctionBegin; 2959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bsymp->update)); 296e49d4f37SHong Zhang PetscFunctionReturn(0); 297e49d4f37SHong Zhang } 298e49d4f37SHong Zhang 299e49d4f37SHong Zhang static PetscErrorCode TSDestroy_BasicSymplectic(TS ts) 300e49d4f37SHong Zhang { 301e49d4f37SHong Zhang PetscFunctionBegin; 3029566063dSJacob Faibussowitsch PetscCall(TSReset_BasicSymplectic(ts)); 303*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticSetType_C",NULL)); 304*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticGetType_C",NULL)); 3059566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 306e49d4f37SHong Zhang PetscFunctionReturn(0); 307e49d4f37SHong Zhang } 308e49d4f37SHong Zhang 309e49d4f37SHong Zhang static PetscErrorCode TSSetFromOptions_BasicSymplectic(PetscOptionItems *PetscOptionsObject,TS ts) 310e49d4f37SHong Zhang { 311e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 312e49d4f37SHong Zhang 313e49d4f37SHong Zhang PetscFunctionBegin; 314d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Basic symplectic integrator options"); 315e49d4f37SHong Zhang { 316e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 317e49d4f37SHong Zhang PetscInt count,choice; 318e49d4f37SHong Zhang PetscBool flg; 319e49d4f37SHong Zhang const char **namelist; 320e49d4f37SHong Zhang 321e49d4f37SHong Zhang for (link=BasicSymplecticSchemeList,count=0; link; link=link->next,count++) ; 3229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count,(char***)&namelist)); 323e49d4f37SHong Zhang for (link=BasicSymplecticSchemeList,count=0; link; link=link->next,count++) namelist[count] = link->sch.name; 3249566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_basicsymplectic_type","Family of basic symplectic integration method","TSBasicSymplecticSetType",(const char*const*)namelist,count,bsymp->scheme->name,&choice,&flg)); 3259566063dSJacob Faibussowitsch if (flg) PetscCall(TSBasicSymplecticSetType(ts,namelist[choice])); 3269566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 327e49d4f37SHong Zhang } 328d0609cedSBarry Smith PetscOptionsHeadEnd(); 329e49d4f37SHong Zhang PetscFunctionReturn(0); 330e49d4f37SHong Zhang } 331e49d4f37SHong Zhang 332e49d4f37SHong Zhang static PetscErrorCode TSView_BasicSymplectic(TS ts,PetscViewer viewer) 333e49d4f37SHong Zhang { 334e49d4f37SHong Zhang PetscFunctionBegin; 335e49d4f37SHong Zhang PetscFunctionReturn(0); 336e49d4f37SHong Zhang } 337e49d4f37SHong Zhang 338e49d4f37SHong Zhang static PetscErrorCode TSInterpolate_BasicSymplectic(TS ts,PetscReal t,Vec X) 339e49d4f37SHong Zhang { 340e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 341e49d4f37SHong Zhang Vec update = bsymp->update; 342e49d4f37SHong Zhang PetscReal alpha = (ts->ptime - t)/ts->time_step; 343e49d4f37SHong Zhang 344e49d4f37SHong Zhang PetscFunctionBegin; 3459566063dSJacob Faibussowitsch PetscCall(VecWAXPY(X,-ts->time_step,update,ts->vec_sol)); 3469566063dSJacob Faibussowitsch PetscCall(VecAXPBY(X,1.0-alpha,alpha,ts->vec_sol)); 347e49d4f37SHong Zhang PetscFunctionReturn(0); 348e49d4f37SHong Zhang } 349e49d4f37SHong Zhang 350e49d4f37SHong Zhang static PetscErrorCode TSComputeLinearStability_BasicSymplectic(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 351e49d4f37SHong Zhang { 352e49d4f37SHong Zhang PetscFunctionBegin; 353e49d4f37SHong Zhang *yr = 1.0 + xr; 354e49d4f37SHong Zhang *yi = xi; 355e49d4f37SHong Zhang PetscFunctionReturn(0); 356e49d4f37SHong Zhang } 357e49d4f37SHong Zhang 358e49d4f37SHong Zhang /*@C 359e49d4f37SHong Zhang TSBasicSymplecticSetType - Set the type of the basic symplectic method 360e49d4f37SHong Zhang 361e49d4f37SHong Zhang Logically Collective on TS 362e49d4f37SHong Zhang 363d8d19677SJose E. Roman Input Parameters: 364e49d4f37SHong Zhang + ts - timestepping context 365e49d4f37SHong Zhang - bsymptype - type of the symplectic scheme 366e49d4f37SHong Zhang 367e49d4f37SHong Zhang Options Database: 368147403d9SBarry Smith . -ts_basicsymplectic_type <scheme> - select the scheme 369e49d4f37SHong Zhang 370e49d4f37SHong Zhang Notes: 371147403d9SBarry 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 372147403d9SBarry Smith that is intended to store the user-provided RHS function. 373e49d4f37SHong Zhang 374e49d4f37SHong Zhang Level: intermediate 375e49d4f37SHong Zhang @*/ 376e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticSetType(TS ts,TSBasicSymplecticType bsymptype) 377e49d4f37SHong Zhang { 378e49d4f37SHong Zhang PetscFunctionBegin; 379e49d4f37SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 380cac4c232SBarry Smith PetscTryMethod(ts,"TSBasicSymplecticSetType_C",(TS,TSBasicSymplecticType),(ts,bsymptype)); 381e49d4f37SHong Zhang PetscFunctionReturn(0); 382e49d4f37SHong Zhang } 383e49d4f37SHong Zhang 384e49d4f37SHong Zhang /*@C 385e49d4f37SHong Zhang TSBasicSymplecticGetType - Get the type of the basic symplectic method 386e49d4f37SHong Zhang 387e49d4f37SHong Zhang Logically Collective on TS 388e49d4f37SHong Zhang 389d8d19677SJose E. Roman Input Parameters: 390e49d4f37SHong Zhang + ts - timestepping context 391e49d4f37SHong Zhang - bsymptype - type of the basic symplectic scheme 392e49d4f37SHong Zhang 393e49d4f37SHong Zhang Level: intermediate 394e49d4f37SHong Zhang @*/ 395e49d4f37SHong Zhang PetscErrorCode TSBasicSymplecticGetType(TS ts,TSBasicSymplecticType *bsymptype) 396e49d4f37SHong Zhang { 397e49d4f37SHong Zhang PetscFunctionBegin; 398e49d4f37SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 399cac4c232SBarry Smith PetscUseMethod(ts,"TSBasicSymplecticGetType_C",(TS,TSBasicSymplecticType*),(ts,bsymptype)); 400e49d4f37SHong Zhang PetscFunctionReturn(0); 401e49d4f37SHong Zhang } 402e49d4f37SHong Zhang 403e49d4f37SHong Zhang static PetscErrorCode TSBasicSymplecticSetType_BasicSymplectic(TS ts,TSBasicSymplecticType bsymptype) 404e49d4f37SHong Zhang { 405e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 406e49d4f37SHong Zhang BasicSymplecticSchemeLink link; 407e49d4f37SHong Zhang PetscBool match; 408e49d4f37SHong Zhang 409e49d4f37SHong Zhang PetscFunctionBegin; 410e49d4f37SHong Zhang if (bsymp->scheme) { 4119566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(bsymp->scheme->name,bsymptype,&match)); 412e49d4f37SHong Zhang if (match) PetscFunctionReturn(0); 413e49d4f37SHong Zhang } 414e49d4f37SHong Zhang for (link = BasicSymplecticSchemeList; link; link=link->next) { 4159566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->sch.name,bsymptype,&match)); 416e49d4f37SHong Zhang if (match) { 417e49d4f37SHong Zhang bsymp->scheme = &link->sch; 418e49d4f37SHong Zhang PetscFunctionReturn(0); 419e49d4f37SHong Zhang } 420e49d4f37SHong Zhang } 42198921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",bsymptype); 422e49d4f37SHong Zhang } 423e49d4f37SHong Zhang 424e49d4f37SHong Zhang static PetscErrorCode TSBasicSymplecticGetType_BasicSymplectic(TS ts,TSBasicSymplecticType *bsymptype) 425e49d4f37SHong Zhang { 426e49d4f37SHong Zhang TS_BasicSymplectic *bsymp = (TS_BasicSymplectic*)ts->data; 427e49d4f37SHong Zhang 428e49d4f37SHong Zhang PetscFunctionBegin; 429e49d4f37SHong Zhang *bsymptype = bsymp->scheme->name; 430e49d4f37SHong Zhang PetscFunctionReturn(0); 431e49d4f37SHong Zhang } 432e49d4f37SHong Zhang 433e49d4f37SHong Zhang /*MC 434e49d4f37SHong Zhang TSBasicSymplectic - ODE solver using basic symplectic integration schemes 435e49d4f37SHong Zhang 436e49d4f37SHong Zhang These methods are intened for separable Hamiltonian systems 437e49d4f37SHong Zhang 438e49d4f37SHong Zhang $ qdot = dH(q,p,t)/dp 439e49d4f37SHong Zhang $ pdot = -dH(q,p,t)/dq 440e49d4f37SHong Zhang 441e49d4f37SHong Zhang where the Hamiltonian can be split into the sum of kinetic energy and potential energy 442e49d4f37SHong Zhang 443e49d4f37SHong Zhang $ H(q,p,t) = T(p,t) + V(q,t). 444e49d4f37SHong Zhang 445e49d4f37SHong Zhang As a result, the system can be genearlly represented by 446e49d4f37SHong Zhang 447e49d4f37SHong Zhang $ qdot = f(p,t) = dT(p,t)/dp 448e49d4f37SHong Zhang $ pdot = g(q,t) = -dV(q,t)/dq 449e49d4f37SHong Zhang 450e49d4f37SHong Zhang and solved iteratively with 451e49d4f37SHong Zhang 452e49d4f37SHong Zhang $ q_new = q_old + d_i*h*f(p_old,t_old) 453e49d4f37SHong Zhang $ t_new = t_old + d_i*h 454e49d4f37SHong Zhang $ p_new = p_old + c_i*h*g(p_new,t_new) 455e49d4f37SHong Zhang $ i=0,1,...,n. 456e49d4f37SHong Zhang 457e49d4f37SHong 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. 458e49d4f37SHong 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. 459e49d4f37SHong Zhang 460e49d4f37SHong Zhang Reference: wikipedia (https://en.wikipedia.org/wiki/Symplectic_integrator) 461e49d4f37SHong Zhang 462e49d4f37SHong Zhang Level: beginner 463e49d4f37SHong Zhang 464db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetType()`, `TSRHSSplitSetIS()`, `TSRHSSplitSetRHSFunction()` 465e49d4f37SHong Zhang 466e49d4f37SHong Zhang M*/ 467e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_BasicSymplectic(TS ts) 468e49d4f37SHong Zhang { 469e49d4f37SHong Zhang TS_BasicSymplectic *bsymp; 470e49d4f37SHong Zhang 471e49d4f37SHong Zhang PetscFunctionBegin; 4729566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticInitializePackage()); 4739566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&bsymp)); 474e49d4f37SHong Zhang ts->data = (void*)bsymp; 475e49d4f37SHong Zhang 476e49d4f37SHong Zhang ts->ops->setup = TSSetUp_BasicSymplectic; 477e49d4f37SHong Zhang ts->ops->step = TSStep_BasicSymplectic; 478e49d4f37SHong Zhang ts->ops->reset = TSReset_BasicSymplectic; 479e49d4f37SHong Zhang ts->ops->destroy = TSDestroy_BasicSymplectic; 480e49d4f37SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_BasicSymplectic; 481e49d4f37SHong Zhang ts->ops->view = TSView_BasicSymplectic; 482e49d4f37SHong Zhang ts->ops->interpolate = TSInterpolate_BasicSymplectic; 483e49d4f37SHong Zhang ts->ops->linearstability = TSComputeLinearStability_BasicSymplectic; 484e49d4f37SHong Zhang 4859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticSetType_C",TSBasicSymplecticSetType_BasicSymplectic)); 4869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBasicSymplecticGetType_C",TSBasicSymplecticGetType_BasicSymplectic)); 487e49d4f37SHong Zhang 4889566063dSJacob Faibussowitsch PetscCall(TSBasicSymplecticSetType(ts,TSBasicSymplecticDefault)); 489e49d4f37SHong Zhang PetscFunctionReturn(0); 490e49d4f37SHong Zhang } 491