File: runge.cpp

package info (click to toggle)
pluto-find-orb 0.0~git20180227-2
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 2,668 kB
  • sloc: cpp: 30,743; makefile: 263
file content (1635 lines) | stat: -rw-r--r-- 64,770 bytes parent folder | download
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
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
/* runge.cpp: numerical integration code, mostly force model stuff

Copyright (C) 2010, Project Pluto

This program is free software; you can redistribute it and/or
modify it under the terms of the GNU General Public License
as published by the Free Software Foundation; either version 2
of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
GNU General Public License for more details.

You should have received a copy of the GNU General Public License
along with this program; if not, write to the Free Software
Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
02110-1301, USA.    */

#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <stdarg.h>
#include <math.h>
#include <assert.h>
#ifndef __GNUC__
// #include <conio.h>
// #include <io.h>
#else
#include <unistd.h>
#endif
#include "watdefs.h"
#include "lunar.h"
#include "afuncs.h"
#include "comets.h"
#include "afuncs.h"
#include "mpc_obs.h"

#define PI 3.1415926535897932384626433832795028841971693993751058209749445923
#define J2000 2451545.
#define GAUSS_K .01720209895
#define SOLAR_GM (GAUSS_K * GAUSS_K)

#define ldouble long double

#ifdef __WATCOMC__
#define sqrtl sqrt
#endif

/* perturbers
   excluded_perturbers = -1;
   general_relativity_factor
   integration_method
   integration_tolerance
   n_extra_params
   planet_mass[],  sort of
   solar_pressure
   best_fit_planet;
   best_fit_planet_dist;
   j2_multiplier
   debug_level
   approx_planet_orientation
   object_mass
   (implicitly) planet_posn cache
*/

double object_mass = 0.;
double j2_multiplier = 1.;
extern unsigned perturbers;

extern int debug_level;
int debug_printf( const char *format, ...)                 /* runge.cpp */
#ifdef __GNUC__
         __attribute__ (( format( printf, 1, 2)))
#endif
;
int planet_posn( const int planet_no, const double jd, double *vect_2000);
                                                /* mpc_obs.cpp */
int earth_lunar_posn( const double jd, double FAR *earth_loc, double FAR *lunar_loc);
const char *get_environment_ptr( const char *env_ptr);     /* mpc_obs.cpp */
double take_rk_step( const double jd, ELEMENTS *ref_orbit,
                 const double *ival, double *ovals,
                 const int n_vals, const double step);      /* runge.cpp */
double take_pd89_step( const double jd, ELEMENTS *ref_orbit,
                 const double *ival, double *ovals,
                 const int n_vals, const double step);      /* runge.cpp */
int symplectic_6( double jd, ELEMENTS *ref_orbit, double *vect,
                                          const double dt);
static int get_planet_posn_vel( const double jd, const int planet_no,
                     double *posn, double *vel);         /* runge.cpp */
int find_best_fit_planet( const double jd, const double *ivect,
                                 double *rel_vect);      /* runge.cpp */
int detect_perturbers( const double jd, const double * __restrict xyz,
                       double *accel);          /* bc405.cpp */
int find_relative_orbit( const double jd, const double *ivect,
               ELEMENTS *elements, const int ref_planet);     /* runge.cpp */
int parallax_to_lat_alt( const double rho_cos_phi, const double rho_sin_phi,
       double *lat, double *ht_in_meters, const int planet_idx); /* ephem0.c */
void calc_approx_planet_orientation( const int planet,        /* runge.cpp */
         const int system_number, const double jde, double *matrix);
double geo_potential_in_au( const double x, const double y, const double z,
                 double *derivs, const int n_terms);    /* geo_pot.c */

#define N_PERTURB 19
#define IDX_MERCURY    1
#define IDX_VENUS      2
#define IDX_EARTH      3
#define IDX_MARS       4
#define IDX_JUPITER    5
#define IDX_SATURN     6
#define IDX_URANUS     7
#define IDX_NEPTUNE    8
#define IDX_PLUTO      9
#define IDX_MOON      10
#define IDX_IO        11
#define IDX_EUROPA    12
#define IDX_GANYMEDE  13
#define IDX_CALLISTO  14
#define IDX_TETHYS    15
#define IDX_DIONE     16
#define IDX_RHEA      17
#define IDX_TITAN     18
#define IDX_IAPETUS   19
#define IDX_ASTEROIDS 20

      /* The following value for various J2s come from the _Explanatory  */
      /* Supplement_,  p 697.  They're in units of planetary radii,  and */
      /* must be converted to AU squared. */

#define J2_IN_EARTH_UNITS (1.08263e-3)
#define J3_IN_EARTH_UNITS (-2.54e-6)
#define J4_IN_EARTH_UNITS (-1.61e-6)
#define EARTH_R (6378.140 / AU_IN_KM)
#define EARTH_R2 (EARTH_R * EARTH_R)
#define EARTH_J2 (J2_IN_EARTH_UNITS * EARTH_R2)
#define EARTH_J3 (J3_IN_EARTH_UNITS * EARTH_R2 * EARTH_R)
#define EARTH_J4 (J4_IN_EARTH_UNITS * EARTH_R2 * EARTH_R2)

#define J2_IN_MARS_UNITS (1.964e-3)
#define MARS_R (3397.2 / AU_IN_KM)
#define MARS_J2 (J2_IN_MARS_UNITS*MARS_R * MARS_R)
#define MARS_J3 0.
#define MARS_J4 0.

#define SUN_R     (696000. / AU_IN_KM)
#define MERCURY_R   (2439. / AU_IN_KM)
#define VENUS_R     (6051. / AU_IN_KM)
#define PLUTO_R     (1500. / AU_IN_KM)
#define MOON_R      (1738. / AU_IN_KM)

/* Start considering atmospheric drag if within 500 km of earth: */

#define ATMOSPHERIC_LIMIT (EARTH_R + 500. / AU_IN_KM)
// #define ATMOSPHERIC_LIMIT 0

/* 2016 Jul 25:  updated all the following gas giant J2/3/4 data using */
/* data from http://ssd.jpl.nasa.gov/?gravity_fields_op */

#define J2_IN_JUPITER_UNITS (.01469562)
#define J3_IN_JUPITER_UNITS 0.
#define J4_IN_JUPITER_UNITS (-5.9131e-4)
#ifdef NOT_USING_ANYTHING_PAST_J4_YET
   #define J6_IN_JUPITER_UNITS   20.78e-6
   #define J6_IN_SATURN_UNITS    86.14e-6
   #define J8_IN_SATURN_UNITS   -10.e-6
#endif
#define JUPITER_R (71492. / AU_IN_KM)
#define JUPITER_R2 (JUPITER_R * JUPITER_R)
#define JUPITER_J2 (J2_IN_JUPITER_UNITS * JUPITER_R2)
#define JUPITER_J3 (J3_IN_JUPITER_UNITS * JUPITER_R2 * JUPITER_R)
#define JUPITER_J4 (J4_IN_JUPITER_UNITS * JUPITER_R2 * JUPITER_R2)

#define J2_IN_SATURN_UNITS (.01629071)
#define J4_IN_SATURN_UNITS (-.00093583)
#define SATURN_R (60330. / AU_IN_KM)
#define SATURN_R2 (SATURN_R * SATURN_R)
#define SATURN_J2 (J2_IN_SATURN_UNITS * SATURN_R2)
#define SATURN_J3 0.
#define SATURN_J4 (J4_IN_SATURN_UNITS * SATURN_R2 * SATURN_R2)

#define J2_IN_URANUS_UNITS  3510.68e-6
#define J4_IN_URANUS_UNITS   -34.17e-6
#define URANUS_R (25559. / AU_IN_KM)
#define URANUS_R2 (URANUS_R * URANUS_R)
#define URANUS_J2 (J2_IN_URANUS_UNITS * URANUS_R * URANUS_R)
#define URANUS_J3 0.
#define URANUS_J4 (J4_IN_URANUS_UNITS * URANUS_R2 * URANUS_R2)

#define J2_IN_NEPTUNE_UNITS  3408.43e-6
#define J4_IN_NEPTUNE_UNITS   -33.40e-6
#define NEPTUNE_R (25225. / AU_IN_KM)
#define NEPTUNE_R2 (NEPTUNE_R * NEPTUNE_R)
#define NEPTUNE_J2 (J2_IN_NEPTUNE_UNITS * NEPTUNE_R * NEPTUNE_R)
#define NEPTUNE_J3 0.
#define NEPTUNE_J4 (J4_IN_NEPTUNE_UNITS * NEPTUNE_R2 * NEPTUNE_R2)

   /* One can compute the acceleration due to the first (J2) oblateness
      term analytically,  and the following function does that.  However,
      the expressions for computing the accelerations due to J3 and J4
      get ugly,  and it becomes simpler to write a function for the _potential_
      and get the accelerations via numerical differentiation.  I'm keeping
      the numerical J2 acceleration version around just for reference. */

