Step 4: Higher-order algorithms
In the last step, we implemented a simple Euler method to integrate the solar system for 200 years. The poor accuracy is expected, because we used a simple first order algorithm with global error \(\mathcal{O}(\Delta t)\). In this step, we will implement 3 new algorithms to improve our simulation.
Relative energy error
Before implementing new algorithms, we need a way to measure the accuracy of the simulation. By conservation of energy, we can assume that better algorithms will have a better conservation on the total energy of the system. Computing the total energy is trivial. All we need is the solution array.
Then, we can compute the relative energy error as
This gives the following code:
def compute_rel_energy_error(
sol_x: np.ndarray, sol_v: np.ndarray, system: System
) -> np.ndarray:
"""
Compute the relative energy error of the simulation.
Parameters
----------
sol_x : np.ndarray
Solution position array with shape (N_steps, num_particles, 3).
sol_v : np.ndarray
Solution velocity array with shape (N_steps, num_particles, 3).
system : System
System object.
Returns
-------
energy_error : np.ndarray
Relative energy error of the simulation, with shape (N_steps,).
"""
# Allocate memory and initialize arrays
n_steps = sol_x.shape[0]
num_particles = system.num_particles
m = system.m
G = system.G
rel_energy_error = np.zeros(n_steps)
# Compute the total energy (KE + PE)
for count in range(n_steps):
x = sol_x[count]
v = sol_v[count]
for i in range(num_particles):
# KE
rel_energy_error[count] += 0.5 * m[i] * np.linalg.norm(v[i]) ** 2
# PE
for j in range(i + 1, num_particles):
rel_energy_error[count] -= G * m[i] * m[j] / np.linalg.norm(x[i] - x[j])
# Compute the relative energy error
initial_energy = rel_energy_error[0]
rel_energy_error = (rel_energy_error - initial_energy) / initial_energy
rel_energy_error = np.abs(rel_energy_error)
return rel_energy_error
Tip
If you found it too slow, you are encouraged to vectorize it, just like what we did in step 2
Now, we can plot the relative energy error with the following code, where y-axis is in log scale.
def plot_rel_energy_error(rel_energy_error: np.ndarray, sol_t: np.ndarray) -> None:
"""
Plot the relative energy error.
Parameters
----------
rel_energy_error : np.ndarray
Relative energy error of the simulation, with shape (N_steps,).
sol_t : np.ndarray
Solution time array with shape (N_steps,).
"""
plt.figure()
plt.plot(sol_t, rel_energy_error)
plt.yscale("log")
plt.xlabel("Time step")
plt.ylabel("Relative Energy Error")
plt.title("Relative Energy Error vs Time Step")
plt.show()
We copy step3.py
to step4.py
and add the following code to the
end of the main
function:
def main() -> None:
...
# Compute and plot relative energy error
rel_energy_error = compute_rel_energy_error(sol_x, sol_v, system)
print(f"Relative energy error: {rel_energy_error[-1]:.3g}")
plot_rel_energy_error(rel_energy_error, sol_t / 365.24)
Euler method
Let us run the simulation again. We have the following plots:
The final relative energy error \(\sim 10^{-1}\), which is not very good.
Euler-Cromer method
Now, we begin implementing our first new algorithm. The Euler-Cromer method, also known as semi-implicit Euler method, is a simple modification of the Euler method,
Notice that the position update is done using the updated velocity \(\mathbf{v}_{n+1}\) instead of \(\mathbf{v}_n\). Therefore, it is an first order semi-implicit method. Although the global error is still \(\mathcal{O}(\Delta t)\), it is a symplectic method, which implies that the energy error over time is bounded. This is a very useful property for long time-scale simulations. We have the following code:
def euler_cromer(a: np.ndarray, system: System, dt: float) -> None:
"""
Advance one step with the Euler-Cromer method.
Parameters
----------
a : np.ndarray
Gravitational accelerations array with shape (N, 3).
system : System
System object.
dt : float
Time step.
"""
acceleration(a, system)
system.v += a * dt
system.x += system.v * dt
Running the simulation again, we have the following plots:
The trajectory looks a lot better than the Euler method! Also, the relative energy error is bounded at \(\sim 10^{-4}\). This provides long term stability for the simulation.
Energy conservation \(\iff\) Higher accuracy?
Although it seems like energy conservation implies higher accuracy, this is not neccessarily true. Recall that the global error for Euler-Cromer method is still \(\mathcal{O}(\Delta t)\), even though its energy error is bounded over time. Therefore, we can only state the opposite direction: Higher accuracy \(\implies\) Energy conservation.
Runge-Kutta method
The Runge-Kutta method is a family of ODE solvers. First, let us look at a second order variation called the midpoint method.
-
In the Euler method, we are using the slope evaluated at the beginning of the interval to update the position and velocity.
-
In the Euler-Cromer method, for the position update, we are using the slope evaluated at the end instead.
What about using the slope evaluated at the center of the interval? We can approximate the position and velocity at the midpoint using one Euler step, then perform the updates using the midpoint values. This gives us the following updates:
Using a more general notation, we can write the midpoint method as
A simple taylor expansion will shows that the local truncation error is indeed \(\mathcal{O}(\Delta t^3)\), which means that it is a second order method. A more commonly used Runge-Kutta method is the 4th order Runge-Kutta method (RK4), which provides a good balance between accuracy and computational cost. It is given by
In our code, we can write out the computation of each term explicitly. However, I prefer a more general
approach by defining the coeff
and weights
arrays instead. The coeff
array is given as
for the computation of \(k_2\), \(k_3\) and \(k_4\). The weights
array is given as
for the final update. The code is given as follows:
def rk4(a: np.ndarray, system: System, dt: float) -> None:
"""
Advance one step with the RK4 method.
Parameters
----------
a : np.ndarray
Gravitational accelerations array with shape (N, 3).
system : System
System object.
dt : float
Time step.
"""
num_stages = 4
coeff = np.array([0.5, 0.5, 1.0])
weights = np.array([1.0, 2.0, 2.0, 1.0]) / 6.0
# Allocate memory and initialize arrays
x0 = system.x.copy()
v0 = system.v.copy()
xk = np.zeros((num_stages, system.num_particles, 3))
vk = np.zeros((num_stages, system.num_particles, 3))
# Initial stage
acceleration(a, system)
xk[0] = v0
vk[0] = a
# Compute the stages
for stage in range(1, num_stages):
# Compute acceleration
system.x = x0 + dt * coeff[stage - 1] * xk[stage - 1]
acceleration(a, system)
# Compute xk and vk
xk[stage] = v0 + dt * coeff[stage - 1] * vk[stage - 1]
vk[stage] = a
# Advance step
dx = 0.0
dv = 0.0
for stage in range(num_stages):
dx += weights[stage] * xk[stage]
dv += weights[stage] * vk[stage]
system.x = x0 + dt * dx
system.v = v0 + dt * dv
Tip
The final loop can be vectorized to improve performance:
Let's run the simulation again. We only show the relative energy error plot:
The final error is in the order of \(10^{-6}\), which is quite nice. However, the error is growing over time, so this may not be a good choice for long term simulations.
Leapfrog method
Our final algorithm is the leapfrog method, which is a second-order symplectic method. Similar to the Euler-Cromer method, it conserves energy over time. We will implement the Kick-Drift-Kick (KDK) variant, which is given by a velocity kick for half time step,
a position drift for a full time step,
and a final velocity kick for half time step,
Tip
For optimization, you can combine the final velocity kick from the last time step with the first velocity kick. However, you need to be careful because now the velocity and position is not synchronized. There is also a synchronized version of the leapfrog method called the velocity Verlet method.
The implementation is simple:
def leapfrog(a: np.ndarray, system: System, dt: float) -> None:
"""
Advance one step with the LeapFrog method.
Parameters
----------
a : np.ndarray
Gravitational accelerations array with shape (N, 3).
system : System
System object.
dt : float
Time step.
"""
# Velocity kick (v_1/2)
acceleration(a, system)
system.v += a * 0.5 * dt
# Position drift (x_1)
system.x += system.v * dt
# Velocity kick (v_1)
acceleration(a, system)
system.v += a * 0.5 * dt
Running the simulation again, we have the relative energy error plot:
The relative energy error is bounded at \(\sim 10^{-6}\), which is better than the Euler-Cromer method!
Which algorithm should I choose?
To choose between RK4 and LeapFrog, below are some factors that you may consider:
- Time scale: If you are simulating a long time scale, LeapFrog should be better as it conserves energy.
- Accuracy: If accuracy matters to you, then RK4 is a better choice since it is higher order.
- Computational cost: RK4 is very expensive as it requires 4 acceleration evaluations per step. In contrast, LeapFrog only requires 1 acceleration evaluation although it is second order.
-
Experiment: If you are not sure which one to choose, you can always experiment and compare the results. Below is a relative energy error plot I made for the solar system simulation using both RK4 and LeapFrog. I chose two very small dt while ensuring that the runtime is the same for both methods. In terms of the results, seems like RK4 is better in this case, but note that I used some techniques to remove the rounding error. (See Reducing round off error)
Summary
In this step, we have implemented 3 new algorithms: Euler-Cromer, RK4 and Leapfrog. RK4 and Leapfrog are both very popular algorithms for N-body simulations. All algorithms we implemented so far are fixed time step methods, which may not be very flexible for chaotic systems with close encounters. It is also a headache to tune the time step. In the next step, we will implement an adaptive time-stepping method which uses a tolerance parameter to control the time step instead.
Full scripts
The full scripts are available at 5_steps_to_n_body_simulation/python/
,
or https://github.com/alvinng4/grav_sim/blob/main/5_steps_to_n_body_simulation/python/
step4.py (Click to expand)
common.py (Click to expand)
5_steps_to_n_body_simulation/python/common.py | |
---|---|
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 |
|