Python · C · C++ · Fortran 四个引擎的蛙跳 (Velocity-Verlet) 算法与边界条件实现对比
最后更新:2026-05-20(力开关系统、万有引力、逐自由度约束)
本框架实现了四套语言的计算引擎。各引擎共享相同的物理模型:
gravity_fielddamping_forceelastic_forcegravity_interactionfix_x fix_y fix_z 列控制各自由度是否参与迭代Python 引擎支持 4 种积分算法,C/C++/Fortran 引擎于 2026-05-20 重写,现支持全部 4 种算法,且边界条件与 Python 完全一致。经验证,四种引擎对于同一输入(case01, NT=1000, leapfrog, B=0)输出完全一致(至双精度浮点精度)。
| 引擎 | 源文件 | 语言标准 | 算法支持 |
|---|---|---|---|
| Python | compute.py | Python 3.8+ | 显式欧拉、隐式欧拉、中点法、蛙跳法 |
| C | engines/c/main.c | C11 | 同上(2026-05-20 新增 3 种) |
| C++ | engines/cpp/main.cpp | C++17 | 同上(2026-05-20 新增 3 种) |
| Fortran | engines/fortran/main.f90 | Fortran 90 | 同上(2026-05-20 新增 3 种) |
compute.py
Python 的 Leapfrog_Method 实现了标准 velocity-Verlet 格式。三步法对应于:
Python917def Leapfrog_Method(x, y, z, vx, vy, vz, dt, m, g, b):
918 ax, ay, az = compute_acceleration(x, y, z, vx, vy, vz, m, g, b)
919 vx_half = vx + 0.5 * ax * dt
920 vy_half = vy + 0.5 * ay * dt
921 vz_half = vz + 0.5 * az * dt
922
923 x = x + vx_half * dt
924 y = y + vy_half * dt
925 z = z + vz_half * dt
926
927 gamma_x = b[0] / m
928 gamma_y = b[1] / m
929 gamma_z = b[2] / m
930 vx_next = (vx_half + 0.5 * g[0] * dt) / (1.0 + 0.5 * gamma_x * dt)
931 vy_next = (vy_half + 0.5 * g[1] * dt) / (1.0 + 0.5 * gamma_y * dt)
932 vz_next = (vz_half + 0.5 * g[2] * dt) / (1.0 + 0.5 * gamma_z * dt)
933 ax_next, ay_next, az_next = compute_acceleration(
934 x, y, z, vx_next, vy_next, vz_next, m, g, b)
935 vx = vx_half + 0.5 * ax_next * dt
936 vy = vy_half + 0.5 * ay_next * dt
937 vz = vz_half + 0.5 * az_next * dt
938 return x, y, z, vx, vy, vz
Python 的蛙跳法中,阻尼项用隐式格式处理。gamma = B / m,然后 v_next = (v_half + ½·g·Δt) / (1 + ½·γ·Δt)。这个 v_next 是阻尼修正后的速度估计值,被传给第二次 compute_acceleration 作为入参。当 B = 0 时退化为显式 v_next = v_half + ½·g·Δt。
注意:vx_next 是独立变量,不覆盖 vx_half。最终速度 vx = vx_half + ½·a_next·dt 使用的是原始的半步速度。这是正确的做法。
Python940def Limit_in_box(a, amin, amax, va):
941 """限制物体在边界内,发生碰撞时反弹。"""
942 over = a > amax
943 under = a < amin
944 a = np.where(over, amax, a)
945 a = np.where(under, amin, a)
946 va = np.where(over | under, -va, va)
947 return a, va
Python949def apply_motion_update(x, y, z, vx, vy, vz, dt, m, g, b):
957 elif METHOD == "leapfrog":
958 x, y, z, vx, vy, vz = Leapfrog_Method(x, y, z, vx, vy, vz, dt, m, g, b)
962 x, vx = Limit_in_box(x, X_MIN, X_MAX, vx)
963 y, vy = Limit_in_box(y, Y_MIN, Y_MAX, vy)
964 z, vz = Limit_in_box(z, Z_MIN, Z_MAX, vz)
965 return x, y, z, vx, vy, vz
Limit_in_boxLimit_in_box:clamp 位置到边界 + 越界方向速度取反(弹性碰撞模型)2026-05-20 添加了独立的力开关系统,每种力可通过 input.txt 中的 0/1 开关独立控制。各引擎的受力计算均基于这些开关条件累加。
| 开关名 | 默认值 | 控制的力 | 物理公式 |
|---|---|---|---|
gravity_field | 1 | 均匀重力场(常数矢量 g) | F = m · g |
gravity_interaction | 0 | 原子间牛顿万有引力 | F = G · mi · mj / r² |
elastic_force | 1 | 弹簧键胡克力 | F = k(r − r0) |
damping_force | 0 | 线性阻尼 | F = −B · v |
Python compute_force() / C/C++ compute_acceleration() / Fortran accel() 统一采用以下模式:
伪代码// 1. 清零
ax = 0; ay = 0; az = 0
// 2. 按开关累加
if (gravity_field): ax += G[0]; ay += G[1]; az += G[2]
if (damping_force): ax -= B[0] * vx / m; ...
if (elastic_force): ax += bond_loop(...)
if (gravity_interaction): ax += pair_gravity_loop(...)
每个力都包裹在 if (开关) 条件中。开关关闭时,对应力的计算代码完全跳过,不影响运行效率。
绘图时 dynamics.py 的 run_case() 根据力开关状态动态计算对应的势能分量:
| 力开关 | 势能公式 | 关=0 时 |
|---|---|---|
gravity_field=1 | Ug = −m · G · r | Ug 归零 |
gravity_interaction=1 | Ugg = −∑i<j G · mi · mj / rij | Ugg 归零 |
elastic_force=1 | Us = ½ ∑ k(r − r0)² | Us 归零 |
动能=蓝色、均匀重力势能=绿色、弹性势能=橙色、万有引力势能=紫色、总能量=红色虚线。
engines/c/main.c
注意:2026-05-20 重写。原实现为 Symplectic Euler(非 Leapfrog),现已修正为 Velocity-Verlet。
C420/* ── 蛙跳法(Velocity-Verlet)── */
421static void leapfrog_step(
422 int n, double *x, double *y, double *z,
423 double *vx, double *vy, double *vz,
424 const double *m, const double G[3], const double B[3],
425 const BondData *bonds, const int *fixed, double dt)
426{
427 double *ax = (double*)alloca(n * sizeof(double));
428 double *ay = (double*)alloca(n * sizeof(double));
429 double *az = (double*)alloca(n * sizeof(double));
430 compute_acceleration(n, x, y, z, vx, vy, vz, m, G, B, bonds, ax, ay, az);
431
432 /* 半推速度:v_half = v + 0.5*a*dt */
433 for (int i = 0; i < n; i++) {
434 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
435 vx[i] += ax[i] * dt * 0.5;
436 vy[i] += ay[i] * dt * 0.5;
437 vz[i] += az[i] * dt * 0.5;
438 }
439
440 /* 全推位置(不含边界)*/
441 for (int i = 0; i < n; i++) {
442 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
443 x[i] += vx[i] * dt; /* vx 此时是 v_half */
444 y[i] += vy[i] * dt;
445 z[i] += vz[i] * dt;
446 }
447
448 /* 隐式阻尼处理:v_next = (v_half + 0.5*g*dt) / (1 + 0.5*gamma*dt)
449 不覆盖 vx/vy/vz,用临时数组存入 */
450 double *dmp_vx = (double*)alloca(n * sizeof(double));
451 double *dmp_vy = (double*)alloca(n * sizeof(double));
452 double *dmp_vz = (double*)alloca(n * sizeof(double));
453 for (int i = 0; i < n; i++) {
454 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
455 double gx = B[0] / m[i], gy = B[1] / m[i], gz = B[2] / m[i];
456 dmp_vx[i] = (vx[i] + 0.5 * G[0] * dt) / (1.0 + 0.5 * gx * dt);
457 dmp_vy[i] = (vy[i] + 0.5 * G[1] * dt) / (1.0 + 0.5 * gy * dt);
458 dmp_vz[i] = (vz[i] + 0.5 * G[2] * dt) / (1.0 + 0.5 * gz * dt);
459 }
460
461 /* 用新位置 + 阻尼处理后的速度重算加速度 */
462 compute_acceleration(n, x, y, z, dmp_vx, dmp_vy, dmp_vz, m, G, B, bonds, ax, ay, az);
463
464 /* 速度后半步:v = v_half + 0.5*a_next*dt
465 vx 仍为 v_half(未被覆盖)*/
466 for (int i = 0; i < n; i++) {
467 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
468 vx[i] += ax[i] * dt * 0.5;
469 vy[i] += ay[i] * dt * 0.5;
470 vz[i] += az[i] * dt * 0.5;
471 }
472}
C 引擎的 leapfrog_step 与 Python 的 Leapfrog_Method 步骤完全相同:
dmp_vx(第448-459行)关键:阻尼值存入 dmp_vx 而不覆盖 vx,与 Python 的 vx_next 独立变量做法一致。
C309/* 边界条件:clamp 位置 + 速度反转 ——与 Python Limit_in_box 一致 */
310static void limit_in_box(double *pos, double *vel, double lo, double hi) {
311 if (*pos > hi) { *pos = hi; *vel = -*vel; }
312 if (*pos < lo) { *pos = lo; *vel = -*vel; }
313}
C473/* ── 分发器:调用对应积分方法 + 边界条件 ── */
488 } else if (strcmp(method, "leapfrog") == 0) {
489 leapfrog_step(n, x, y, z, vx, vy, vz, m, G, B, bonds, fixed, dt);
490 }
495 /* 边界条件(与 Python Limit_in_box 一致) */
496 for (int i = 0; i < n; i++) {
497 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
498 limit_in_box(&x[i], &vx[i], -box_a, box_a);
499 limit_in_box(&y[i], &vy[i], -box_a, box_a);
500 limit_in_box(&z[i], &vz[i], -box_a, box_a);
501 }
502}
原 C 引擎仅有 clamp(只钳位置不反转速度),现已改为 limit_in_box(钳位置 + 反转速度),与 Python 一致。
engines/cpp/main.cpp
C++360/* ── 蛙跳法(Velocity-Verlet)——与 Python Leapfrog_Method 一致 ── */
361static void leapfrog_full_step(
362 int n, double *x, double *y, double *z,
363 double *vx, double *vy, double *vz,
364 const double *m, const double G[3], const double B[3],
365 const BondData &bonds, const int *fixed, double dt)
366{
367 // 第一次加速度
368 std::vector ax(n), ay(n), az(n);
369 compute_acceleration(n, x, y, z, vx, vy, vz, m, G, B, bonds, ax.data(), ay.data(), az.data());
370
371 // 半推速度:v_half = v + 0.5*a*dt (存入 vx, vy, vz)
372 for (int i = 0; i < n; i++) {
373 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
374 vx[i] += ax[i] * dt * 0.5;
375 vy[i] += ay[i] * dt * 0.5;
376 vz[i] += az[i] * dt * 0.5;
377 }
378
379 // 全推位置(不含边界,边界在外层统一处理)
380 for (int i = 0; i < n; i++) {
381 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
382 x[i] += vx[i] * dt; // vx 此时是 v_half
383 y[i] += vy[i] * dt;
384 z[i] += vz[i] * dt;
385 }
386
387 // 隐式阻尼处理得到 v_next(不覆盖 vx/vy/vz,Python 第929-931行 vx_next)
388 std::vector dmp_vx(n), dmp_vy(n), dmp_vz(n);
389 for (int i = 0; i < n; i++) {
390 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
391 double gamma_x = B[0] / m[i];
392 double gamma_y = B[1] / m[i];
393 double gamma_z = B[2] / m[i];
394 dmp_vx[i] = (vx[i] + 0.5 * G[0] * dt) / (1.0 + 0.5 * gamma_x * dt);
395 dmp_vy[i] = (vy[i] + 0.5 * G[1] * dt) / (1.0 + 0.5 * gamma_y * dt);
396 dmp_vz[i] = (vz[i] + 0.5 * G[2] * dt) / (1.0 + 0.5 * gamma_z * dt);
397 }
398
399 // 用新位置 + 阻尼处理后的速度重算加速度
400 compute_acceleration(n, x, y, z, dmp_vx.data(), dmp_vy.data(), dmp_vz.data(),
401 m, G, B, bonds, ax.data(), ay.data(), az.data());
402
403 // 速度后半步:v = v_half + 0.5*a_next*dt
404 // vx 仍为 v_half(未被覆盖),直接加上 0.5*a_next*dt
405 for (int i = 0; i < n; i++) {
406 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
407 vx[i] += ax[i] * dt * 0.5;
408 vy[i] += ay[i] * dt * 0.5;
409 vz[i] += az[i] * dt * 0.5;
410 }
411}
C++265/* 边界条件:clamp 位置 + 速度反转 ——与 Python Limit_in_box 一致 */
266static void limit_in_box(double &pos, double &vel, double lo, double hi) {
267 if (pos > hi) { pos = hi; vel = -vel; }
268 if (pos < lo) { pos = lo; vel = -vel; }
269}
C++413/* ── 分发器:调用对应积分方法 + 边界条件 ── */
429 } else if (method == "leapfrog") {
430 leapfrog_full_step(n, x, y, z, vx, vy, vz, m, G, B, bonds, fixed, dt);
431 }
435 // 边界条件(与 Python Limit_in_box 一致)
436 for (int i = 0; i < n; i++) {
437 if (fixed[i*3+0] && fixed[i*3+1] && fixed[i*3+2]) continue;
438 limit_in_box(x[i], vx[i], -box_a, box_a);
439 limit_in_box(y[i], vy[i], -box_a, box_a);
440 limit_in_box(z[i], vz[i], -box_a, box_a);
441 }
442}
C++ 引擎的 leapfrog_full_step 在重写后已与 Python 完全对齐,阻尼使用 std::vector<double> dmp_vx 临时数组而非原处覆盖。边界条件新增了速度反转。
engines/fortran/main.f90
Fortran498! ── 蛙跳法(Velocity-Verlet)──
499subroutine leapfrog_full(n, x, y, z, vx, vy, vz, m, G, B, &
500 nb, bp, bk, br, fixed, dt)
501 integer, intent(in) :: n, nb, bp(nb, 2), fixed(n, 3)
502 double precision, intent(inout) :: x(n), y(n), z(n), vx(n), vy(n), vz(n)
503 double precision, intent(in) :: m(n), G(3), B(3), bk(nb), br(nb), dt
504 double precision :: ax(n), ay(n), az(n)
505 double precision :: dmp_vx(n), dmp_vy(n), dmp_vz(n)
506 double precision :: gx, gy, gz
507 integer :: i
508
509 call accel(n, x, y, z, vx, vy, vz, m, G, B, nb, bp, bk, br, ax, ay, az)
510
511 ! 速度半步推
512 do i = 1, n
513 if (fixed(i,1) /= 0 .and. fixed(i,2) /= 0 .and. fixed(i,3) /= 0) cycle
514 vx(i) = vx(i) + ax(i) * dt * 0.5d0
515 vy(i) = vy(i) + ay(i) * dt * 0.5d0
516 vz(i) = vz(i) + az(i) * dt * 0.5d0
517 end do
518
519 ! 全推位置(不含边界)
520 do i = 1, n
521 if (fixed(i,1) /= 0 .and. fixed(i,2) /= 0 .and. fixed(i,3) /= 0) cycle
522 x(i) = x(i) + vx(i) * dt
523 y(i) = y(i) + vy(i) * dt
524 z(i) = z(i) + vz(i) * dt
525 end do
526
527 ! 隐式阻尼处理(用临时数组 dmp_v,不覆盖 vx/vy/vz)
528 do i = 1, n
529 if (fixed(i,1) /= 0 .and. fixed(i,2) /= 0 .and. fixed(i,3) /= 0) then
530 dmp_vx(i) = 0; dmp_vy(i) = 0; dmp_vz(i) = 0; cycle
531 end if
532 gx = B(1) / m(i); gy = B(2) / m(i); gz = B(3) / m(i)
533 dmp_vx(i) = (vx(i) + 0.5d0 * G(1) * dt) / (1.0d0 + 0.5d0 * gx * dt)
534 dmp_vy(i) = (vy(i) + 0.5d0 * G(2) * dt) / (1.0d0 + 0.5d0 * gy * dt)
535 dmp_vz(i) = (vz(i) + 0.5d0 * G(3) * dt) / (1.0d0 + 0.5d0 * gz * dt)
536 end do
537
538 ! 用新位置 + 阻尼速度重算加速度
539 call accel(n, x, y, z, dmp_vx, dmp_vy, dmp_vz, m, G, B, nb, bp, bk, br, ax, ay, az)
540
541 ! 速度后半步:v = v_half + 0.5*a_next*dt(vx 仍为 v_half)
542 do i = 1, n
543 if (fixed(i,1) /= 0 .and. fixed(i,2) /= 0 .and. fixed(i,3) /= 0) cycle
544 vx(i) = vx(i) + ax(i) * dt * 0.5d0
545 vy(i) = vy(i) + ay(i) * dt * 0.5d0
546 vz(i) = vz(i) + az(i) * dt * 0.5d0
547 end do
548end subroutine leapfrog_full
Fortran395! 边界条件:clamp 位置 + 速度反转
396subroutine limit_in_box(pos, vel, lo, hi)
397 double precision, intent(inout) :: pos, vel
398 double precision, intent(in) :: lo, hi
399 if (pos > hi) then
400 pos = hi; vel = -vel
401 else if (pos < lo) then
402 pos = lo; vel = -vel
403 end if
404end subroutine limit_in_box
Fortran571 else if (trim(method) == 'leapfrog') then
572 call leapfrog_full(n, x, y, z, vx, vy, vz, mm, G, B, &
573 nb, bp, bk, br, fixed, dt)
574 end if
579 ! 边界条件
580 do i = 1, n
581 if (fixed(i,1) /= 0 .and. fixed(i,2) /= 0 .and. fixed(i,3) /= 0) cycle
582 call limit_in_box(x(i), vx(i), -box_a, box_a)
583 call limit_in_box(y(i), vy(i), -box_a, box_a)
584 call limit_in_box(z(i), vz(i), -box_a, box_a)
585 end do
586end subroutine apply_step
Fortran 引擎的 leapfrog_full 在重写中修复了阻尼覆盖问题(使用 dmp_vx 临时数组),边界也改为速度反转。
以下列出了 2026-05-20 三次迭代中所有修复和改进:
| # | 问题描述 | 影响引擎 | 修复方式 |
|---|---|---|---|
| 1 | C 引擎蛙跳算法错误 原 leapfrog_step 实现的是 Symplectic Euler:v += a·Δt; x += v·Δt(一次加速度,使用新速度推位置)。这比 Velocity-Verlet 的截断误差大一阶(O(Δt) vs O(Δt²)),且只用了一次加速度计算。 |
C | 重写为 velocity-Verlet 三步法: v½ → x → anew → v |
| 2 | 阻尼值覆盖速度变量 蛙跳法的隐式阻尼处理直接写回 vx[i],导致后面的速度后半步 vx[i] += ½·a·dt 在覆盖后的值上叠加。等效于重力加速度被多算了一次,引入 ½·g·Δt 的系统偏差(Fortran)或 g·Δt 的偏差(C)。 |
C, Fortran | 用临时数组 dmp_vx 存储阻尼估计值,不覆盖原始速度变量 |
| 3 | bond.txt 表头未跳过 C/C++ 引擎的 read_bonds 用 while (fb >> bn >> bk >> br) 读取 bond.txt,但第一行是表头 bond_name k rest_length。解析 k 时尝试将字符串 "k" 转 double 失败,使用默认值 1.0(应为 100.0),弹簧劲度系数差 100 倍,完全改变了轨迹。Python 和 Fortran 正确跳过了表头。 |
C, C++ | 在数据读取循环前用 std::getline(fb, header)(C++)或 fgets(header, ...)(C)跳过表头行 |
| 4 | C++ JSON 输出精度不足 C++ 的 JSON 输出用默认 << 运算符,只有 6 位有效数字,而 C 的 %.15g 和 Fortran 的 g0 都用全精度。导致 C++ 与其他引擎的输出到第 6 位之后截断。 |
C++ | 添加 std::setprecision(15) 到输出流 |
| 5 | C 引擎 record_steps 忽略 warmup_stepsC 引擎的轨迹缓冲区大小写死为 params.NT 而非 params.NT - params.warmup_steps。Python/C++/Fortran 正确。 |
C | 改为 int record_steps = params.NT - params.warmup_steps; |
| 6 | 力开关系统 新增 4 个独立力开关 gravity_field/gravity_interaction/elastic_force/damping_force,通过 input.txt 的 0/1 控制各力是否参与计算。所有引擎统一实现:清零 → 按条件累加。 |
四引擎 | Python: 重写 compute_force();C/C++: 重写 compute_acceleration();Fortran: 重写 accel() |
| 7 | 万有引力引擎 原子间牛顿万有引力 F = G · mi · mj / r²,O(N²) 双循环对计算, gravity_interaction=0 时零开销。 |
四引擎 | 所有引擎的受力计算函数中添加多原子对循环 |
| 8 | 逐自由度固定约束 coord.txt 末尾增加 fix_x fix_y fix_z 三列,每列 0/1 控制该方向是否锁定。锁定后位置回复初始值、速度置零。 |
四引擎 | C/C++/Fortran: apply_step() 新增 pos_0 参数,边界条件后添加约束;Python: 已有 apply_fixed_constraints() 完美支持 |
| 9 | 能量图随力开关联动 绘图时根据力开关状态动态计算势能分量。仅 gravity_interaction=1 时计算万有引力势能并以紫色曲线绘制。 |
Python | dynamics.py 能量计算添加 ug_grav 和开关判断 |
四引擎的边界条件在重写后完全一致:
| 引擎 | 函数名 | 钳位置 | 反转速度 | 调用时机 |
|---|---|---|---|---|
| Python | Limit_in_box | ✅ | ✅ | 蛙跳步之后(后处理) |
| C | limit_in_box | ✅ | ✅ | 积分步之后(apply_step) |
| C++ | limit_in_box | ✅ | ✅ | 积分步之后(apply_step) |
| Fortran | limit_in_box | ✅ | ✅ | 积分步之后(apply_step) |
速度反转是必须的。 硬壁反射的物理模型要求速度在碰撞时反向以保持运动方向的一致性。没有速度反转时,粒子到达边界后会被持续钳制在边界上("黏壁"现象),这是功能性缺陷。
2026-05-20 新增了逐自由度固定约束。在 coord.txt 末尾添加 fix_x fix_y fix_z 三列:
coord.txt 格式n mass radius x y z vx vy vz fix_x fix_y fix_z
1 1 0.28 -1 0 0 0 0 0 1 0 0 ← 仅固定 x 方向
2 1 0.28 1 0 1 0 0 0 0 1 1 ← 固定 y 和 z 方向
Python 的 apply_fixed_constraints() 使用 NumPy 逐元素操作:
Pythonfixed = ATOM_FIXED != 0 # (n_atoms, 3) 布尔矩阵
positions = np.where(fixed, ATOM_POSITIONS, positions) # 逐度约束
velocities = np.where(fixed, 0.0, velocities)
C/C++/Fortran 引擎在 apply_step() 的边界条件之后添加了相同的约束逻辑:
C/C++/Fortran 伪代码// 边界条件后,逐自由度固定约束
for i = 1..n:
if fixed_x[i]: x[i] = x0[i]; vx[i] = 0
if fixed_y[i]: y[i] = y0[i]; vy[i] = 0
if fixed_z[i]: z[i] = z0[i]; vz[i] = 0
约束在积分器完成位置/速度更新、边界钳位之后执行,确保固定自由度不参与任何物理演化。
Velocity-Verlet 是 symplectic 积分器,在平滑势场中能量几乎精确守恒(前 20 步漂移仅 -0.000078%)。但硬壁反射是一个非光滑势(无限高势垒),每个碰撞步的工作流程如下:
在典型 case01 参数下(100,000 步,83 次碰底),累积能量漂移约 36%。这是硬壁边界的固有限制,不是代码 Bug。
用光滑的排斥势替换硬壁 clamp 可保持 symplectic 守恒性:
| 维度 | Python | C | C++ | Fortran |
|---|---|---|---|---|
| 积分器类型 | Velocity-Verlet | Velocity-Verlet 原为 Symplectic Euler |
Velocity-Verlet | Velocity-Verlet |
| 算法数量 | 4 种 | 4 种 新增 3 种 |
4 种 新增 3 种 |
4 种 新增 3 种 |
| 阻尼处理 | 隐式 + 显式混合 | 同 Python | 同 Python | 同 Python |
| 边界:钳位置 | ✅ | ✅ | ✅ | ✅ |
| 边界:反转速度 | ✅ | ✅ 原缺失 |
✅ 原缺失 |
✅ 原缺失 |
| 输出精度(JSON) | 双精度 | %.15g | setprecision(15) 原默认 6 位 |
g0 |
| 算法选择方式 | YAML config | param.json | param.json | param.json |
| 引擎一致性 | 四个引擎输出完全一致(双精度验证通过) | |||
| 力开关系统 | ✅ | ✅ | ✅ | ✅ |
| 万有引力 | ✅ | ✅ | ✅ | ✅ |
| 逐自由度约束 | ✅ np.where |
✅ pos_0 参数 |
✅ pos_0 参数 |
✅ pos_0 参数 |