static double jn_potential( const double *loc, const double j3,
                                        const double j4)
{
   const double r2 = loc[0] * loc[0] + loc[1] * loc[1] + loc[2] * loc[2];
   const double r = sqrt( r2);
   const double mu = loc[2] / r;
   const double mu2 = mu * mu;
   const double p3 = mu * (2.5 * mu2 - 1.5);
   const double p4 = (35. * mu2 * mu2 - 30. * mu2 + 3.) / 8.;
// const double p2 = 1.5 * mu2 - .5;     /* Danby, p. 115 */

//    return( (j2 * p2 + j3 * p3 / r + j4 * p4 / r2) / (r * r2));
      return( (          j3 * p3 / r + j4 * p4 / r2) / (r * r2));
}

/* For an input planetocentric location in AU,  and a GM in AU^3/day^2,
computes the planetocentric acceleration due to J2,  J3,  and J4,  in
AU/day^2.  For the Earth,  the GGM03 model can be used;  see 'geo_pot.cpp'. */

static void numerical_gradient( double *grad, const double *loc,
                   const double planet_gm,
                   const double j2, const double j3, const double j4)
{
   const double r2 = loc[0] * loc[0] + loc[1] * loc[1] + loc[2] * loc[2];
   const double r = sqrt( r2), delta = r * 0.01;
   const double r7 = r2 * r2 * r2 * sqrt( r2);
   const double tval = (1.5 * r2 - 7.5 * loc[2] * loc[2]) / r7;
   int i;

   if( j3 == EARTH_J3)
      {
      static int n_terms = -9999;
                                    /* height in units of earth radii */
      if( n_terms == -9999)
         {
         n_terms = atoi( get_environment_ptr( "GEO_TERMS"));
         if( !n_terms)
            n_terms = 3;
         }
      if( n_terms > 0)     /* n_term <= 0 -> use "usual" J2 & J3 & J4 */
         {
         double ht_above_ground = (r / EARTH_R) - 1.;

         if( ht_above_ground < 0.04)         /* less than about 250 km */
            ht_above_ground = 0.04;
         geo_potential_in_au( loc[0], loc[1], loc[2], grad,
                      (int)( (double)n_terms / ht_above_ground) + 3);
         return;
         }
      }
   grad[0] = j2 * loc[0] * tval;   /* compute J2 accel analytically */
   grad[1] = j2 * loc[1] * tval;
   grad[2] = j2 * loc[2] * (3. * r2 / r7 + tval);
   for( i = 0; i < 3; i++)         /* then do J3 & J4 numerically */
      {
      double tval, tloc[3];
      int j;

      for( j = 0; j < 3; j++)
         tloc[j] = loc[j];
      tloc[i] += delta;
      tval = jn_potential( tloc, j3, j4);
      tloc[i] -= delta + delta;
      tval -= jn_potential( tloc, j3, j4);

      grad[i] += tval / (delta * 2.);
      grad[i] *= planet_gm;
      }
#ifdef OLD_DEBUGGING_CODE
   if( j3 == EARTH_J3)
      {
      double grad2[3];

      debug_printf( "Loc (in earth radii): %e %e %e; radius %e\n",
               loc[0] / EARTH_R,
               loc[1] / EARTH_R,
               loc[2] / EARTH_R, r / EARTH_R);
      debug_printf( "Gradient (numerical): %e %e %e\n",
                  grad[0], grad[1], grad[2]);

      geo_potential_in_au( loc[0], loc[1], loc[2], grad2,
                  atoi( get_environment_ptr( "GEO_TERMS")));

      debug_printf( "Gradient (analytica): %e %e %e\n",
                  grad2[0], grad2[1], grad2[2]);
      debug_printf( "Ratios                %e %e %e\n",
                  grad[0] / grad2[0],
                  grad[1] / grad2[1],
                  grad[2] / grad2[2]);
      }
#endif
}

   /* For testing purposes,  one can multiply the relativistic effect by */
   /* a constant 'general_relativity_factor'.  Set to zero,  this disables */
   /* GR;  set to one,  it mirrors the actual universe.  I've used this    */
   /* just to see how much GR matters for a given object.   */

double general_relativity_factor = 1.;

   /* Somewhat artificially,  we set the relativistic effect to be zero */
   /* when an object goes inside the sun (r < radius of the sun).  This */
   /* just avoids a situation wherein the relativistic acceleration     */
   /* explodes to infinity,  and the program drags to a crawl.  Similar */
   /* trickery takes place for J2, J3, and J4 when an object runs inside */
   /* a planet. */

static void set_relativistic_accel( double *accel, const double *posnvel)
{
   int i;
   const double c = AU_PER_DAY;           /* speed of light in AU per day */
   const double r_squared = posnvel[0] * posnvel[0] + posnvel[1] * posnvel[1]
                                                    + posnvel[2] * posnvel[2];
   const double v_squared = posnvel[3] * posnvel[3] + posnvel[4] * posnvel[4]
                                                    + posnvel[5] * posnvel[5];
   const double v_dot_r   = posnvel[0] * posnvel[3] + posnvel[1] * posnvel[4]
                                                    + posnvel[2] * posnvel[5];
   const double r = sqrt( r_squared), r_cubed_c_squared = r_squared * r * c * c;
#ifndef PREVIOUS_EQUATION
   const double r_component =
                  (4. * SOLAR_GM / r - v_squared) / r_cubed_c_squared;
   const double v_component = 4. * v_dot_r / r_cubed_c_squared;
#else
   const double v_component = 3. * v_dot_r / r_cubed_c_squared;
   const double r_component = 0.;
#endif

   if( r > SUN_R)
      for( i = 0; i < 3; i++)
         {
         accel[i] = r_component * posnvel[i] + v_component * posnvel[i + 3];
         accel[i] *= general_relativity_factor;
         }
   else                             /* shut off relativity inside the sun */
      for( i = 0; i < 3; i++)
         accel[i] = 0.;
}

/* This is explained on page 3 of

http://www.lpi.usra.edu/books/CometsII/7009.pdf

   'Cometary Orbit Determination and Nongravitational Forces',  D. K. Yeomans,
P. W. Chodas, G. Sitarski, S. Szutowicz, M. Krolikowska, "Comets II".

   Essentially,  the idea is that the following expression gives a decent
approximation to the magnitude of comet non-gravitational accelerations.
This 'g' function is then multiplied by a comet-specific,  fitted parameter
A1 to give the radial acceleration;  by A2 to give the transverse acceleration;
and by A3 to give the "out-of-orbit" component of the acceleration,  with A1
through A3 to be fitted parameters.  The resulting accelerations should
be in AU/day^2,  which fortunately matches the units used in this program.
(Note that MPC uses 10^-8 AU/day^2.  Patrick Rocher and JPL use AU/day^2 in
their comet elements.  Find_Orb should have a toggle for the units...)

   Frequently,  only A1 and A2 are used;  A3 often doesn't appear unless
one has three apparitions.  A3 has been added to the console version of
Find_Orb,  but it appears to be just about unnecessary.

   A small bit of trickery:  for the Marsden-Sekanina formula, g(1) = 1,
requiring a normalization constant alpha.  We get that on the first call
by setting alpha = 1,  determining what value we get for g(1),  and setting
alpha to the inverse of that.
*/

static inline double comet_g_func( const double r)
{
   if( object_type == OBJECT_TYPE_COMET)
      {                      /* default, Marsden/Sekanina formula */
      static double r0 = 2.808;         /* AU */
      static double m = 2.15;
      static double n = 5.093;
      static double k = 4.6142;
      static double alpha = 0.;
      double r_over_r0;

      if( !alpha)
         {
         const char *comet_params = get_environment_ptr( "COMET_CONSTANTS");

         sscanf( comet_params, "%lf,%lf,%lf,%lf", &r0, &m, &n, &k);
         alpha = 1.;
         alpha = 1. / comet_g_func( 1.);
         }
      r_over_r0 = r / r0;
      return( alpha * pow( r_over_r0, -m) * pow( 1. + pow( r_over_r0, n), -k));
      }
   else        /* just an inverse-square force,  used for 2009 BD */
      return( 1. / (r * r));
}

