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