Skip to content

Commit

Permalink
Missing Inertial<->DHC transform in TRACE (#772)
Browse files Browse the repository at this point in the history
* Missing Inertial<->DHC transform

* added pericenter test

---------

Co-authored-by: Tiger Lu <tigerlu@Tigers-MBP-2.wireless.yale.internal>
  • Loading branch information
tigerchenlu98 and Tiger Lu committed May 7, 2024
1 parent 7ebc08c commit 5b5d4db
Show file tree
Hide file tree
Showing 3 changed files with 9 additions and 7 deletions.
4 changes: 2 additions & 2 deletions docs/index.md
Original file line number Diff line number Diff line change
Expand Up @@ -19,8 +19,8 @@ Abstract submission and registration are now open at [https://hannorein.github.i
* No dependencies on external libraries.
* Runs natively on Linux, MacOS, and Windows.
* Symplectic integrators ([WHFast](integrators/#whfast), [SEI](integrators/#sei), [LEAPFROG](integrators/#leapfrog), [EOS](integrators/#embedded-operator-splitting-method-eos))
* Hybrid symplectic integrators for planetary dynamics with close encounters ([MERCURIUS](integrators/#mercurius)
* Hybrid reversible integrators for planetary dynamics with arbitrary close encounters [TRACE](integrators/#trace))
* Hybrid symplectic integrators for planetary dynamics with close encounters ([MERCURIUS](integrators/#mercurius))
* Hybrid reversible integrators for planetary dynamics with arbitrary close encounters ([TRACE](integrators/#trace))
* High order symplectic integrators for integrating planetary systems ([SABA](integrators/#saba), WH Kernel methods)
* High accuracy non-symplectic integrator with adaptive time-stepping ([IAS15](integrators/#ias15))
* Can integrate arbitrary user-defined ODEs that are coupled to N-body dynamics for tides, spin, etc
Expand Down
10 changes: 5 additions & 5 deletions rebound/tests/test_trace.py
Original file line number Diff line number Diff line change
Expand Up @@ -442,27 +442,27 @@ def jacobi(sim):

self.assertLess(dE_trace, 1e-6) # reasonable precision for trace
self.assertLess(time_trace,2.0*time_ias15) # not much slower than ias15
'''

def test_pericenter(self):

sim = pericenter_sim()
sim.integrator = "ias15"
start_ias15=datetime.now()
sim.integrate(3000.)
sim.integrate(10 * 2 * math.pi * 29.4)
time_ias15 = (datetime.now()-start_ias15).total_seconds()

sim = pericenter_sim()
sim.integrator = "trace"
sim.dt = 0.15 * 2 * math.pi
E0 = sim.energy()
start_trace=datetime.now()
sim.integrate(3000.)
sim.integrate(10 * 2 * math.pi * 29.4)
time_trace = (datetime.now()-start_trace).total_seconds()
dE_trace = abs((sim.energy() - E0)/E0)

self.assertLess(dE_trace,1e-3) # reasonable precision for trace
self.assertLess(dE_trace,1e-4) # reasonable precision for trace
self.assertLess(time_trace,time_ias15) # faster than ias15
'''
def test_trace_simulationarchive(self):
sim = chaotic_exchange_sim()
sim.integrator = "trace"
Expand Down
2 changes: 2 additions & 0 deletions src/integrator_trace.c
Original file line number Diff line number Diff line change
Expand Up @@ -717,6 +717,7 @@ static void reb_integrator_trace_step(struct reb_simulation* const r){
const double old_t = r->t;
r->gravity = REB_GRAVITY_BASIC;
r->ri_trace.mode = REB_TRACE_MODE_FULL; // for collision search
reb_integrator_trace_dh_to_inertial(r);
switch (r->ri_trace.peri_mode){
case REB_TRACE_PERI_FULL_IAS15:
// Run default IAS15 integration
Expand Down Expand Up @@ -779,6 +780,7 @@ static void reb_integrator_trace_step(struct reb_simulation* const r){
r->gravity = REB_GRAVITY_TRACE;
r->t = old_t; // final time will be set later
r->dt = old_dt;
reb_integrator_trace_inertial_to_dh(r);
}
}

Expand Down

0 comments on commit 5b5d4db

Please sign in to comment.