/* I wrote a little code to dump the following constants from DE-432.
They were given in AU^3/day^2;  from that,  I got the following,  using
the first column to fill most of the 'planet_mass' array.  Note that for
EMB through Uranus,  these are "system" masses including moons.

       mass(obj)/mass(sun)   mass(sun)/mass(obj)    GM (km^3/s^2)
Merc 1.660114153054348e-07 6.023682155592479e+06 2.203177969072598e+04
Venu 2.447838287784771e-06 4.085237186582997e+05 3.248585874397545e+05
EMB  3.040432648022641e-06 3.289005598102476e+05 4.035032298380295e+05
Mars 3.227156037554996e-07 3.098703590290707e+06 4.282837461279101e+04
Jupi 9.547919101886966e-04 1.047348630972762e+03 1.267127623546989e+08
Satu 2.858856727222416e-04 3.497901767786634e+03 3.794058466740400e+07
Uran 4.366249662744965e-05 2.290295052370693e+04 5.794556384409937e+06
Nept 5.151383772628673e-05 1.941225977597307e+04 6.836527004611366e+06
Plut 7.350487833457740e-09 1.360453921776768e+08 9.755011621830380e+02
Sun  1.000000000000000e+00 1.000000000000000e+00 1.327124381789709e+11
Eart 3.003489614792921e-06 3.329460488475656e+05 3.986004298243866e+05
Moon 3.694303322972000e-08 2.706870315119437e+07 4.902800013642883e+03
Cere 4.725582914451237e-10 2.116141052021147e+09 6.271436303937109e+01
Pall 1.018229468217943e-10 9.820968958501590e+09 1.351317153528802e+01
Juno 1.371898570800420e-11 7.289168611179182e+10 1.820680042651692e+00
Vest 1.302666122601934e-10 7.676564106868835e+09 1.728799972636488e+01
*/

#define MASS_EARTH_PLUS_MOON 3.040432648022641e-06
#define EARTH_MOON_MASS_RATIO 81.300568800524701
#define MASS_MOON (MASS_EARTH_PLUS_MOON / (1. + EARTH_MOON_MASS_RATIO))
#define MASS_EARTH (MASS_EARTH_PLUS_MOON - MASS_MOON)

#define MASS_OF_SUN_IN_KILOGRAMS 1.9891E+30

#define MASS_IO               (893.3E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_EUROPA           (479.7E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_GANYMEDE         (1482E+20  / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_CALLISTO         (1076E+20  / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_GALILEANS  (MASS_IO + MASS_EUROPA + MASS_GANYMEDE + MASS_CALLISTO)
#define MASS_JUPITER_SYSTEM  9.547919101886966e-04
#define MASS_JUPITER   (MASS_JUPITER_SYSTEM - MASS_GALILEANS)

   /* 2015 Jun 4: revised masses for Saturn's moons to Cassini values, */
   /* and added masses for Mimas,  Enceladus,  Hyperion,  and Phoebe.  */
   /* We aren't using those yet,  but if we do,  we now have them.     */
#define MASS_SATURN_SYSTEM    2.858856727222416e-04
#define MASS_MIMAS             (0.37493E+20  / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_ENCELADUS         (1.08022E+20  / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_TETHYS            (6.17449E+20  / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_DIONE             (10.95452E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_RHEA              (23.06518E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_TITAN            (1345.2E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_HYPERION          (0.056199E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_IAPETUS           (18.05635E+20 / MASS_OF_SUN_IN_KILOGRAMS)
#define MASS_PHOEBE            (0.08292E+20  / MASS_OF_SUN_IN_KILOGRAMS)
   /* "MASS_SATURN_SATS" = mass of those we're treating as perturbers: */
#define MASS_SATURN_SATS  (MASS_TETHYS + MASS_DIONE + MASS_RHEA + MASS_TITAN + MASS_IAPETUS)
#define MASS_SATURN          (MASS_SATURN_SYSTEM - MASS_SATURN_SATS)

      /* 22 Oct 2001:  replaced Uranus,  Neptune masses w/DE-405 values */
#define EARTH_MOON_RATIO 81.30056

extern const double planet_mass[N_PERTURB + 1] = { 1.,             /*  0 */
        1.660114153054348e-07,                /* mercury */        /*  1 */
        2.447838287784771e-06,                /* venus */          /*  2 */
        MASS_EARTH,                           /* Earth */          /*  3 */
        3.227156037554996e-07,                /* Mars */           /*  4 */
        MASS_JUPITER,                         /* Jupiter */        /*  5 */
        MASS_SATURN,                          /* saturn */         /*  6 */
        4.366249662744965e-05,                /* Uranus */         /*  7 */
        5.151383772628673e-05,                /* Neptune */        /*  8 */
        7.350487833457740e-09,                /* Pluto */          /*  9 */
        MASS_MOON,                            /* Moon */           /* 10 */

         MASS_IO,                                                    /* 11 */
         MASS_EUROPA,                                                /* 12 */
         MASS_GANYMEDE,                                              /* 13 */
         MASS_CALLISTO,                                              /* 14 */

         MASS_TETHYS,                                                /* 15 */
         MASS_DIONE,                                                 /* 16 */
         MASS_RHEA,                                                  /* 17 */
         MASS_TITAN,                                                 /* 18 */
         MASS_IAPETUS };                                             /* 19 */

static double planetary_system_mass( const int planet_no)
{
   double rval;
// const int bc405_start = 100;

   if( planet_no == 5)
      rval = MASS_JUPITER_SYSTEM;
   else if( planet_no == 6)
      rval = MASS_SATURN_SYSTEM;
// else if( planet_no == bc405_start)     /* Ceres: special fixed mass  */
//    rval = 4.76e-10;                    /* used for Dawn computations */
   else
      rval = planet_mass[planet_no];
   return( rval);
}

/* "approx_planet_orientation" returns the planet orientation matrix,  as
described in 'cospar.cpp',  for a time no more than range/2 days away.
(Which,  at present,  means half a day.)  The idea is that the orientation
of the planet's pole doesn't change much over that time scale,
and we don't want to do a lot of expensive recomputes of planet orientation
that won't have any noticeable effect.  (It's particularly expensive for
the earth,  where we do a full-meal-deal precession and nutation matrix.)

   For the earth,  we then rotate the planet to the time jd.  This will be
needed for handling non-zonal (tesseral) harmonics.  For the other planets,
we are (at least thus far) just dealing with the zonal terms J2, J3, J4,
which are spherically symmetrical around the planet's poles;  we just need
the polar orientation and can ignore rotation.  Though this may not be
true forever... at some point,  I may want to add tesseral harmonics for,
say,  the Moon or Mars.  If that happens,  this code will need revising,
perhaps with a table of omega values. */

void calc_approx_planet_orientation( const int planet,
         const int system_number, const double jde, double *matrix)
{
   static double cached_matrix[9];
   static double cached_jde = 0.;
   static int cached_planet = -1;
   const double range = 1.;
   const double new_jde = floor( jde / range + .5) * range;

   if( cached_planet != planet || new_jde != cached_jde)
      {
      const double ut = new_jde - td_minus_ut( new_jde) / seconds_per_day;

      calc_planet_orientation( planet, system_number, ut, cached_matrix);
      cached_planet = planet;
      cached_jde = new_jde;
      }
   memcpy( matrix, cached_matrix, 9 * sizeof( double));
   if( planet == 3)        /* rotate matrix */
      {
      const double omega = 360.9856235 * PI / 180.;
      const double theta = omega * (jde - cached_jde);

      spin_matrix( matrix, matrix + 3, theta);
      }
}

/* The idea of the following is as follows.  If our distance from the sun
is more than 120% of the planet's semimajor axis (given in the 'radii'
table),  we throw its entire mass into the sun.  (Assuming it's not
turned on as a perturber.)  If our distance is less than 100%
of the semimajor axis,  we don't add in anything.  In between,  we ramp
up the mass thrown into the sun linearly.  The reason for this is to avoid
a discontinuity,  and to instead have a gradual,  linear increase in the
mass we use for the sun.         */

static double include_thrown_in_planets( const double r)
{
   const int n_radii = 9;
   const double radii[9] = { 0., .38709927, .72333566, 1.00000261,
               1.52371034, 5.20288799, 9.53667594,  19.18916464,  30.06992276};
   const double fraction = .2;
   double rval = 1.;
   int i;

   for( i = 1; i < n_radii && r > radii[i]; i++)
      if( ((perturbers >> i) & 1) == 0)
         {
         if( r > radii[i] * (1. + fraction))
            rval += planet_mass[i];
         else
            rval += planet_mass[i] * ((1. + fraction) - r / radii[i]) / fraction;
         }
   return( rval);
}

/* Returns a passable estimate of the atmospheric density,  in kg/m^3,
as a function of height above sea level.  This just does an interpolation
within a table giving the atmospheric density,  at ten-kilometer intervals,
from 0 to 1000 km (extrapolating off the ends).  The table is from the
"Unofficial Australian Standard Atmosphere 2000",

http://www.sworld.com.au/steven/space/atmosphere/

   which is an update to the US Standard Atmosphere of 1976.  So far, only
the earth's atmophere is handled;  but note that similar tables are
available for Mars,  Venus,  Saturn,  and Titan :

https://solarsystem.nasa.gov/docs/03b_Justus_EDLatmospheres.pdf

   The interpolation is actually done within a table of the (natural)
logs of the original table.  That eliminates taking a logarithm at
runtime and makes the code slightly simpler.

   The line marked '(1)' ensures that,  for negative ht_in_km,  the
returned value extrapolates correctly,  contining to increase for a few
kilometers at roughly the correct rate (to handle objects hitting the Dead
Sea),  but then levels off and drops to a near-vacuum inside the earth.
That was the easiest way of avoiding problems with the density climbing to
absurd heights as one approached the center of the earth... doing this
avoids discontinuities that it would be awkward to code around.   */

static double atmospheric_density( const double ht_in_km)
{
#ifdef ORIGINAL_TABLE_FOR_REFERENCE_ONLY
   const double rho[] = {          1.22500e+0,  4.13510e-1,  8.89099e-2,
         1.84102e-2,  3.99568e-3,  1.02688e-3,  3.09678e-4,  8.28286e-5,
         1.84580e-5,  3.40454e-6,  5.60400e-7,  9.39781e-8,  2.22485e-8,
         8.15200e-9,  3.81812e-9,  2.07600e-9,  1.23365e-9, 7.81486e-10,
        5.19431e-10, 3.58062e-10, 2.54100e-10, 1.84452e-10, 1.36522e-10,
        1.02784e-10, 7.85246e-11, 6.07300e-11, 4.74445e-11, 3.73990e-11,
        2.97195e-11, 2.37875e-11, 1.91600e-11, 1.55181e-11, 1.26324e-11,
        1.03320e-11, 8.48753e-12, 7.00038e-12, 5.79500e-12, 4.81307e-12,
        4.00937e-12, 3.34858e-12, 2.80300e-12, 2.35089e-12, 1.97535e-12,
        1.66279e-12, 1.40213e-12, 1.18435e-12, 1.00205e-12, 8.49163e-13,
        7.20722e-13, 6.12626e-13, 5.21500e-13, 4.44559e-13, 3.79525e-13,
        3.24498e-13, 2.77890e-13, 2.38368e-13, 2.04815e-13, 1.76297e-13,
        1.52027e-13, 1.31346e-13, 1.13700e-13, 9.86238e-14, 8.57305e-14,
        7.46934e-14, 6.52351e-14, 5.71206e-14, 5.01507e-14, 4.41564e-14,
        3.89945e-14, 3.45434e-14, 3.07000e-14, 2.73763e-14, 2.44951e-14,
        2.19914e-14, 1.98103e-14, 1.79058e-14, 1.62389e-14, 1.47768e-14,
        1.34915e-14, 1.23593e-14, 1.13600e-14, 1.04763e-14, 9.69241e-15,
        8.99505e-15, 8.37279e-15, 7.81590e-15, 7.31609e-15, 6.86621e-15,
        6.46015e-15, 6.09261e-15, 5.75900e-15, 5.45526e-15, 5.17752e-15,
        4.92236e-15, 4.68680e-15, 4.46826e-15, 4.26448e-15, 4.07349e-15,
        3.89356e-15, 3.72317e-15, 3.56100e-15  };
#endif
   const double log_rho[] = { 0.2029, -0.8831, -2.4201, -3.9949,
         -5.5225, -6.8812, -8.0800, -9.3987, -10.9000, -12.5904,
        -14.3946, -16.1802, -17.6210, -18.6250, -19.3835, -19.9928,
        -20.5133, -20.9698, -21.3783, -21.7503, -22.0933, -22.4136,
        -22.7145, -22.9984, -23.2676, -23.5246, -23.7715, -24.0094,
        -24.2392, -24.4619, -24.6782, -24.8890, -25.0948, -25.2958,
        -25.4924, -25.6851, -25.8740, -26.0597, -26.2424, -26.4225,
        -26.6003, -26.7762, -26.9503, -27.1225, -27.2930, -27.4618,
        -27.6290, -27.7945, -27.9585, -28.1210, -28.2821, -28.4417,
        -28.5999, -28.7565, -28.9116, -29.0650, -29.2167, -29.3666,
        -29.5147, -29.6609, -29.8052, -29.9475, -30.0876, -30.2254,
        -30.3608, -30.4936, -30.6237, -30.7510, -30.8754, -30.9966,
        -31.1145, -31.2291, -31.3403, -31.4481, -31.5526, -31.6537,
        -31.7514, -31.8457, -31.9367, -32.0244, -32.1087, -32.1897,
        -32.2674, -32.3421, -32.4138, -32.4826, -32.5487, -32.6122,
        -32.6731, -32.7317, -32.7880, -32.8422, -32.8945, -32.9450,
        -32.9940, -33.0418, -33.0885, -33.1343, -33.1795, -33.2242,
        -33.2687 };
   double fraction, log_rval;
   int range = (int)( ht_in_km / 10.);
   const int table_size = sizeof( log_rho) / sizeof( log_rho[0]);  /* = 101 */

   if( range < 0)
      range = 0;
   if( range > table_size - 2)
      range = table_size - 2;
   fraction = ht_in_km / 10. - (double)range;
   log_rval = log_rho[range] + fraction * (log_rho[range + 1] - log_rho[range])
                - (fraction < 0. ? fraction * fraction : 0.);    /* (1) */
   if( !range && log_rval < -4.)
      log_rval = -4.;
   return( exp( log_rval));
}


      /* If we're closer than .1 AU,  include Galileans separately: */

#define GALILEAN_LIMIT .03
#define TITAN_LIMIT .03

unsigned excluded_perturbers = (unsigned)-1;
int best_fit_planet;
double best_fit_planet_dist;

/* The Earth and Moon pose a special problem in the following function.
The way we want things to work is this:  if the earth's perturbations are
included,  but not the moon's,  then we pretend that there is one object,
with the combined mass of the two,  at the earth-moon barycenter.  We've
got functions (JPL DE ephemerides and VSOP) to compute the EMB location,
so this is actually pretty straightforward.

   If Earth and Moon are handled as separate objects (always the case for
objects really close to us),  things are a little stickier.  When it comes
time for the Earth's position to be computed,  we call earth_lunar_posn(),
which actually computes both the earth _and_ lunar position simultaneously.
We store the latter,  and use it when it's time to compute lunar perturbations.

   Another peculiarity that should be explained:  in the real universe,
there are limits as to how close you can get to an object (except for a
black hole),  so accelerations cannot climb toward infinity.  But objects
do hit planets,  and in the preliminary steps of orbit determination,  you
may have orbits passing through planets.  To get around this,  I added
a fictional "compute_accel_multiplier".  If you're outside the planet,
this multiplier is 1,  i.e.,  no effect.  If you're within a fraction r0
radius of the center,  deep inside the planet (or sun or moon),  this
multiplier is zero,  i.e.,  there is no force (instead of acceleration
growing to infinity).

   In between,  a cubic spline linking the boundary conditions is used.
This ensures that both the accelerations and their derivatives will
be continuous.
*/

static double compute_accel_multiplier( double fraction)
{
   const double r0 = .8;  /* acceleration drops to zero at 80% of planet radius */
   double rval;

   assert( fraction >= 0. && fraction <= 1.);
   if( fraction < r0)
      rval = 0.;
   else
      {
      fraction = (fraction - r0) / (1. - r0);
      assert( fraction >= 0. && fraction <= 1.);
      rval =  fraction * fraction * (3. - 2. * fraction);
      }
   return( rval);
}

#define FUDGE_FACTOR .9

int planet_hit = -1;

int calc_derivatives( const double jd, const double *ival, double *oval,
                           const int reference_planet)
{
   double r, r2 = 0., solar_accel = 1. + object_mass;
   int i, j;
   unsigned local_perturbers = perturbers;
   double lunar_loc[3], jupiter_loc[3], saturn_loc[3];
   double relativistic_accel[3];
   extern int n_extra_params;
   static const double sphere_of_influence_radius[10] = {
            10000., 0.00075, 0.00412, 0.00618,   /* sun, mer, ven, ear */
            0.00386, 0.32229, 0.36466, 0.34606,  /* mar, jup, sat, ura */
            0.57928, 0.02208 };                  /* nep, plu */
   static const double planet_radius[11] = {
                        SUN_R * FUDGE_FACTOR, MERCURY_R * FUDGE_FACTOR,
                        VENUS_R * FUDGE_FACTOR, EARTH_R * FUDGE_FACTOR,
                        MARS_R * FUDGE_FACTOR, JUPITER_R * FUDGE_FACTOR,
                        SATURN_R * FUDGE_FACTOR, URANUS_R * FUDGE_FACTOR,
                        NEPTUNE_R * FUDGE_FACTOR, PLUTO_R * FUDGE_FACTOR,
                        MOON_R * FUDGE_FACTOR };

   assert( fabs( jd) < 1e+9);
   oval[0] = ival[3];
   oval[1] = ival[4];
   oval[2] = ival[5];
   best_fit_planet = 0;
   planet_hit = -1;
   for( i = 0; i < 3; i++)
      r2 += ival[i] * ival[i];
   r = sqrt( r2);
   if( n_extra_params == 1)  /* straightforward radiation pressure */
      {
      extern double solar_pressure[];

      solar_accel -= solar_pressure[0];
      }
   if( r < planet_radius[0])     /* special fudge to keep acceleration from reaching */
      {                          /* infinity inside the sun;  see above notes        */
      solar_accel *= compute_accel_multiplier( r / planet_radius[0]);
      planet_hit = 0;
      if( debug_level)
         debug_printf( "Inside the sun: %f km\n", r * AU_IN_KM);
      if( !solar_accel)
         {
         for( i = 3; i < 6; i++)
            oval[i] = 0.;
         return( 0);
         }
      }

   solar_accel *= -SOLAR_GM / (r2 * r);

   if( perturbers)
      set_relativistic_accel( relativistic_accel, ival);
   else                           /* shut off relativity if no perturbers */
      for( i = 0; i < 3; i++)
         relativistic_accel[i] = 0.;

   solar_accel *= include_thrown_in_planets( r);

   for( i = 0; i < 3; i++)
      oval[i + 3] = solar_accel * ival[i]
                 + SOLAR_GM * relativistic_accel[i];

   if( (local_perturbers >> IDX_ASTEROIDS) & 1)
      if( r < 11.5 && r > 1.)
         detect_perturbers( jd, ival, oval);        /* bc405.cpp */


   if( n_extra_params == 2 || n_extra_params == 3)
      {                  /* Marsden & Sekanina comet formula */
      const double g = comet_g_func( r);
      extern double solar_pressure[];
      double transverse[3], dot_prod = 0.;

      memcpy( transverse, ival + 3, 3 * sizeof( double));
      for( i = 0; i < 3; i++)
         dot_prod += transverse[i] * ival[i];
      for( i = 0; i < 3; i++)
         transverse[i] -= ival[i] * dot_prod / r;
      dot_prod = vector3_length( transverse);
      for( i = 0; i < 3; i++)
         oval[i + 3] += g * (solar_pressure[0] * ival[i] / r
                     + solar_pressure[1] * transverse[i] / dot_prod);
      if( n_extra_params == 3)
         {
         double out_of_plane[3];

         vector_cross_product( out_of_plane, ival, transverse);
         dot_prod = vector3_length( out_of_plane);
         for( i = 0; i < 3; i++)
            oval[i + 3] += g * solar_pressure[2] * out_of_plane[i] / dot_prod;
         }
      }

   if( perturbers)
      for( i = 1; i < N_PERTURB + 1; i++)
         if( ((local_perturbers >> i) & 1)
                   && !((excluded_perturbers >> i) & 1))
            {
            double planet_loc[15], accel[3], mass_to_use = planet_mass[i];
            double accel_multiplier = 1.;

            r = r2 = 0.;
            if( i >= IDX_IO)       /* Galileans,  Titan */
               {
               double matrix[10], sat_loc[15];
               const double t_years = (jd - J2000) / 365.25;

               if( i >= IDX_TETHYS)         /* Saturnian satell */
                  calc_ssat_loc( jd, sat_loc,
                               ((i == IDX_IAPETUS) ? 7 : i - 13), 0L);
               else
                  {
                  calc_jsat_loc( jd, sat_loc, 1 << (i - IDX_IO), 0L);
                  memmove( sat_loc, sat_loc + (i - IDX_IO) * 3,
                                                      3 * sizeof( double));
                  }
                                 /* turn ecliptic of date to equatorial: */
               rotate_vector( sat_loc, mean_obliquity( t_years / 100.), 0);
                                 /* then to equatorial J2000: */
               setup_precession( matrix, 2000. + t_years, 2000.);
               precess_vector( matrix, sat_loc, planet_loc + 12);
                                 /* then to ecliptic J2000: */
               equatorial_to_ecliptic( planet_loc + 12);
               for( j = 0; j < 3; j++)
                  {
                  double coord;

                  if( i >= IDX_TETHYS)         /* Saturnian */
                     coord = saturn_loc[j] + planet_loc[12 + j];
                  else
                     coord = jupiter_loc[j] + planet_loc[12 + j] * JUPITER_R;
                  r2 += coord * coord;
                  planet_loc[12 + j] = coord;
                  }
               planet_loc[2] = sqrt( r2);
               }
            else
               {
               if( local_perturbers & 1024)   /* if the moon is included */
                  {
                  if( i == 3)
                     earth_lunar_posn( jd, planet_loc, lunar_loc);
                  else if( i == 10)
                     memcpy( planet_loc, lunar_loc, 3 * sizeof( double));
                  else
                     planet_posn( i, jd, planet_loc);
                  }
               else
                  planet_posn( i, jd, planet_loc);

               for( j = 0; j < 3; j++)
                  r2 += planet_loc[j] * planet_loc[j];
               memcpy( planet_loc + 12, planet_loc, 3 * sizeof( double));
               planet_loc[2] = sqrt( r2);
               }

            for( j = 0; j < 3; j++)
               {
               accel[j] = ival[j] - planet_loc[12 + j];
               r += accel[j] * accel[j];
               }
            r = sqrt( r);

            if( i == IDX_JUPITER)
               {
               if( r < GALILEAN_LIMIT)
                  {
                  memcpy( jupiter_loc, planet_loc + 12, 3 * sizeof( double));
                  local_perturbers |= (15 << 11);
                  }
               else     /* "throw" Galileans into Jupiter: */
                  mass_to_use = MASS_JUPITER_SYSTEM;
               }

            if( i == IDX_SATURN)
               {
               if( r < TITAN_LIMIT)
                  {
                  memcpy( saturn_loc, planet_loc + 12, 3 * sizeof( double));
                  local_perturbers |= (31 << 15);
                  }
               else        /* "throw" saturn's satellites into the primary: */
                  mass_to_use = MASS_SATURN_SYSTEM;
               }

            if( r < planet_radius[i])
               {
               accel_multiplier = compute_accel_multiplier( r / planet_radius[i]);
               planet_hit = i;
               }
            if( i >= IDX_EARTH && i <= IDX_NEPTUNE && r < .015 && j2_multiplier
                                 && accel_multiplier)
               {          /* Within .015 AU,  we take J2 into account: */
               double grad[3], delta_j2000[3], matrix[10], delta_planet[3];
               const double j2[6] = { EARTH_J2, MARS_J2, JUPITER_J2,
                        SATURN_J2, URANUS_J2, NEPTUNE_J2 };
               const double j3[6] = { EARTH_J3, MARS_J3, JUPITER_J3,
                        SATURN_J3, URANUS_J3, NEPTUNE_J3 };
               const double j4[6] = { EARTH_J4, MARS_J4, JUPITER_J4,
                        SATURN_J4, URANUS_J4, NEPTUNE_J4 };
               const double total_j_mul = j2_multiplier * accel_multiplier;

               calc_approx_planet_orientation( i, 0, jd, matrix);
                           /* Remembering the 'accels' are 'deltas' now... */
               memcpy( delta_j2000, accel, 3 * sizeof( double));
                           /* Cvt ecliptic to equatorial 2000...: */
               ecliptic_to_equatorial( delta_j2000);
                           /* ...then to planet-centered coords: */
               precess_vector( matrix, delta_j2000, delta_planet);

               if( total_j_mul)
                  numerical_gradient( grad, delta_planet,
                        planet_mass[i] * SOLAR_GM,
                        j2[i - 3], j3[i - 3], j4[i - 3]);
               else     /* inside the planet */
                  for( j = 0; j < 3; j++)
                     grad[j] = 0.;
                        /* Cvt gradient from planet-centric to equatorial: */
               deprecess_vector( matrix, grad, delta_j2000);
                           /* Cvt equatorial to ecliptic: */
               equatorial_to_ecliptic( delta_j2000);
                           /* And add 'em to the output acceleration: */
               for( j = 0; j < 3; j++)
                  oval[j + 3] -= total_j_mul * delta_j2000[j];
               if( i == IDX_EARTH && r < ATMOSPHERIC_LIMIT && n_extra_params == 1
                           && !*get_environment_ptr( "DRAG_SHUTOFF"))
                  {
                  extern double solar_pressure[];
                  const double SRP1AU = 2.3e-7;   /* kg*AU^3 / (m^2*d^2) */
                  const double amr_drag = solar_pressure[0] * SOLAR_GM / SRP1AU;
                  const double rho_cos_phi = sqrt( delta_planet[0] * delta_planet[0]
                              + delta_planet[1] * delta_planet[1]) / EARTH_R;
                  const double rho_sin_phi = delta_planet[2] / EARTH_R;
                  double ht_in_meters, rho;
                  const double meters_per_km = 1000.;
                  const double dt = .5 / minutes_per_day;  /* half a minute */
                  const double Cd = 0.47;  /* drag coeff for a sphere */
                  double earth_loc[3];  /* 2nd location of earth, used for getting vel */
                  double vel[3];       /* velocity of obj relative to the earth */
                  double speed;        /* magnitude of the vel[] vector */
                  double drag[3];     /* drag acceleration,  in m/s^2 */
                  double accel_coeff;

                  parallax_to_lat_alt( rho_cos_phi, rho_sin_phi, NULL,
                                    &ht_in_meters, i);
//                debug_printf( "Parallax: %f %f; ht in meters %f\n",
//                               rho_cos_phi,  rho_sin_phi, ht_in_meters);
                  rho = atmospheric_density( ht_in_meters / meters_per_km);
                  earth_lunar_posn( jd + dt, earth_loc, NULL);
                  for( j = 0; j < 3; j++)
                     {
                                    /* in degrees/day,  sidereal... */
                     const double earth_rotation_rate = 360.98564736629;
                                    /* ...and converted to radians/day:    */
                     const double earth_omega = earth_rotation_rate * (PI / 180.);
                     const double earth_vel = (earth_loc[j] - planet_loc[j + 12]) / dt;
                     const double topo_vel = earth_omega *
                             (delta_planet[0] * matrix[j + 3] - delta_planet[1] * matrix[j]);

                     vel[j] = earth_vel - oval[j];       /* in AU/day */
                     vel[j] -= topo_vel * accel_multiplier;
                     vel[j] *= AU_IN_METERS / seconds_per_day;  /* cvt to m/s */
                     }
                  speed = vector3_length( vel);    /* also in m/s */
                  accel_coeff = rho * Cd * speed * amr_drag / 2.;
                  for( j = 0; j < 3; j++)
                     drag[j] = vel[j] * accel_coeff;
//                debug_printf( "Height %7.2f km; speed %7.1f m/s; accel %6.1f m/s^2\n",
//                            ht_in_meters / meters_per_km, speed, vector3_length( drag));
                  for( j = 0; j < 3; j++)
                     oval[j + 3] += drag[j] * seconds_per_day * seconds_per_day / AU_IN_METERS;
                  }
               }

                     /* If we're including the earth,  but the moon isn't */
                     /* being included separately,  add the mass of the moon:*/
            if( i == IDX_EARTH)
               if( !((local_perturbers >> IDX_MOON) & 1) ||
                    ((excluded_perturbers >> IDX_MOON) & 1))
                  mass_to_use += planet_mass[IDX_MOON];


            if( i < 10 && r < sphere_of_influence_radius[i]
                                && reference_planet != -1)
               {
               best_fit_planet = i;
               best_fit_planet_dist = r;
               }

            if( accel_multiplier)
               {
               const double accel_factor = -SOLAR_GM * accel_multiplier
                                         * mass_to_use / (r * r * r);

               for( j = 0; j < 3; j++)
                  oval[j + 3] += accel_factor * accel[j];
               }
            if( i != reference_planet)
               {
               if( reference_planet >= 0)
                  {
                  double planet_posn[3];

                  get_planet_posn_vel( jd, reference_planet, planet_posn, NULL);
                  for( j = 0; j < 3; j++)
                     planet_loc[j + 12] -= planet_posn[j];
                  r = vector3_length( planet_loc + 12);
                  }
               else
                  r = planet_loc[2];
               if( r)
                   r = -SOLAR_GM * mass_to_use / (r * r * r);
               for( j = 0; j < 3; j++)
                  oval[j + 3] += r * planet_loc[j + 12];
               }
            else        /* subtract sun's attraction to reference planet */
               {
               r = planet_loc[2];
               if( r)
                   r = SOLAR_GM / (r * r * r);
               for( j = 0; j < 3; j++)
                  oval[j + 3] += r * planet_loc[j + 12];
               }
            }
   return( planet_hit);
}

int calc_derivativesl( const ldouble jd, const ldouble *ival, ldouble *oval,
                           const int reference_planet)
{
   unsigned i;
   int rval;
   double ival1[6], oval1[6];

   assert( fabs( (double)jd) < 1e+9);
   for( i = 0; i < 6; i++)
      ival1[i] = (double)ival[i];
   rval = calc_derivatives( (double)jd, ival1, oval1, reference_planet);
   for( i = 0; i < 6; i++)
      oval[i] = (ldouble)oval1[i];
   return( rval);
}

static int get_planet_posn_vel( const double jd, const int planet_no,
                     double *posn, double *vel)
{
   assert( fabs( jd) < 1e+9);
   if( posn)
      {
      if( !planet_no)       /* sun doesn't move in the heliocentric frame */
         memset( posn, 0, 3 * sizeof( double));
      else if( planet_no == 3)
         earth_lunar_posn( jd, posn, NULL);
      else if( planet_no == 10)
         earth_lunar_posn( jd, NULL, posn);
      else
         planet_posn( planet_no, jd, posn);
      }
   if( vel)
      {
      if( !planet_no)       /* sun doesn't move in the heliocentric frame */
         memset( vel, 0, 3 * sizeof( double));
      else
         {
         int i;
         double loc1[3], loc2[3];
         const double delta = 1. / minutes_per_day;    /* one minute delta... */

         get_planet_posn_vel( jd + delta, planet_no, loc2, NULL);
         get_planet_posn_vel( jd - delta, planet_no, loc1, NULL);
         for( i = 0; i < 3; i++)
            vel[i] = (loc2[i] - loc1[i]) / (2. * delta);
         }
      }
   return( 0);
}

int find_relative_orbit( const double jd, const double *ivect,
               ELEMENTS *elements, const int ref_planet)
{
   double local_rel_vect[6];

   assert( ref_planet >= 0);
   memcpy( local_rel_vect, ivect, 6 * sizeof( double));
   if( ref_planet)
      {
      double planet_state[6];
      int i;

      get_planet_posn_vel( jd, ref_planet, planet_state, planet_state + 3);
      for( i = 0; i < 6; i++)
         local_rel_vect[i] -= planet_state[i];
      }
// if( rel_vect)
//    memcpy( rel_vect, local_rel_vect, 6 * sizeof( double));
   if( elements)
      {
      elements->gm = SOLAR_GM * planetary_system_mass( ref_planet);
      calc_classical_elements( elements, local_rel_vect, jd, 1);
      elements->epoch = jd;
      elements->central_obj = ref_planet;
      }
   return( 0);
}

int check_for_perturbers( const double t_cen, const double *vect); /* sm_vsop*/

int find_best_fit_planet( const double jd, const double *ivect,
                                 double *rel_vect)
{
   int i, j, rval = 0;
   double best_fit = 0., vel[3];
   int included = 1;          /* the Sun is always included in this */
   const int possible_perturber = check_for_perturbers( (jd - J2000) / 36525.,
                                             ivect);

   assert( ivect);
   assert( rel_vect);
   included |= (1 << possible_perturber);
   if( possible_perturber == 3)
      included |= (1 << 10);
   for( i = 0; i < 11; i++)
      if( (included >> i) & 1)
         {
         double planet_loc[3], delta[3], dist2 = 0., curr_fit;

         get_planet_posn_vel( jd, i, planet_loc, NULL);
         for( j = 0; j < 3; j++)
            {
            delta[j] = ivect[j] - planet_loc[j];
            dist2 += delta[j] * delta[j];
            }

         curr_fit = planet_mass[i] / exp( log( dist2) * 1.25);
         if( !i || curr_fit > best_fit)
            {
            rval = i;
            best_fit = curr_fit;
            memcpy( rel_vect, delta, 3 * sizeof( double));
            }
         }

   get_planet_posn_vel( jd, rval, NULL, vel);
   for( j = 0; j < 3; j++)
      rel_vect[j + 3] = ivect[j + 3] - vel[j];
   return( rval);
}

static void compute_ref_state( ELEMENTS *ref_orbit, double *ref_state,
                                          const double jd)
{
   double r2 = 0., accel;
   int i;

   if( ref_orbit->central_obj < 0)      /* using the method of Cowell */
      {
      memset( ref_state, 0, 9 * sizeof( double));
      return;
      }

   comet_posn_and_vel( ref_orbit, jd, ref_state, ref_state + 3);
   for( i = 0; i < 3; i++)
      r2 += ref_state[i] * ref_state[i];
   accel = -SOLAR_GM * planetary_system_mass( ref_orbit->central_obj) / (r2 * sqrt( r2));
   for( i = 0; i < 3; i++)
      ref_state[i + 6] = accel * ref_state[i];

   if( ref_orbit->central_obj)
      {
      double planet_state[6];

      get_planet_posn_vel( jd, ref_orbit->central_obj, planet_state, planet_state + 3);
      for( i = 0; i < 6; i++)
         ref_state[i] += planet_state[i];
      }
}

#define B_2_1     .0555555555555555555555555555556
#define B_3_1     .0208333333333333333333333333333
#define B_3_2     .0625
#define B_4_1     .03125
#define B_4_2    0.
#define B_4_3     .09375
#define B_5_1     .3125
#define B_5_2    0.
#define B_5_3   -1.171875
#define B_5_4    1.171875
#define B_6_1     .0375
#define B_6_2    0.
#define B_6_3    0.
#define B_6_4     .1875
#define B_6_5     .15
#define B_7_1     .0479101371111111111111111111111
#define B_7_2    0.
#define B_7_3    0.
#define B_7_4     .112248712777777777777777777778
#define B_7_5    -.0255056737777777777777777777778
#define B_7_6     .0128468238888888888888888888889
#define B_8_1        .016917989787292281181431107136
#define B_8_2       0.
#define B_8_3       0.
#define B_8_4       .387848278486043169526545744159
#define B_8_5       .0359773698515003278967008896348
#define B_8_6       .196970214215666060156715256072
#define B_8_7      -.172713852340501838761392997002
#define B_9_1       .0690957533591923006485645489846
#define B_9_2      0.
#define B_9_3      0.
#define B_9_4      -.634247976728854151882807874972
#define B_9_5      -.161197575224604080366876923982
#define B_9_6       .138650309458825255419866950133
#define B_9_7       .94092861403575626972423968413
#define B_9_8       .211636326481943981855372117132
#define B_10_1      .183556996839045385489806023537
#define B_10_2     0.
#define B_10_3     0.
#define B_10_4    -2.46876808431559245274431575997
#define B_10_5     -.291286887816300456388002572804
#define B_10_6     -.026473020233117375688439799466
#define B_10_7     2.84783876419280044916451825422
#define B_10_8      .281387331469849792539403641827
#define B_10_9      .123744899863314657627030212664
#define B_11_1    -1.21542481739588805916051052503
#define B_11_2     0.
#define B_11_3     0.
#define B_11_4    16.6726086659457724322804132886
#define B_11_5      .915741828416817960595718650451
#define B_11_6    -6.05660580435747094755450554309
#define B_11_7   -16.0035735941561781118417064101
#define B_11_8    14.849303086297662557545391898
#define B_11_9   -13.3715757352898493182930413962
#define B_11_10    5.13418264817963793317325361166
#define B_12_1      .258860916438264283815730932232
#define B_12_2     0.
#define B_12_3     0.
#define B_12_4    -4.77448578548920511231011750971
#define B_12_5     -.43509301377703250944070041181
#define B_12_6    -3.04948333207224150956051286631
#define B_12_7     5.57792003993609911742367663447
#define B_12_8     6.15583158986104009733868912669
#define B_12_9    -5.06210458673693837007740643391
#define B_12_10    2.19392617318067906127491429047
#define B_12_11     .134627998659334941535726237887
#define B_13_1      .822427599626507477963168204773
#define B_13_2     0.
#define B_13_3     0.
#define B_13_4   -11.6586732572776642839765530355
#define B_13_5     -.757622116690936195881116154088
#define B_13_6      .713973588159581527978269282765
#define B_13_7    12.0757749868900567395661704486
#define B_13_8    -2.12765911392040265639082085897
#define B_13_9     1.99016620704895541832807169835
#define B_13_10    -.234286471544040292660294691857
#define B_13_11     .17589857770794226507310510589
#define B_13_12    0.

/*  The coefficients BHAT(*) refer to the formula used to advance the
   integration, here the one of order 8.  The coefficients B(*) refer
   to the other formula, here the one of order 7.  */

#define CHAT_1     .0417474911415302462220859284685
#define CHAT_2    0.
#define CHAT_3    0.
#define CHAT_4    0.
#define CHAT_5    0.
#define CHAT_6    -.0554523286112393089615218946547
#define CHAT_7     .239312807201180097046747354249
#define CHAT_8     .70351066940344302305804641089
#define CHAT_9    -.759759613814460929884487677085
#define CHAT_10    .660563030922286341461378594838
#define CHAT_11    .158187482510123335529614838601
#define CHAT_12   -.238109538752862804471863555306
#define CHAT_13    .25

#define C_1      .029553213676353496981964883112
#define C_2     0.
#define C_3     0.
#define C_4     0.
#define C_5     0.
#define C_6     -.828606276487797039766805612689
#define C_7      .311240900051118327929913751627
#define C_8     2.46734519059988698196468570407
#define C_9    -2.54694165184190873912738007542
#define C_10    1.44354858367677524030187495069
#define C_11     .0794155958811272872713019541622
#define C_12     .0444444444444444444444444444445
#define C_13    0.

#define A_1     0.
#define A_2      .0555555555555555555555555555556
#define A_3      .0833333333333333333333333333334
#define A_4      .125
#define A_5      .3125
#define A_6      .375
#define A_7      .1475
#define A_8      .465
#define A_9      .564865451382259575398358501426
#define A_10     .65
#define A_11     .924656277640504446745013574318
#define A_12    1.
#define A_13      A_12

#define N_EVALS 13
#define N_EVALS_PLUS_ONE 14

double take_pd89_step( const double jd, ELEMENTS *ref_orbit,
                 const double *ival, double *ovals,
                 const int n_vals, const double step)
{
   double *ivals[N_EVALS_PLUS_ONE], *ivals_p[N_EVALS], rval = 0.;
   int i, j, k;
   const double bvals[91] = { B_2_1,
       B_3_1, B_3_2,
       B_4_1, B_4_2, B_4_3,
       B_5_1, B_5_2, B_5_3, B_5_4,
       B_6_1, B_6_2, B_6_3, B_6_4, B_6_5,
       B_7_1, B_7_2, B_7_3, B_7_4, B_7_5, B_7_6,
       B_8_1, B_8_2, B_8_3, B_8_4, B_8_5, B_8_6, B_8_7,
       B_9_1, B_9_2, B_9_3, B_9_4, B_9_5, B_9_6, B_9_7, B_9_8,
       B_10_1, B_10_2, B_10_3, B_10_4, B_10_5, B_10_6, B_10_7, B_10_8, B_10_9,
       B_11_1, B_11_2, B_11_3, B_11_4, B_11_5, B_11_6, B_11_7, B_11_8, B_11_9, B_11_10,
       B_12_1, B_12_2, B_12_3, B_12_4, B_12_5, B_12_6, B_12_7, B_12_8, B_12_9, B_12_10, B_12_11,
       B_13_1, B_13_2, B_13_3, B_13_4, B_13_5, B_13_6, B_13_7, B_13_8, B_13_9, B_13_10, B_13_11, B_13_12,
       CHAT_1, CHAT_2, CHAT_3, CHAT_4, CHAT_5, CHAT_6, CHAT_7, CHAT_8, CHAT_9, CHAT_10, CHAT_11, CHAT_12, CHAT_13 };
   const double avals[N_EVALS_PLUS_ONE] = { 0, A_1, A_2, A_3, A_4, A_5,
             A_6, A_7, A_8, A_9, A_10, A_11, A_12, A_13 };
   const double *bptr = bvals;

   ivals[0] = (double *)calloc( (2 * N_EVALS + 1) * n_vals, sizeof( double));
   assert( ivals[0]);
   if( !ivals[0])
      return( 0.);
   for( i = 0; i < N_EVALS; i++)
      {
      ivals[i + 1] = ivals[0] + (i + 1) * n_vals;
      ivals_p[i] = ivals[0] + (i + N_EVALS + 1) * n_vals;
      }

   for( j = 0; j <= N_EVALS; j++)
      {
      double ref_state_j[9], state_j[6];
      const double jd_j = jd + step * avals[j];

      compute_ref_state( ref_orbit, ref_state_j, jd_j);
      if( !j)
         {
         memcpy( state_j, ival, 6 * sizeof( double));
               /* subtract the analytic posn/vel from the numeric: */
         for( i = 0; i < n_vals; i++)
            ivals[0][i] = ival[i] - ref_state_j[i];
         }
      else
         for( i = 0; i < n_vals; i++)
            {
            double tval = 0.;

            for( k = 0; k < j; k++)
               tval += bptr[k] * ivals_p[k][i];
            ivals[j][i] = tval * step + ivals[0][i];
            state_j[i] = ivals[j][i] + ref_state_j[i];
            }
      bptr += j;
      if( j != N_EVALS)
         {
         assert( fabs( jd_j) < 1e+9);
         calc_derivatives( jd_j, state_j, ivals_p[j], ref_orbit->central_obj);
         for( k = 0; k < 6; k++)
            ivals_p[j][k] -= ref_state_j[k + 3];
         }
      else     /* on last iteration,  we have our answer: */
         memcpy( ovals, state_j, n_vals * sizeof( double));
      }

   for( i = 0; i < n_vals; i++)
      {
      double tval = 0.;
      const double err_coeff[N_EVALS] = { CHAT_1 - C_1, CHAT_2 - C_2,
               CHAT_3  -  C_3, CHAT_4  -  C_4, CHAT_5  -  C_5, CHAT_6 - C_6,
               CHAT_7  -  C_7, CHAT_8  -  C_8, CHAT_9  -  C_9, CHAT_10 - C_10,
               CHAT_11 - C_11, CHAT_12 - C_12, CHAT_13 - C_13 };

      for( k = 0; k < N_EVALS; k++)
         tval += err_coeff[k] * ivals_p[k][i];
      rval += tval * tval;
      }
   return( sqrt( rval * step * step));
}

#define ORIGINAL_FEHLBERG_CONSTANTS

         /* These "original" constants can be found in Danby, p. 298.  */
         /* The ones actually used are from _Numerical Recipes_,  and  */
         /* are for the Cash-Karp variant of Runge-Kutta.  NR says      */
         /* these constants are slightly better.  (Must admit,  I've not */
         /* done a really careful comparison!)  Dormand-Prince's method */
         /* (the 5/4 order version,  not the 9/8 one given above) may also */
         /* be worth a look.   */
#ifdef ORIGINAL_FEHLBERG_CONSTANTS
#define RKF_B21 2. / 9.
#define RKF_B31 1. / 12.
#define RKF_B32 1. / 4.
#define RKF_B41 69. / 128.
#define RKF_B42 -243. / 128.
#define RKF_B43 135. / 64.
#define RKF_B51 -17. / 12.
#define RKF_B52 27. / 4.
#define RKF_B53 -27. / 5.
#define RKF_B54 16. / 15.
#define RKF_B61 65. / 432.
#define RKF_B62 -5. / 16.
#define RKF_B63 13 / 16.
#define RKF_B64 4 / 27.
#define RKF_B65 5. / 144.
#define RKF_CHAT1 47. / 450.
#define RKF_CHAT2 0.
#define RKF_CHAT3 12 / 25.
#define RKF_CHAT4 32. / 225.
#define RKF_CHAT5 1. / 30.
#define RKF_CHAT6 6. / 25.
#define RKF_C1  1. / 9.
#define RKF_C2 0.
#define RKF_C3 9. / 20.
#define RKF_C4 16. / 45.
#define RKF_C5 1. / 12.
#define RKF_C6 0.
#define RKF_A1       0.
#define RKF_A2       2. / 9.
#define RKF_A3       1. / 3.
#define RKF_A4       .75
#define RKF_A5       1.
#define RKF_A6       5. / 6.
#else
#define RKF_B21     1. / 5.
#define RKF_B31     3. / 40.
#define RKF_B32     9. / 40.
#define RKF_B41     3. / 10.
#define RKF_B42    -9. / 10.
#define RKF_B43     6. / 5.
#define RKF_B51   -11. / 54.
#define RKF_B52     5. / 2.
#define RKF_B53   -70. / 27.
#define RKF_B54    35. / 27.
#define RKF_B61  1631. / 55296
#define RKF_B62   175. / 512.
#define RKF_B63   575. / 13824.
#define RKF_B64 44275. / 110592.
#define RKF_B65   253. / 4096.
#define RKF_CHAT1  2825. / 27648.
#define RKF_CHAT2 0.
#define RKF_CHAT3 18575. / 48384.
#define RKF_CHAT4 13525. / 55296.
#define RKF_CHAT5 277. / 14336.
#define RKF_CHAT6 .25
#define RKF_C1    37. / 378.
#define RKF_C2            0.
#define RKF_C3   250. / 621.
#define RKF_C4   125. / 594.
#define RKF_C5            0.
#define RKF_C6  512. / 1771.
#define RKF_A1           0.
#define RKF_A2            .2
#define RKF_A3            .3
#define RKF_A4            .6
#define RKF_A5           1.
#define RKF_A6            .875
#endif
/*   Butcher tableau looks like this :
A1 |
A2 | B21
A3 | B31 B32
A4 | B41 B42 B43
A5 | B51 B52 B53 B54
A6 | B61 B62 B63 B64 B65
---+---------------------
   | C1  C2  C3  C4  C5  C6
   | C^1 C^2 C^3 C^4 C^5 C^6          */

ldouble take_rk_stepl( const ldouble jd, ELEMENTS *ref_orbit,
                 const ldouble *ival, ldouble *ovals,
                 const int n_vals, const ldouble step)
{
   ldouble *ivals[7], *ivals_p[6], rval = 0.;
   int i, j, k;
            /* Revised values from _Numerical Recipes_: */
   const ldouble bvals[21] = { RKF_B21,
            RKF_B31, RKF_B32,
            RKF_B41, RKF_B42, RKF_B43,
            RKF_B51, RKF_B52, RKF_B53, RKF_B54,
            RKF_B61, RKF_B62, RKF_B63, RKF_B64, RKF_B65,
            RKF_CHAT1, RKF_CHAT2, RKF_CHAT3,
            RKF_CHAT4, RKF_CHAT5, RKF_CHAT6 };

   const ldouble avals[7] = { RKF_A1, RKF_A2, RKF_A3, RKF_A4, RKF_A5, RKF_A6, 1.};
   const ldouble *bptr = bvals;
   ldouble temp_ivals[78];

   if( n_vals > 6)
      ivals[0] = (ldouble *)calloc( 13 * n_vals, sizeof( ldouble));
   else
      ivals[0] = temp_ivals;
   assert( ivals[0]);
   if( !ivals[0])
      return( 0.);
   for( i = 0; i < 6; i++)
      {
      ivals[i + 1] = ivals[0] + (i + 1) * n_vals;
      ivals_p[i] = ivals[0] + (i + 7) * n_vals;
      }

   for( j = 0; j < 7; j++)
      {
      ldouble ref_state_j[9], state_j[6];
      const ldouble jd_j = jd + step * avals[j];
      double temp_array[9];

      compute_ref_state( ref_orbit, temp_array, (double)jd_j);
      for( i = 0; i < 9; i++)
         ref_state_j[i] = (ldouble)temp_array[i];
      if( !j)
         {
         memcpy( state_j, ival, 6 * sizeof( ldouble));
               /* subtract the analytic posn/vel from the numeric: */
         for( i = 0; i < n_vals; i++)
            ivals[0][i] = ival[i] - ref_state_j[i];
         }
      else
         for( i = 0; i < n_vals; i++)
            {
            ldouble tval = 0.;

            for( k = 0; k < j; k++)
               tval += bptr[k] * ivals_p[k][i];
            ivals[j][i] = tval * step + ivals[0][i];
            state_j[i] = ivals[j][i] + ref_state_j[i];
            }
      bptr += j;
      if( j != 6)
         {
#ifndef __WATCOMC__
         assert( fabsl( jd_j) < 1e+9);
#endif
         calc_derivativesl( jd_j, state_j, ivals_p[j], ref_orbit->central_obj);
         for( k = 0; k < 6; k++)
            ivals_p[j][k] -= ref_state_j[k + 3];
         }
      else     /* on last iteration,  we have our answer: */
         memcpy( ovals, state_j, n_vals * sizeof( ldouble));
      }

   for( i = 0; i < n_vals; i++)
      {
      ldouble tval = 0.;
      static const ldouble err_coeffs[6] = {
            RKF_CHAT1 - RKF_C1, RKF_CHAT2 - RKF_C2, RKF_CHAT3 - RKF_C3,
            RKF_CHAT4 - RKF_C4, RKF_CHAT5 - RKF_C5, RKF_CHAT6 - RKF_C6 };

      for( k = 0; k < 6; k++)
         tval += err_coeffs[k] * ivals_p[k][i];
      rval += tval * tval;
      }

   if( n_vals > 6)
      free( ivals[0]);
   return( sqrtl( rval * step * step));
}

double take_rk_step( const double jd, ELEMENTS *ref_orbit,
                 const double *ival, double *ovals,
                 const int n_vals, const double step)
{
   ldouble ivall[6], ovalsl[6], rvall;
   unsigned i;

   for( i = 0; i < 6; i++)
      ivall[i] = (ldouble)ival[i];
   rvall = take_rk_stepl( (long double)jd, ref_orbit, ivall,
               ovalsl, n_vals, (long double)step);
   for( i = 0; i < 6; i++)
      ovals[i] = (double)ovalsl[i];
   return( (double)rvall);
}

int symplectic_6( double jd, ELEMENTS *ref_orbit, double *vect,
                                          const double dt)
{
   int i, j;
#ifdef FOR_REFERENCE_ONLY
         /* Some compilers object to mathematically defined consts,  so  */
         /* I had to replace these lines with explicit numerical consts: */
   const double w1 = -0.117767998417887E1;
   const double w2 = 0.235573213359357E0;
   const double w3 = 0.784513610477560E0;
   const double w0 = (1-2*(w1+w2+w3));
   const double d6[7] = { w3, w2, w1, w0, w1, w2, w3 };
   const double c6[8] = { w3/2, (w3+w2)/2, (w2+w1)/2, (w1+w0)/2,
                         (w1+w0)/2, (w2+w1)/2, (w3+w2)/2, w3/2 };
#endif
   static const double d6[7] = { 0.7845136104775600,  0.2355732133593570,
            -1.1776799841788700, 1.3151863206839060, -1.1776799841788700,
             0.2355732133593570, 0.7845136104775600 };
   static const double c6[8] = { 0.3922568052387800,  0.5100434119184585,
            -0.4710533854097566, 0.0687531682525180,  0.0687531682525180,
            -0.4710533854097566, 0.5100434119184585,  0.3922568052387800 };

   for( i = 0; i < 8; i++)
      {
      double deriv[6];
      const double step = dt * c6[i];

      for( j = 0; j < 3; j++)
         vect[j] += step * vect[j + 3];
      jd += step;
      if( i != 7)
         {
         assert( fabs( jd) < 1e+9);
         calc_derivatives( jd, vect, deriv, ref_orbit->central_obj);
         for( j = 3; j < 6; j++)
            vect[j] += dt * d6[i] * deriv[j];
         }
      }
   return( 0);
}