source: proiecte/NBody/NBI/NBIserial.c @ 124

Last change on this file since 124 was 124, checked in by (none), 14 years ago

added NBI

File size: 58.7 KB
Line 
1
2
3#define CODENAME "NBI"
4#define VERSION    "1"
5
6/* ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
7
8                        NBI
9  A set of numerical integrators for the gravitational N-body problem.
10
11        (C) Copyright 1996 Ferenc Varadi.  All rights reserved.
12  This copyright notice is intended to prevent the commercial use of this code.
13  The code is free and can be freely distributed.
14  Do not remove this copyright notice unless you make changes in the code.
15  If you make changes in the code, add your name to the copyright notice.
16
17  The usual disclaimers apply, i.e., use the code at your own risk.
18  The author does not guarantee anything.
19
20  This code was developed by F. Varadi in collaboration
21  with M. Ghil, W. I. Newman, W. M. Kaula, K. Grazier, D. Goldstein
22  and M. Lessnick. The partial financial support of NSF Grant
23  ATM95-23787 (M. Ghil and F. Varadi) is gratefully acknowledged.
24
25  If you have any questions, suggestions etc., contact F. Varadi,
26  e-mail: varadi@ucla.edu
27  WWW: http://www.atmos.ucla.edu/~varadi
28
29 
30  For the impatient:
31        Save this file as NBI.c
32        Compile the code as
33        cc -o NBI NBI.c -lm
34
35        Download the sample input file NBIsampleinput and save it as
36        NBIsampleinput.  It is an input file for the outer Solar System, with
37        the planets from Jupiter through Neptune based on the Jet Propulsion
38        Laboratory's DE403 (Digital Ephemeris 403).  This was obtained
39        from the ftp site navigator.jpl.nasa.gov, updates can be accessed
40        there.
41
42        Run the code:
43        NBI NBIsampleinput
44        Take a look at the output in NBIsampleoutput and the log file
45        NBIsampleoutput.log Change the MAXNBODIES macro, i.e., the
46        maximum number of particles, to suit your needs.  It is set to
47        2100 in this version.
48
491) Description:
50        Four numerical integrators for the gravitational (nonrelativistic)
51        N-body problem are bundled into a single source code. The N bodies
52        are assumed to be point masses. Particles with zero mass are called
53        test particles; they are assumed not to influence each other and
54        the particles with nonzero mass.
55
56        The integrators are:
57
58        a) A Wisdom-Holman-type mapping for hierarchical N-body systems.
59                Code name: HWH, maximum global order = 4
60           This can deal with more general arrangements than the standard
61           Wisdom-Holman mapping.  Hierarchical N-body systems are described
62           by Roy (1988). These are represented in the code by binary trees.
63           The integrator takes advantage of certain properties of
64           symplectic forms and Jacobi coordinates, these will be described
65           in a paper (in preparation). We also use singularly weighted
66           symplectic forms, these are discussed by Varadi et al. (1995).
67
68        b) A modified Cowell-Stormer integrator, with modifications by
69           W. I. Newman and his students.
70                Code name: CSN, max global order = 15
71           Note: (global order) = (local order) - 2 (Goldstein, 1996)
72
73        c) A Gragg-Bulirsch-Stoer integrator, as it is described by
74           Hairer et al., (1993).
75                Code name: GBS, max order = 20 (max stage = 10)
76        IMPORTANT: The user has to specify the number of stages
77        in the input file and not the order. Order = 2*(number of stages),
78        see the Notes below.
79
80        d) A symplectic mapping with Kinetic-Potential energy splitting
81                Code name: SKP, max global order = 4
82           This is included for the sake of completeness.
83
84       
852) References:
86       
87        Goldstein, D.: 1996, The Near-Optimality of Stormer
88        Methods for Long Time Integrations of y''=g(y),
89        Ph.D. Dissertation, Univ. of California, Los Angeles,
90        Dept. of Mathematics
91
92        Hairer, E., Norsett, S. P. and Wanner, G.: 1993,
93        Solving Ordinary Differential Equations I.  Nonstiff Problems.
94        Second Revised Edition
95       
96        Lessnick, M.: 1996, Stability Analysis of Symplectic
97        Integration Schemes, Ph.D. Dissertation, Univ. of California,
98        Los Angeles, Dept. of Mathematics
99
100        Roy, A. E.: 1988, Orbital Motion, Institute of Physics
101        Publishing, Bristol
102       
103        Varadi, F. De la Barre, Kaula, W. M. and Ghil, M.: 1995,
104        Singularly Weighted Symplectic Forms and Applications to
105        Asteroid Motion, Celestial Mechanics and Dynamical Astronomy,
106        vol. 62, pp. 23-41
107
108        Yoshida, H.: 1993, Recent Progress in The Theory and
109        Application of Symplectic Integrators,
110        Celestial Mechanics and Dynamical Astronomy, vol. 56, pp. 27-43
111       
112       
113
1143) Notes
115       
116        All long-term integrations suffer from roundoff errors.
117        Symplectic schemes guarantee symplecticity for short
118        integrations but they are not strictly symplectic
119        for very long integrations.  The implementation of
120        the Wisdom-Holman mapping in this code was aimed at
121        reducing the effects of roundoff errors and thus make
122        the scheme as symplectic as possible.
123       
124        The Cowell-Stormer integrator, being a multistep
125        integrator, has to be intialized.  In theory, the starting
126        points should be computed with machine precision. This
127        can be done using quadrupole precision arithmetic but
128        the Gragg-Bulirsch-Stoer scheme in double precision
129        appers to be quite adequate. It should be noted that
130        the stability of the Cowell-Stormer scheme ensures
131        that initial errors do not grow due to computational
132        limitations (Goldstein, 1996).
133
134        In the case of the Gragg-Bulirsch-Stoer scheme one
135        specifies the number of stages k in the extrapolation
136        scheme. The order of the method is then 2*k.
137        This factor of 2 comes from the fact that the
138        underlying method (Gragg's) has only even powers
139        of the step size in the expansion of the local error
140        (cf. Hairer et al., 1993, especially equation 9.22).
141
142        All integrators use only G*mass, there is no need for
143        the masses themselves nor the value of G.  Any set of units
144        can be used:
145        distance_unit for position,
146        time_unit for time and step size,
147        distance_unit/time_unit for velocity.
148        E.g., one can use AU, day and AU/day as units.
149
150        None of the codes is regularized in any sense. We intend
151        to provide regularized code in the future. For the Wisdom-Holman
152        mapping this involves changing the tree of Jacobian reference
153        frames and dealing with the transition region accurately.
154
155        In the case of the Wisdom-Holman mapping, the computation of
156        acceleration differences between Jacobi and inertial coordinates
157        could be improved in order to further reduce roundoff error.
158        These differences can be represented by multi-pole expansions
159        of the potential but we have not explored the details yet.
160
161        The present version of the integrator does not monitor error.
162
163        The choice of step size is up to the user. In most cases
164        one has to adjust the step size to the shortest orbital
165        period T in the system. The Wisdom-Holman mapping at second
166        order is usually used with step size of T/10-T/100.
167        We recommend T/1000 for Cowell-Stormer integrator
168        since this renders the calculation EXACT to double
169        precision, even for high (e=0.5) eccentricities.
170        The cost of using this short time step is substantially
171        compensated by reduced computational complexity (as
172        compared to the Wisdom-Holman mapping).  For much larger
173        step sizes the integrator is not stable.
174
175        The Gragg-Bulirsch-Stoer method appears to be stable
176        for large orders and step sizes. For order 20 one can
177        even use T/5 as the step size and still get accurate results.
178        Since this is an accurate but compuationally very expensive
179        integrator, it is mainly used by us for short integrations.
180        We have not investigated the propagation of roundoff
181        error in the GBS method.
182       
183
1844) How to compile:
185        Since there is only a single source file, it is very easy to
186        compile the code with whatever C compiler is available. The
187        code was compiled and tested on Sun SPARCstations (Solaris 2.5.1
188        and earlier) and DEC Alpha workstations (Digital UNIX V3.2).
189        For instance, on the SPARCstations,
190       
191        cc -fast -xO5 -xdepend -o NBI NBI.c -lm
192
193        should work with Sun's C compiler.  The code uses only the standard
194        C library, more exactly: I/O functions, exit, malloc and free for
195        the tree, sqrt, cbrt, sin, cos, sinh, cosh.
196
197        Caution: optimization may introduce problems with some compilers
198        and it is prudent to check the correctness of the optimized code
199        by compiling without optimization, e.g.,
200        cc -g -o NBI_no_op NBI.c -lm
201        By comparing the output from NBI and NBI_no_op one can detect if
202        the optimization leads to any problems.  (It should not but in
203        practice the situation is not so clear.)
204
205        One also can play with various compiler switches which control
206        roundoff but we have not explored these.
207
2085) How to run
209
210        NBI inputfile
211
212        (Everything has to be specified in inputfile.)
213
2146) On hierarchical N-body systems
215
216        We can give only a brief introduction into these; a detailed technical
217        description is in preparation.  Our approach is somewhat different
218        from that of Roy (1988).
219
220        One starts with the concept of joining subsystems and applies this
221        concept recursively to form new systems. When two subsystems are
222        joined, the two subsystems are assumed to behave as point masses
223        and thus define the appropriate two-body problem for the Wisdom-Holman
224        mapping. The subsystems are, in general, not point masses and one
225        takes this into account in the perturbation step of the Wisdom-Holman
226        mapping.
227
228        Instead of a simple chain of Jacobi coordinates, the code uses a
229        tree of Jacobi coordinates. This makes it possible to use the code
230        e.g., in situations when the massive particles have sattelites, either
231        massive ones or not. It also reduces, when carefully specified, the
232        integration error since more of the perturbations are absorbed into
233        two-body interactions. The tree has to be specified in the input file. 
234        One can define this tree formally as it is processed by the code which
235        works like a simple finite automaton but we did not use lex, yacc or
236        any other compiler-compiler tools.
237
238 Tree description rules:
239 Rule 1) particles are indexed by consecutive natural numbers, from 1 to N,
240 Rule 2) particles are subsystems,
241 Rule 3) (x y) : join the subsystems x and y to create a new subsystem.,
242 Rule 4) [i1 - i2] denotes a set of test particles, from index i1 to index i2,
243        which are treated the same way,
244        Note: the space before and after - is necessary!
245 Rule 5) ([i1 - i2] i3) is forbidden, use instead (i3 [i1 - i2])
246
247 Example: 2000 main belt asteroids, from  6 to 2006, and
248        Sun, Jupiter, Saturn, Uranus and Neptune (1 - 5)
249        can described as
250 (((((1 [6 - 2006]) 2) 3) 4) 5)
251
252 Another example for the tree specification when test particles between
253 Jupiter and Saturn are added with indeces from 6 to 1000 and another set
254 of test particles is added between Saturn and Uranus with indices from
255 1001 to 2000. The indices for the massive bodies are:
256 Sun 1, Jupiter, 2, Saturn 3, Uranus 4, Neptune 5.
257((((((1 2) [6 - 1000]) 3) [1001 - 2000]) 4) 5)
258
259 Yet another example: the Sun, Earth, Moon, Jupiter and Io with indeces
260 1,2,3,4 and 5, respectively, for which the tree is
261 ((1 (2 3)) (4 5))
262
263
2647) Format of the input file
265Id method order step_size  total_integration_time printing_time
266        number_of_particles
267
268 tree_description
269
270GMs of particles, possibly zeros
271
272Initial conditions: x,y,z, vx, vy, vz
273
274Explanation:
275
276        Id is an arbitrary string which specifies the name of the experiment.
277
278        Method can be SKP, HWH, CSN and GBS.
279
280        Tree_description is used only when the method is HWH or the HWH is used
281        to initialize CSN. One can put a dummy tree description of the form
282        (1 2)
283        in the input file when HWH is not used.  This will satisfy the parser
284        which processes the tree description but will not play any role in the
285        integration.
286
287
288Example:
289 For Sun, Jupiter, Saturn, Uranus and Neptune,
290 with indeces 1,2,3,4 and 5, in heliocentric coordinates.
291
292SomeId HWH 4 100.0 365.0e+4 1.0e+4
2935
294((((1 2) 3) 4) 5)
295
2960.295912208285591095e-03   
2970.282534590952422643e-06  0.845971518568065874e-07 
2980.129202491678196939e-07   0.152435890078427628e-07 
299
3000.0 0.0 0.0 0.0 0.0 0.0
301
302-0.538420940678015203e+01   -0.831247035699103631e+00 -0.225095848657298592e+00 
3030.109236311953937086e-02   -0.652329334256263223e-02 -0.282301443780311710e-02 
304
3050.788989169798583756e+01   0.459570785756998301e+01  0.155843220515231762e+01
306-0.321720261683698921e-02   0.433063255278440598e-02  0.192641782724487106e-02
307
308-0.182698980946940850e+02   -0.116271406901560681e+01 -0.250371938307287545e+00 
3090.221540447608191936e-03   -0.376765389967669917e-02 -0.165324449144407023e-02
310
311-0.160595376822070470e+02   -0.239429495784517528e+02 -0.940042953888094424e+01 
3120.264312302715172384e-02   -0.150349219713705076e-02 -0.681271071793930938e-03
313
314
315There are 2 output files:
316        1) SomeId contains the actual integration data
317        2) SomeId.log contains a copy of the input data and error messages.
318        Some error messages might be sent to the standard output and the user
319        has to redirect it in order to capture them.
320
321In the case of the Cowell-Stormer integrator, an additional file, SomeId.lerr
322is created. This contains the LOCAL truncation error of the integration,
323computed for each step but only printed for the last step in the printing
324cycle. This also should help to determine the appropriate step size.
325
326The output is in the original coordinate system, i.e., inertial and not Jacobian.
327Output is produced when time is an integer multiple of the specified printing
328time. The following is printed:
329time    relative_change_in_total_energy relative_change_in_ang.mom.z.comp.
330after which
331x y z
332dx dy dz
333for each particle.
334
335
3369) The uset can change the code, of course. The section MAIN LOOP is where
337   this can be done easily. You can insert code to do more things, e.g.,
338   detecting close approaches.
339
34010) MACROS TO GIVE SOME CONTROL OVER THINGS.
341Changes in the code that you can make HERE:
342        CHANGE the DIM macro to 2 for planar motions.
343        (Note: a variables dim is also used, easy to change
344                it everywhere to DIM. It is still better to keep
345                the Kepler solver in the more general 3-dimensional form.)
346        CHANGE the MAXNBODIES macro to whatever you want,
347                including test particles.
348
349**************** END OF DESCRIPTION ********************** */
350
351/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++
352        USER-CHANGEABLE MACROS    */
353
354/* dimension of physical space */
355#define DIM                     3
356/* max number of bodies, obviously */
357#define MAXNBODIES              2100
358
359/* TOLERANCES */
360
361/* the extrapolation scheme stops when
362        Tj,k - T(j-1),k becomes smaller than this */
363#define GBSTOL                  5.0e-32
364
365/* The solution of Kepler's equation is accepted when
366   when iterates change by an amount smaller than this  */
367#define KEPLERSOLVERTOLERANCE   5.0e-16
368
369/* what method and order to use to initialize the CSN integrator, GBS or HWH  */
370#define CSNINITGBS              1
371#define CSNINITGBSORDER         4       
372#define CSNINITHWH              2
373#define CSNINITHWHORDER         4
374int csninit = CSNINITGBS;
375
376
377/* which extrapolation sequence to use in the GBS method */
378#define EXSEQHARM int extrapolseq[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10};
379#define EXSEQBUL  int extrapolseq[] = { 0, 1, 2, 3, 4, 6, 8, 12, 16, 24, 32};
380
381/* cange this line accordingly */
382EXSEQBUL
383
384
385/* VVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVVV */
386
387/*      The actual code.
388        Change things below this line at your own peril.
389       
390 */
391
392/* MACROS TO DISTINGUISH BETWEEN TEST PARTICLES AND
393   PARTICLES WITH FINITE MASS
394   You can add more types here, these values are stored
395        in tpstore for each particle
396   e.g., you could define stopped test particles,
397   but you have to change the code accordingly.
398   It is perhaps easier to add more integer arrays.
399   IMPORTANT: THE TREE MAINTAINS ITS OWN INTEGER ARRAYS
400   WHICH HOLD TEST PARTICLE INDICES.
401*/
402
403#define NOTTESTPARTICLE         1
404#define TESTPARTICLE            0
405
406/* the header files we need, only standard stuff */
407#include <stdio.h>
408#include <stdlib.h>
409#include <string.h>
410#include <math.h>
411
412/* Macro definitions */
413#define NUMMETHODS      4
414#define SKP 1
415#define HWH 2
416#define CSN 3
417#define GBS 4
418
419/* it is really the number of stages for GBS */
420char *methodnames[] = { "   ", "SKP", "HWH", "CSN", "GBS" };
421#define MAXCSN  15
422#define MAXGBS  10
423int MAXORDERS[] = { 0, 4, 4, MAXCSN, MAXGBS };
424
425/*
426   MACROS WHICH DEFINE WHAT IS STORED FOR THE PARTICLES
427        from 1 to 6: x y z vx vy vz
428                  7: G*mass
429        from 8 to 10:  accelerations in x, y and z
430   e.g., particle i's y velocity is
431         parts[i][5]
432*/
433
434#define TMPSIZE      (2*DIM+1)
435/* Note: we need only the G*mass for the actual integrations */ 
436#define GMASSCO      (7)
437#define ACC          (7)
438/* you can add more data for particles here
439   but keep NPARTDATA to specify the total number of
440   coordinates and other components  */
441#define NPARTDATA    (ACC+3)
442/* this is used from index 1 in both indeces */
443#define DECLARR [MAXNBODIES+1][NPARTDATA+1]
444
445/* Cowell-Stormer-Newman array */
446#define DECLSTORMER [MAXNBODIES+1][DIM+1][MAXCSN+2]
447/* Gragg-Bulirsch-Stoer array */
448#define DECLBS [MAXNBODIES+1][2*DIM+1][2][MAXGBS+2]
449
450#define DIM2         (2*DIM)
451
452/* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
453                GLOBAL DATA
454        ARRAYS THAT HOLD THE PARTICLES AND AUXILIARY DATA
455*/
456
457/* ARRAYS THAT HOLD THE PARTICLES' DATA.
458        This can be moved into main() without
459        any adverse effect.  */
460double parts DECLARR;
461
462/* This holds test particle indices.
463   It seems easier to keep this global. It is used, e.g.,
464   by PlainAccel to avoid the computation of force between
465   test particles.  */
466int tpstore[MAXNBODIES+1];
467
468/* this is needed for error messages from deep within the integrator */
469FILE *FLOG;
470
471/* Local error monitoring for the Cowell-Stormer integrator */
472double CSNerror = 0;
473
474
475/* FUNCTION AND STRUCT DECLARATIONS */
476
477/* the tree to store Jacobi coordinates */
478struct tree {
479        double tmass; /* total G*mass */
480        double rm1;   /* total mass / mass1 */
481        double rm2;   /* total mass / mass2 */ 
482        double bc[DIM2+1]; /* barycenter: position and velocity
483                                for actual particle these  are pos. and vel. */
484        double acc[DIM+1];
485/* nonzero particle index indicates that this is one with finite mass
486 or it is an array of test particles which can be processed in a simple loop,
487   zero means that this is a node in the tree
488*/
489        int particleindex;  /* nonzero indicates actual particle(s) */
490        struct tree *p1;
491        struct tree *p2;
492        int *testparticles; /* nonzero indicates array of test particles */
493        };
494
495
496int CStep(int ord, double h, int nump, double parts DECLARR, 
497                struct tree *treein);
498
499int GBSStep(int ord, double h, int nump, double parts DECLARR);
500
501/* makes a simple chain, usual planetary motion case */
502struct tree *MakeSimpleTree(int nump);
503/* makes an arbitrary tree from a textual description, i.e.,
504  parses the text and builds the tree */
505struct tree *MakeTree(FILE *f);
506
507/* Fill tree with barycenters, used to initialize the tree */
508int FillCMinTree(struct tree *tree, double p DECLARR);
509/* copy particle coordinates from tree;
510   since the tree is the primary storage for the coordinates
511   of particles with finite mass, one needs to call this
512   before accessing data in the "parts" array
513 */
514int ParticlesFromTree(struct tree *tree, double p DECLARR);
515
516/* The actual HWH integrator */
517int HWHStep(int order, struct tree *tree, double stepsize, double p DECLARR);
518
519/* Kinetic-potential energy split */
520int SKPStep(int order, double stepsize, int nump, double p DECLARR);
521
522
523/* access integer arrays in trees that hold the indeces
524        of test particles attached to that point,
525        returns pointer to testparticles,
526        included for the sake of completeness, not actually used */
527int **GetTestParticles(struct tree *tree, int ind);
528
529
530/*    VERY SIMPLE IO, you can write your own if you need something more */
531int readdata(FILE *f, int num, double p DECLARR) {
532        int i, j; double tmp;
533        for ( i = 1; i <= num; ++i ) {
534                        fscanf(f, "%le", &tmp); p[i][GMASSCO] = tmp; }
535        for ( i = 1; i <= num; ++i ) {
536                for ( j = 1; j <= 6; ++j ) {
537                        fscanf(f, "%le", &tmp); p[i][j] = tmp; } }
538        return num; }
539
540int printdataA(int num, FILE *f, double p DECLARR) {
541        int i, j;
542        for ( i = 1; i <= num; ++i ) {
543                for ( j = 1; j <= DIM; ++j ) { fprintf(f, "%.14le  ", p[i][j]); } 
544                fprintf(f, "%c", '\n');
545                for ( j = DIM+1; j <= 2*DIM; ++j ) { fprintf(f, "%.14le  ", p[i][j]); } 
546                fprintf(f, "%c", '\n');
547                }
548        return num; }
549
550int printdata(int num, FILE *f, double p DECLARR, int c1, int c2) {
551        int i, j;
552        for ( i = 1; i <= num; ++i ) {
553                for ( j = c1; j <= c2; ++j ) {
554                        fprintf(f, "%.14le  ", p[i][j]); } 
555                fprintf(f, "%c", '\n');
556                }
557        return num; }
558
559/* total energy for particles with finite mass */
560double energy(int nump, int dim, double p DECLARR, int *tps) {
561        double h = 0, tmp, tmp2;
562        int i,j,k;
563        for ( i = 1; i <= nump; ++i ) {
564                if ( tps[i] == 0 ) continue;
565                tmp = 0.0;
566                for ( k = 1; k <= dim; ++k ) {
567                        tmp = tmp + p[i][k+dim]*p[i][k+dim]; }
568                tmp = tmp*p[i][GMASSCO]/2.0;
569                h = h + tmp;
570                for ( j = i+1; j <= nump; ++j ) {
571                        tmp = 0.0;
572                        for ( k = 1; k <= dim; ++k ) {
573                                tmp2 = p[i][k]-p[j][k];
574                                tmp = tmp+tmp2*tmp2; }
575                        tmp = 1/sqrt(tmp);
576                        h = h - tmp*p[i][GMASSCO]*p[j][GMASSCO]; }
577                }
578        return h; }
579
580/* z component of angular momentum */
581double angmomz(int nump, int dim, double p DECLARR, int *tps) {
582        double az = 0.0;
583        int i;
584        for ( i = 1; i <= nump; ++i ) {
585                if ( tps[i] == 0 ) continue;
586                az += (p[i][1]*p[i][2+dim]-p[i][2]*p[i][1+dim])*p[i][GMASSCO]; }
587        return az; }
588
589/* transforms all particles to coordinates relative to barycenter */
590int tobarycenter(int nump, int dim, double p DECLARR) {
591        double bc[10];
592        int i, j;
593        for (j = 0; j <= 2*dim; ++j ) bc[j] = 0.0;
594        for ( i = nump; i >= 1; --i ) {
595                bc[0] += p[i][GMASSCO];
596                for (j = 1; j <= 2*dim; ++j ) bc[j] += p[i][j]*p[i][GMASSCO]; }
597        for (j = 1; j <= 2*dim; ++j ) bc[j] = bc[j]/bc[0];
598        for ( i = nump; i >= 1; --i ) {
599                for (j = 1; j <= 2*dim; ++j ) p[i][j] -= bc[j]; }
600        return 0; }
601
602
603
604
605/* """"""""""""""""""""""""""""""""""""""""""""""""""""""""""
606        MAIN
607*/
608
609int main(int argc, char **argv) {
610FILE *f, *fin, *flog, *flocerr;
611char sicode[100];
612int icode; /* SKP = 1 HWH = 2 CSN = 3 GBS = 4 */
613struct tree *tree; int num, num2;
614int ii, i, order = 1, nin, nps;
615double stepsize;
616double time = 0.0, ftime, ptime;
617char fnamin[1000]; char id[1000];
618double az, az0, toten, toten0;
619
620/* get command-line argument */
621if ( argc < 2 ) { exit(1); }
622sscanf(argv[1], "%s", fnamin); 
623
624/* read input file */
625fin = fopen(fnamin, "r");
626fscanf(fin, "%s", id); strcpy(fnamin, id);
627/* create output and log file */
628f = fopen(id, "w"); strcat(fnamin, ".log"); flog = fopen(fnamin, "w"); 
629FLOG = flog;
630
631#define FAILEDIN { fclose(fin); fclose(f); fclose(flog); exit(1); }
632
633fprintf(flog, "ID = %s\n", id);
634fscanf(fin, "%s", sicode); 
635fscanf(fin, "%d", &order); fprintf(flog, "ORDER = %2d\n", order);
636icode = 0;
637for ( icode = 1; icode <= NUMMETHODS; ++icode ) {
638        if ( strcmp(sicode, methodnames[icode] ) == 0 ) { 
639                if ( order < 1 || order > MAXORDERS[icode] ) {
640fprintf(flog, " Order %2d for method %s is not implemented \n", order, sicode);
641        FAILEDIN }
642                break; } }
643if ( icode > NUMMETHODS ) {
644fprintf(flog, " Invalid method code %s\n", sicode);
645        FAILEDIN }
646
647
648fscanf(fin, "%le", &stepsize); fprintf(flog, "STEPSIZE = %.12f\n", stepsize);
649fscanf(fin, "%le", &ftime); fprintf(flog, "FINAL TIME = %.12f\n", ftime);
650fscanf(fin, "%le", &ptime); fprintf(flog, "PRINTING TIME = %.12f\n", ptime);
651fscanf(fin, "%d", &num); fprintf(flog, "NUMBER OF PARTS = %2d\n", num);
652if ( num > MAXNBODIES ) {
653fprintf(flog, " Too many particles, change MAXNBODIES macro in code ");
654        FAILEDIN }
655
656
657/* initialize tpstore */
658for ( i = 1; i <= num+2; ++i ) { tpstore[i] = NOTTESTPARTICLE; }
659
660/* this will TESTPARTICLE for test particles in tpstore */
661tree = MakeTree(fin);
662
663/* For CSN create error file */
664if ( icode == CSN ) {
665strcpy(fnamin, id); strcat(fnamin, ".lerr"); flocerr = fopen(fnamin, "w"); }
666
667num2 = readdata(fin, num, parts);
668/* DONE READING ALL INPUTS */
669fclose(fin);
670
671/* zero out test particles' tpstore, again, just in case */
672for ( i = 1; i <= num; ++i ) { 
673        if ( parts[i][GMASSCO] == 0.0 ) { tpstore[i] = TESTPARTICLE; } }
674
675if ( num != num2 ) exit(1);
676fprintf(flog, "GMs\n");
677 printdata(num, flog, parts, GMASSCO, GMASSCO); 
678fprintf(flog, "Inits \n");
679         printdataA(num, flog, parts);
680fflush(flog);
681
682/* WATCH OUT FOR THIS
683   The transformation to barycenter
684   tends to improve things.
685*/
686tobarycenter(num, DIM, parts);
687fprintf(flog, "Actual Barycentric Inits \n"); printdataA(num, flog, parts);
688
689fprintf(flog, "Actual integrator is %s \n", sicode); 
690fprintf(flog, "CODE: %s \n", CODENAME); 
691fprintf(flog, "VERSION: %s \n", VERSION); 
692fflush(flog);
693
694
695/* z component of total angular momentum,
696  it is used to monitor accumulation of roundoff */
697az0 = angmomz(num, DIM, parts, tpstore);
698toten0 = energy(num, DIM, parts, tpstore);
699
700
701
702/* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
703        ACTUAL INTEGRATION
704*/
705
706nin = (int) (fabs(ptime)/fabs(stepsize));  /* the length of the inner loop, i.e.,
707                                        between printing */
708if ( nin == 0 ) {
709        fprintf(flog, "%s\n", 
710        " Stepsize is larger than printing time, aborted");
711        exit(1); }
712nps = (int) (fabs(ftime)/fabs(ptime));  /* the lenght of the outer loop,
713                        i.e., printings from beginning to end */
714
715/* WARNING: VERY LONG INTEGRATIONS CAN PRODUCE VERY LARGE
716        INTEGERS AND ONE HAS TO MAKE SURE THAT THEY ARE STILL
717        SMALLER THAN THE LARGEST REPRESENTABLE INTEGER
718  ANOTHER WARNING: ROUNDOFF EFFECTS MAY CAUSE nin AND/OR nps TO BE
719  DIFFERENT FROM WHAT THE USER EXPECTS
720*/
721
722
723/* print the first data point, i.e., the initial data */ 
724#define PRINTDATA fprintf(f, "%.12f   ", time);\
725        az = angmomz(num, DIM, parts, tpstore); \
726        toten = energy(num, DIM, parts, tpstore); \
727        fprintf(f, "%.16e  ", (toten-toten0)/toten0); \
728        fprintf(f, "%.16e\n", (az-az0)/az0); \
729         printdataA(num, f, parts);
730
731#define PRINTERROR \
732        if ( icode == CSN ) { \
733        fprintf(flocerr, "%.12f    %.12e\n", time, CSNerror); }
734
735PRINTDATA
736PRINTERROR
737
738
739/* MAIN LOOP */
740for ( ii = 1; ii <= nps; ++ii ) {
741for ( i = 1; i <= nin; ++i ) {
742switch(icode) {
743        case SKP:
744          SKPStep(order, stepsize, num, parts); break;
745        case HWH:
746          HWHStep(order, tree, stepsize, parts); break;
747        case CSN:
748          CStep(order, stepsize, num, parts, tree); break;
749        case GBS:
750          GBSStep(order, stepsize, num, parts); break;
751        default: 
752                fprintf(flog, "Not implemented\n");
753                fclose(flog); fclose(f); exit(1); break;
754        } 
755/* DONE ONE TIME STEP
756        ADD MORE CODE HERE IF NECESSARY
757        variables you can use:
758        double time;
759        double parts DECLARR;
760        Note that actual time should be computed as
761        actualtime = time+stepsize*i
762        */
763
764}
765
766 time += stepsize*nin;
767
768 if ( icode == HWH ) { ParticlesFromTree(tree, parts); }
769
770/* CHANGE THIS PART OR ADD YOUR CODE IF NECESSARY */
771       
772        PRINTDATA
773        PRINTERROR
774
775 if ( stepsize > 0.0 && time > ftime ) break;
776 if ( stepsize < 0.0 && time < ftime ) break;
777
778        }
779
780
781fclose(f);
782fprintf(flog, "\n Finished \n");
783fclose(flog);
784if ( icode == CSN ) fclose(flocerr);
785       
786return 0; }
787
788
789/* >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
790        FOR TEST: THE SIMPLE SPLIT, I.E., Kinetic - Potential
791        It does not provide special treatment for test particles.
792*/
793int PlainSplitKineticStep(int num, double stepsize, double p DECLARR) {
794        int i, k;
795        for ( i = 1; i <= num; ++i ) {
796                for ( k = 1; k <= DIM; ++k ) {
797                        p[i][k] += stepsize*p[i][k+DIM]; } }
798        return 0; }
799
800int PlainSplitPotentialStep(int num, double stepsize, double p DECLARR) {
801        int i, j, k;
802        double tmp, tmp2, tmp3;
803        for ( i = 1; i <= num; ++i ) {
804                tmp = 0.0;
805                for ( j = i+1; j <= num; ++j ) {
806                        tmp2 = 0.0;
807                        for ( k = 1; k <= DIM; ++k ) {
808                                tmp3 = p[i][k] - p[j][k];
809                                tmp2 = tmp2 + tmp3*tmp3; }
810                        tmp3 = 1.0/sqrt(tmp2);
811                        tmp2 = tmp3/tmp2;
812                        for ( k = 1; k <= DIM; ++k ) {
813                                tmp3 = p[i][k] - p[j][k];
814                                tmp = -tmp3*tmp2*p[j][GMASSCO]*stepsize;
815                                p[i][k+DIM] += tmp;
816                                tmp = tmp3*tmp2*p[i][GMASSCO]*stepsize;
817                                p[j][k+DIM] += tmp; }
818                        }
819        }
820        return 0; }
821
822
823
824/* >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
825        COWELL-STORMER INTEGRATOR */
826
827
828/* this computes accelerations, used both for CSN and GBS,
829        not to be confused with Accel used for HWH */
830
831int PlainAccel(int num, double p DECLARR) {
832        int i, j, k; int tpi, tpj;
833        double tmp, tmp2, tmp3;
834        double massi, massj;
835        for ( i = 1; i <= num; ++i ) {
836                for ( k = 1; k <= DIM; ++k ) { p[i][ACC+k] = 0.0; } }
837/* for testing things
838for ( i = 1; i <= num; ++i ) {
839        for ( k = 1; k <= DIM; ++k ) { p[i][ACC+k] = 1.0; } }
840return 0;
841*/
842        for ( i = num; i >= 1; --i ) { 
843                tpi = tpstore[i];
844                massi = p[i][GMASSCO];
845                for ( j = i+1; j <= num; ++j ) {
846                        tpj = tpstore[j];
847                        if ( tpi == TESTPARTICLE && tpj == TESTPARTICLE ) continue;
848                        tmp2 = 0.0;
849                        for ( k = 1; k <= DIM; ++k ) {
850                                tmp3 = p[i][k] - p[j][k];
851                                tmp2 = tmp2 + tmp3*tmp3; }
852                        tmp3 = 1.0/sqrt(tmp2);
853                        tmp2 = tmp3/tmp2;
854                        massj = p[j][GMASSCO];
855                        for ( k = 1; k <= DIM; ++k ) {
856                                tmp3 = p[i][k] - p[j][k];
857                                tmp = tmp3*tmp2;
858                        if ( tpj == NOTTESTPARTICLE ) {
859                                p[i][ACC+k] -= (tmp*massj); }
860                        if ( tpi == NOTTESTPARTICLE ) {
861                                p[j][ACC+k] += (tmp*massi); }
862                                }
863                        }
864        }
865        return 0; }
866
867/* Cowell-Stormer-Newman 14(6)th order */
868
869
870/* Newman's 14(6)th order Stormer's coeffs */
871
872double CSalpha[] =  {   0.0,
873                        1.000000000000000000000,
874                        0.000000000000000000000,
875                        0.083333333333333333333,
876                        0.083333333333333333333,
877                        0.079166666666666666666,
878                        0.075000000000000000000,
879                        0.071345899470899470899,
880                        0.068204365079365079365,
881                        0.065495756172839506172,
882                        0.063140432098765432098,
883                        0.061072649861712361712,
884                        0.059240564123376623376,
885                        0.057603625837453135733,
886                        0.056129980884507174190,
887                        0.054794379107071147139
888                                                 };
889 
890 
891double CSbeta[] =  {    0.0,
892                        1.50000000000000000000,
893                        0.33333333333333333333,
894                        0.37500000000000000000,
895                        0.35277777777777777777,
896                        0.33402777777777777777,
897                        0.31924603174603174603,
898                        0.30736607142857142857,
899                        0.29757660934744268077,
900                        0.28933077050264550264,
901                        0.28225737868098979210,
902                        0.27609762576993479771,
903                        0.27066578505957226195,
904                        0.26582499331955247133,
905                        0.26147199790503706399,
906                        0.25752728193649391773
907                                               };
908 
909       
910
911double  sumfrombottom(double *a, int nr, double *coef) {
912        int i;
913        double sum = 0.0;
914        for ( i = nr; i >= 1;  --i ) {
915                sum += a[i]*coef[i]; }
916        return sum; }
917
918int differencing(double *a, int nr, double in) {
919        int i, j, p;
920        double tmp1, tmp2;
921        tmp1 = a[1]; a[1] = in;
922        for ( i = 2; i <= nr; ++i ) {
923                tmp2 = a[i]; a[i] = a[i-1]-tmp1; tmp1 = tmp2; } 
924        return 0; }
925
926
927
928static int CStepCount = 0;
929static int CStepInit = 0;
930static struct tree *initCStree;
931
932int CStep(int ord, double h, int nump, 
933                double parts DECLARR, struct tree *treein) {
934        double *a, asav, vx, sum, tmp; int p,j, i;
935static double df DECLSTORMER;
936/* actual local order to determine how terms to sum */
937if ( CStepInit == 0 ) {
938        initCStree = treein;
939        CStepInit = 1;
940        for ( i = 1; i <= ord+1; ++i ) {
941        for ( p = 1; p <= nump; ++p ) {
942        for ( j = 1; j <= DIM; ++j ) {
943                df[p][j][i] = 0.0; } } }        }
944if ( CStepInit == 1 ) {
945/* initialization stage */
946                for ( p = 1; p <= nump; ++p ) {
947                for ( j = 1; j <= DIM; ++j ) {
948                        df[p][j][0] = -parts[p][j]; } }
949
950                if ( csninit == CSNINITHWH ) {
951                        HWHStep(CSNINITHWHORDER, initCStree, h, parts); 
952                        ParticlesFromTree(initCStree, parts);  }
953                else { GBSStep(CSNINITGBSORDER, h, nump, parts); }
954                if ( CStepCount > ord ) CStepInit = 2;
955                for ( p = 1; p <= nump; ++p ) {
956                for ( j = 1; j <= DIM; ++j ) {
957                        df[p][j][0] += parts[p][j]; } }
958                }
959else {
960/* actual run */
961CSNerror = 0;
962for ( p = 1; p <= nump; ++p ) {
963for ( j = 1; j <= DIM; ++j ) {
964/* formula: x(n+1) - x(n) = (x(n) - x(n-1)) + h^2*sum
965        with X(n) = x(n)-x(n-1)
966        X(n+1) = X(n) + h^2*sum
967        x(n+1) = X(N+1) + x(n)
968X(n) is in a[0]
969*/
970        a = df[p][j];
971        asav = a[0];
972        sum = sumfrombottom(a, ord, CSalpha);
973        sum *= h*h;
974        a[0] += sum;
975        parts[p][j] += a[0];
976/* vx = (1.0/h)*(x(n)-x(n-1)) + h*sum */
977        sum = h*sumfrombottom(a, ord, CSbeta);
978        vx = asav/h + sum;
979        parts[p][DIM+j] = vx;
980/* record error */
981        tmp = fabs(CSalpha[ord]*a[ord]*h*h);
982        if ( CSNerror < tmp ) CSNerror = tmp;
983        } } }
984PlainAccel(nump, parts);
985for ( p = 1; p <= nump; ++p ) {
986for ( j = 1; j <= DIM; ++j ) {
987        differencing(df[p][j], ord, parts[p][ACC+j]);
988        parts[p][ACC+j] = 0.0;
989        } }
990++CStepCount;
991return 0; }
992
993
994
995/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
996        BULIRSCH-STOER
997        actually, GRAGG-BULIRSCH-STOER, following Hairer et al., 1993
998*/
999
1000/* Gragg's symmetric integrator, this is used to generate Tj,1  */
1001int GStep(int nump, double parts DECLARR, double s DECLARR, double h, int len) {
1002int k, j, p, d, knext, kprev, m;
1003double dif, difsum, tmp1, astep;
1004static double y0 DECLARR;
1005static double y1 DECLARR;
1006/* save originial in y0 */
1007for ( p = 1; p <= nump; ++p ) {
1008                y0[p][GMASSCO] = parts[p][GMASSCO];
1009for ( d = 1; d <= DIM; ++d ) { 
1010                y0[p][d] = parts[p][d]; 
1011                y0[p][d+DIM] = parts[p][d+DIM]; } }
1012/* consult Hairer et al. why we need this
1013        in short: we need a method which has only even powers
1014        of h in the expansion of its error, i.e., a symmetric method
1015 */
1016h = h/2.0;
1017PlainAccel(nump, y0);
1018for ( p = 1; p <= nump; ++p ) {
1019/* load masses into y1 */
1020                y1[p][GMASSCO] = y0[p][GMASSCO];
1021for ( d = 1; d <= DIM; ++d ) { 
1022                y1[p][d] = y0[p][d]+h*y0[p][DIM+d]; 
1023                y1[p][d+DIM] = y0[p][d+DIM]+h*y0[p][ACC+d]; } }
1024/* compute y2 (m=1), y3' (m=2), y4 (m=3), y5' (m=4) y(2n+1)' (m=2*n)
1025        n = len
1026        */
1027for ( m = 1; m <= 2*len; ++m ) {
1028        if ( (m%2) ) {
1029        for ( p = 1; p <= nump; ++p ) {
1030        for ( d = 1; d <= DIM; ++d ) { 
1031                y1[p][d] = y0[p][d]+2.0*h*y1[p][d+DIM]; 
1032                s[p][d] = y1[p][d]; y0[p][d] = y1[p][d]; 
1033                }  } }
1034        else {
1035        PlainAccel(nump, y1);
1036        for ( p = 1; p <= nump; ++p ) {
1037        for ( d = 1; d <= DIM; ++d ) { 
1038                y0[p][d+DIM] = y1[p][d+DIM]; 
1039                y1[p][d+DIM] = y0[p][d+DIM]+2.0*h*y1[p][ACC+d]; 
1040                s[p][d+DIM] = (y0[p][d+DIM]+y1[p][d+DIM])/2.0; 
1041                } } 
1042                }
1043        }
1044return 0; }
1045
1046
1047int GBSStep(int bsord, double h, int nump, double parts DECLARR) {
1048int k, j, jfinal, p, d, jnext, jprev, len;
1049static double gbsarr DECLBS;
1050double dif, difsum, tmp1, astep;
1051int *n = extrapolseq;
1052static double s DECLARR;
1053/* for testing GStep
1054        astep = h; len = 1;
1055        GStep(nump, parts, s, astep, len);
1056for ( p = 1; p <= nump; ++p ) {
1057for ( d = 1; d <= 2*DIM; ++d ) { parts[p][d] = s[p][d]; } }
1058return 0;
1059*/
1060jfinal = bsord;
1061for ( j = 1; j <= bsord; ++j ) {
1062        difsum = 0.0;
1063/* map j and j-1 to storage index */
1064        jnext = ((j)%2); 
1065        jprev = ((j-1)%2);
1066/* compute Tj,1 */
1067        len = n[j];
1068        astep = h/len;
1069        GStep(nump, parts, s, astep, len);
1070        for ( p = 1; p <= nump; ++p ) {
1071        for ( d = 1; d <= 2*DIM; ++d ) { gbsarr[p][d][jnext][1] = s[p][d]; } }
1072/* compute Tj,k */
1073for ( k = 1; k <= j-1; ++k ) {
1074/* Aitken-Neville for nonsymmetric underlying method */
1075        tmp1 = ((double) n[j])/((double) n[j-k]); 
1076/* Aitken-Neville for symmetric underlying method (square of previous) */
1077        tmp1 = tmp1*tmp1;
1078for ( p = 1; p <= nump; ++p ) {
1079for ( d = 1; d <= 2*DIM; ++d ) {
1080        dif = (gbsarr[p][d][jnext][k] - gbsarr[p][d][jprev][k]);
1081        if ( k == j-1 ) difsum += fabs(dif);
1082        gbsarr[p][d][jnext][k+1] = gbsarr[p][d][jnext][k] + dif/(tmp1-1.0);
1083        } }
1084        } 
1085/* stop computing more is reached roundoff level */
1086        if ( j > 1 && difsum <= GBSTOL ) { jfinal = j; break; }
1087        }
1088j = jfinal; k = j; jnext = ((j)%2); 
1089for ( p = 1; p <= nump; ++p ) {
1090for ( d = 1; d <= 2*DIM; ++d ) { parts[p][d] = gbsarr[p][d][jnext][k]; } } 
1091return 0; }
1092       
1093
1094
1095
1096
1097
1098/* >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
1099        Actual  HWH integration */
1100int MoveCMandKepler(struct tree *tree, double *cmove, double stepsize, 
1101                double p DECLARR);
1102int KeplerOnTree(struct tree *tree, double stepsize, double p DECLARR);
1103int RecursivePerturbations(struct tree *tree, double stepsize, double p DECLARR);
1104int RecPertAcc(struct tree *tree, double stepsize, double p DECLARR);
1105
1106/* fourth order coeffs */
1107static double Coeff4x0 = 0.0;
1108static double Coeff4x1 = 0.0;
1109
1110static int SympCalled = 0;
1111
1112void initSympCoefs() {
1113        double cbrt2 = pow(2.0, 1.0/3); //cbrt(2.0);
1114        Coeff4x0 = -cbrt2/(2.0-cbrt2);
1115        Coeff4x1 = 1.0/(2.0-cbrt2);
1116        SympCalled = 1;
1117}
1118       
1119
1120static int HWHCalled = 0;
1121
1122int HWHStep(int order, struct tree *tree, double stepsize, double p DECLARR) {
1123        if ( SympCalled == 0 ) { initSympCoefs(); }
1124        if ( HWHCalled == 0) {
1125/* initialize tree */
1126                FillCMinTree(tree, p); 
1127                HWHCalled = 1; }
1128        if ( order == 1 ) {
1129                KeplerOnTree(tree, stepsize, p);
1130                RecursivePerturbations(tree, stepsize, p);
1131                return 0; }
1132        if ( order == 2 ) {
1133                KeplerOnTree(tree, stepsize/2.0, p);
1134                RecursivePerturbations(tree, stepsize, p);
1135                KeplerOnTree(tree, stepsize/2.0, p); 
1136                return 0; }
1137        if ( order == 4 ) {
1138                HWHStep(2, tree, stepsize*Coeff4x1, p);
1139                HWHStep(2, tree, stepsize*Coeff4x0, p);
1140                HWHStep(2, tree, stepsize*Coeff4x1, p);
1141                return 0; }
1142        printf(" Order %2d  is not implemented, Aborted", order);
1143        exit(1); }
1144
1145
1146int SKPStep(int order, double stepsize, int nump, double p DECLARR) {
1147        if ( SympCalled == 0 ) { initSympCoefs(); }
1148        if ( order == 1 ) {
1149                PlainSplitKineticStep(nump, stepsize, p);
1150                PlainSplitPotentialStep(nump, stepsize, p);
1151                return 0; }
1152        if ( order == 2 ) {
1153                PlainSplitKineticStep(nump, stepsize/2.0, p);
1154                PlainSplitPotentialStep(nump, stepsize, p);
1155                PlainSplitKineticStep(nump, stepsize/2.0, p);
1156                return 0; }
1157        if ( order == 4 ) {
1158                SKPStep(2, stepsize*Coeff4x1, nump, p);
1159                SKPStep(2, stepsize*Coeff4x0, nump, p);
1160                SKPStep(2, stepsize*Coeff4x1, nump, p);
1161                return 0; }
1162        printf(" Order %2d  is not implemented, Aborted", order);
1163        exit(1); }
1164
1165
1166
1167/* +++++++++++++++++++++++++++++++++++++++++++++++++++++
1168/       TREE
1169/
1170/
1171*/
1172
1173
1174int inittree(struct tree *tmp) {
1175        int i; tmp->tmass = 0.0; tmp->particleindex = 0;
1176        tmp->p1 = 0; tmp->p2 = 0;
1177        tmp->testparticles = 0; tmp->rm1 = 0; tmp->rm2 = 0;
1178        for ( i = 1; i <= DIM2; ++i ) { tmp->bc[i] = 0; }
1179        for ( i = 1; i <= DIM; ++i ) { tmp->acc[i] = 0; }
1180        return 0; }
1181
1182
1183struct tree *newtree() {
1184        struct tree *tmp = (struct tree*) malloc(sizeof(struct tree)); 
1185        inittree(tmp); return tmp; }
1186
1187void deletetree(struct tree *p) { 
1188        if ( p->p1 ) deletetree(p->p1);
1189        if ( p->p2 ) deletetree(p->p2);
1190        free((void *) p); }
1191
1192
1193struct tree *MakeSimpleTree(int nump) {
1194        int i;
1195        struct tree *root;
1196        struct tree *t1, *p1, *p2;
1197        root = newtree();
1198        t1 = root;
1199        for ( i = nump; i > 1; --i ) {
1200                p1 = newtree();
1201                p2 = newtree();
1202                p1->particleindex = 0;
1203                p2->particleindex = i;
1204                t1->p1 = p1; t1->p2 = p2;
1205                t1 = p1; }
1206        t1->particleindex = 1;
1207        return root; }
1208
1209int eatwhite(FILE *f) {
1210        int c; 
1211        for ( ; ; ) {
1212                c = getc(f);
1213                if ( !isspace(c) ) break; }
1214        ungetc(c, f); 
1215        return c; }
1216
1217
1218/* just an array to hold indeces read in,
1219        portions of this are assigned to
1220        tpa's below */
1221int tmpReadList[MAXNBODIES+1000];
1222/* bookkeeping */
1223int freelist = 0;
1224
1225int *ReadList(FILE *f) {
1226        int i1, i2; int c; int i; int *arr;
1227        c = eatwhite(f);
1228        if ( c != '[' ) {
1229                printf(" Error: bad list, aborted \n"); exit(1); }
1230        c = getc(f);
1231        fscanf(f, "%d - %d", &i1, &i2);
1232        c = eatwhite(f);
1233        if ( c != ']' ) {
1234                printf(" Error: bad list, aborted \n"); exit(1); }
1235        c = getc(f); 
1236        arr = tmpReadList+freelist;
1237/* put zero in tpstore */
1238        for ( i = i1; i <= i2; ++i ) tpstore[i] = TESTPARTICLE;
1239        for ( i = i1; i <= i2; ++i ) arr[i-i1+1] = i;
1240        arr[0] = i2-i1+1;
1241        freelist += i2-i1+2;
1242        return arr; }
1243       
1244               
1245
1246struct tree *MakeTree(FILE *f) {
1247        int i; int c; int pind; int *tpa;
1248        struct tree *root;
1249        struct tree *p1, *p2;
1250        c = eatwhite(f);
1251        if ( c != '(' ) {
1252                if ( isdigit(c) ) {
1253/* single particle */
1254                        fscanf(f,"%d",&pind);
1255                        root = newtree();
1256                        root->particleindex = pind;
1257                        return root; } 
1258/* set of test particles */
1259                tpa = ReadList(f);
1260                root = newtree();
1261                root->particleindex = tpa[1];
1262                root->testparticles = tpa;
1263                return root; }
1264        c = getc(f);
1265        root = newtree();
1266        p1 = MakeTree(f);
1267        p2 = MakeTree(f);
1268        root->p1 = p1; root->p2 = p2;
1269        root->particleindex = 0;
1270        c = eatwhite(f);
1271        if ( c != ')' ) {
1272                printf(" Error: bad tree, aborted \n"); exit(1); }
1273        c = getc(f);
1274        return root; }
1275
1276
1277
1278/* not used now but might come handy */
1279struct tree *getparticleintree(struct tree *tree, int ind) {
1280        struct tree *tmp;
1281        if ( tree->particleindex ) { 
1282                if ( tree->particleindex == ind ) { return tree; }
1283                return 0; }
1284        tmp = getparticleintree(tree->p1, ind);
1285        if ( tmp != 0 ) return tmp;
1286        tmp = getparticleintree(tree->p2, ind);
1287        if ( tmp != 0 ) return tmp;
1288        return 0; }
1289       
1290
1291/* not used now but might come handy */
1292int **GetTestParticles(struct tree *tree, int ind) {
1293        struct tree *tmp = getparticleintree(tree, ind);
1294        if ( tmp == 0 ) return 0;
1295        return &(tmp->testparticles); }
1296       
1297
1298/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1299        ACTUAL NUMERICAL STUFF ON TREE
1300*/
1301
1302/* fills in initial data: masses, mass ratios, barycenters */
1303int FillCMinTree(struct tree *tree, double p DECLARR) {
1304        int j; double *bc = tree->bc;
1305        int i, ind, tp, num, *tpa = tree->testparticles;
1306        for ( j = 1; j <= DIM; ++j ) tree->acc[j] = 0.0;
1307        tree->tmass = 0.0; tree->rm1 = 0.0; tree->rm2 = 0.0;
1308        for ( j = 1; j <= DIM2; ++j ) bc[j] = 0.0;
1309        ind = tree->particleindex;
1310           if ( ind ) {
1311              if ( tpa != 0 ) {  num = tpa[0];
1312                for ( tp = 1; tp <= num; ++tp ) { 
1313                        ind = tpa[tp]; if ( ind < 1 ) continue;
1314                        for ( i = 1; i <= DIM; ++i ) p[ind][ACC+i] = 0.0; } } 
1315               else { tree->tmass = p[ind][GMASSCO];
1316                for ( j = 1; j <= 2*DIM; ++j ) bc[j] = p[ind][j]; }
1317           return 0; } 
1318        FillCMinTree(tree->p1, p);
1319        FillCMinTree(tree->p2, p);
1320        tree->tmass = tree->p1->tmass + tree->p2->tmass;
1321        tree->rm1 = (tree->p1->tmass)/(tree->tmass); 
1322        tree->rm2 = (tree->p2->tmass)/(tree->tmass);
1323        for ( j = 1; j <= DIM2; ++j ) {
1324                bc[j] = tree->p1->bc[j] * tree->rm1
1325                + tree->p2->bc[j] * tree->rm2; } 
1326        return 0; }
1327
1328/* copy from tree into p */
1329int ParticlesFromTree(struct tree *tree, double p DECLARR) {
1330        int j; double *bc = tree->bc;
1331        int i, ind, *tpa = tree->testparticles;
1332        ind = tree->particleindex;
1333        if ( ind ) {
1334           if ( tpa == 0 ) {
1335                tree->tmass = p[ind][GMASSCO];
1336                for ( j = 1; j <= DIM2; ++j ) p[ind][j] = bc[j]; }
1337           return 0; } 
1338        ParticlesFromTree(tree->p1, p);
1339        ParticlesFromTree(tree->p2, p);
1340        return 0; }
1341
1342
1343/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1344        KEPLER STEP ON TREE
1345*/
1346
1347/* computes Stumpff functions several ways, depending on x */
1348int StumpffFunctions(double *vals, double x);
1349/* advances the Kepler motion 1/2*v^2 - mu/|x| for time deltat */
1350double KeplerStumpffStep(int dim, double mu, double *xv, double deltat);
1351/* return only the increment in xv */
1352double IncrementKeplerStumpffStep(int dim, double mu, double *xv, double deltat);
1353
1354/* NOTE: this does not CHANGE centers of masses */
1355int KeplerOnTree(struct tree *tree, double stepsize, double p DECLARR) { 
1356/* move the center of mass */
1357        double cmove[10]; int i;
1358        for ( i = 1; i <= 2*DIM; ++i ) cmove[i] = 0.0;
1359/* move by 1/2*(Tmass)*v^2 */
1360        for ( i = 1; i <= DIM; ++i ) cmove[i] = (tree->bc[DIM+i])*stepsize; 
1361        MoveCMandKepler(tree, cmove, stepsize, p);
1362return 0; }
1363
1364
1365/* Assume that test particles are attached to 2nd */
1366int DoTestPartKepler(struct tree *tree, double *cmove, 
1367                        double stepsize, double p DECLARR) {
1368        int i, tp, num, *tpa, ind;
1369        double mu, tmpkep[TMPSIZE], *bc;
1370        tpa = tree->p2->testparticles; 
1371        num = tpa[0];
1372        if ( num > 0 ) {
1373                bc = tree->p1->bc; mu = tree->tmass;
1374        for ( tp = 1; tp <= num; ++tp ) { ind = tpa[tp]; if ( ind < 1 ) continue;
1375                for ( i = 1; i <= DIM2; ++i )   tmpkep[i] = p[ind][i] - bc[i];
1376                IncrementKeplerStumpffStep(DIM, mu, tmpkep, stepsize); 
1377                for ( i = 1; i <= DIM2; ++i ) p[ind][i] += (cmove[i]+tmpkep[i]);
1378                } }
1379        return 0; }
1380
1381/* move the center of mass and do a Kepler step on relative
1382  coordinates, transmit the new "moves" down */
1383int MoveCMandKepler(struct tree *tree, double *cmove, double stepsize, 
1384                double p DECLARR) {
1385        int i, *tpa, ind;
1386        double tmp[TMPSIZE], tmpkep[TMPSIZE], *xv1, *xv2, mu;
1387        double *bc = tree->bc;
1388        for ( i = 1; i <= DIM2; ++i ) { bc[i] += cmove[i]; }
1389/* do nothing else for actual particles */
1390        if ( tree->particleindex ) { return 0; }
1391        if ( tree->p2->testparticles ) {       
1392/* this is a pair where p2 is an array of test particles */
1393                DoTestPartKepler(tree, cmove, stepsize, p); 
1394                MoveCMandKepler(tree->p1, cmove, stepsize, p); 
1395                return 0; }
1396        mu = tree->tmass;
1397        xv1 = tree->p1->bc; xv2 = tree->p2->bc;
1398        for ( i = 1; i <= DIM2; ++i ) { tmpkep[i] = xv2[i]-xv1[i]; }
1399/* do Kepler step on tmpkep, we get back the increment */
1400        IncrementKeplerStumpffStep(DIM, mu, tmpkep, stepsize); 
1401/* use only the relative change */
1402        for ( i = 1; i <= DIM2; ++i ) { tmp[i] = cmove[i] - tmpkep[i]*(tree->rm2); }
1403        MoveCMandKepler(tree->p1, tmp, stepsize, p);
1404        for ( i = 1; i <= DIM2; ++i ) { tmp[i] = cmove[i] + tmpkep[i]*(tree->rm1); }
1405        MoveCMandKepler(tree->p2, tmp, stepsize, p);
1406        return 0; }
1407       
1408/*  ++++++++++++++++++++++++++++++++++++++++++++++++++++++
1409        PERTURBATION ACCELERATIONS
1410*/
1411
1412/* simple acceleration */
1413int Accel(double *x1, double *x2, double *acc) {
1414        double tmp[10], r = 0.0; int i; 
1415        for ( i = 1; i <= DIM; ++i ) {
1416                tmp[i] = x1[i]-x2[i];
1417                r = r + tmp[i]*tmp[i]; }
1418        r = sqrt(r);
1419        r = r*r*r;
1420        for ( i = 1; i <= DIM; ++i ) { acc[i] = -tmp[i]/r; }
1421        return 0; }
1422
1423int AdvanceTestParticleVel(struct tree *tree, double stepsize, double p DECLARR) {
1424                int i, ind, tp, num, *tpa = tree->testparticles;
1425                if ( tpa == 0 ) return 0; num = tpa[0];
1426        for ( tp = 1; tp <= num; ++tp ) {
1427                ind = tpa[tp]; if ( ind < 1 ) continue;
1428                for ( i = 1; i <= DIM; ++i ) {
1429                        p[ind][DIM+i] += stepsize*p[ind][ACC+i]; 
1430                        p[ind][ACC+i] = 0.0; }
1431                        } 
1432        return 0; }
1433                       
1434       
1435/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1436        PERTURBATIONS ON TREE
1437*/
1438
1439int RecursivePerturbations(struct tree *tree, double stepsize, double p DECLARR) {
1440        int j, ind; double *bc = tree->bc;
1441/* advance attached test particles */
1442        if ( tree->particleindex ) {
1443/* if particle we are ready to advance v */
1444                if ( tree->testparticles ) {
1445                        AdvanceTestParticleVel(tree, stepsize, p); 
1446                        return 0; }
1447/* this is a nonzero mass actual particle, advance v */
1448                for ( j = 1; j <= DIM; ++j ) {
1449                        bc[DIM+j] += stepsize*tree->acc[j]; tree->acc[j] = 0.0; }
1450                return 0; }
1451/* if not, first we take care of accs. for pairs with the two parts 
1452        IMPORTANT TO DO THIS FIRST
1453*/
1454        RecPertAcc(tree, stepsize, p);
1455/* now descend and do the same for both parts */
1456RecursivePerturbations(tree->p2, stepsize, p);
1457RecursivePerturbations(tree->p1, stepsize, p);
1458/* now regenerate center of mass etc. */
1459if ( tree->p2->testparticles ) {
1460        for ( j = 1; j <= DIM2; ++j ) { bc[j] = tree->p1->bc[j]; }
1461        return 0; }
1462for ( j = 1; j <= DIM2; ++j ) {
1463        bc[j] = tree->p1->bc[j] * tree->rm1
1464        + tree->p2->bc[j] * tree->rm2; }
1465return 0; }
1466
1467/* assume test particles in tree2 for dir > 1, reverse signs for dir < 0 */
1468int TestParticleAcc(int dir, double *cacc, struct tree *tree1, struct tree *tree2, 
1469                                double p DECLARR) {
1470                int i, ind, tp, num, *tpa = tree2->testparticles;
1471                double acc[TMPSIZE];
1472                double *bc = tree1->bc;
1473                double tmp, gm1 = tree1->tmass;
1474                num = tpa[0];
1475        for ( tp = 1; tp <= num; ++tp ) {
1476                ind = tpa[tp]; if ( ind < 1 ) continue;
1477/* rewrite this to avoid bad numerical attitude */
1478                        Accel(bc, &p[ind][0], acc);
1479                for ( i = 1; i <= DIM; ++i ) {
1480/* since test particles are in 2nd */
1481                        tmp = acc[i]-cacc[i]*dir;
1482                        p[ind][ACC+i] += -gm1*tmp; }
1483                        } 
1484        return 0; }
1485                       
1486
1487/* c1 and c2 are for possible future use */
1488int RecPertAccSub(struct tree *tree1, 
1489        struct tree *tree2, double *cacc, double *c1, double *c2, double p DECLARR) {
1490/* do the work if both are particles */
1491        int i, ind1, ind2; double acc[TMPSIZE];
1492        double tmp, gm1 = tree1->tmass, gm2 = tree2->tmass;
1493        ind1 = tree1->particleindex;
1494        ind2 = tree2->particleindex;
1495        if ( ind1 && ind2 ) {
1496/* do nothing if both are test particle arrays */
1497                if ( tree1->testparticles  && tree2->testparticles  ) return 0; 
1498/* test particles first */
1499                if ( tree2->testparticles ) {
1500                        TestParticleAcc(1, cacc, tree1, tree2, p); return 0; }
1501                if ( tree1->testparticles ) {
1502                        TestParticleAcc( -1, cacc, tree2, tree1, p); return 0; }
1503/* if actual particles with nonzero mass */
1504                Accel(tree1->bc, tree2->bc, acc);
1505                for ( i = 1; i <= DIM; ++i ) {
1506/* rewrite this to avoid bad numerical attitude */
1507                        tmp = acc[i]-cacc[i];
1508                        tree1->acc[i] += gm2*tmp; tree2->acc[i] += -gm1*tmp; }
1509        return 0; }
1510/* if not, descend but only for one branch, the other is taken
1511                care in turn */
1512        if ( !ind1 ) {
1513                RecPertAccSub(tree1->p2, tree2, cacc, c1, c2, p);
1514                RecPertAccSub(tree1->p1, tree2, cacc, c1, c2, p);
1515                return 0; }
1516        RecPertAccSub(tree1, tree2->p1, cacc, c1, c2, p);
1517        RecPertAccSub(tree1, tree2->p2, cacc, c1, c2, p);
1518        return 0; }
1519
1520/* assume test particles in tree2 */
1521int TestParticleAccSpec(double *cm, struct tree *tree1, 
1522                        struct tree *tree2, double p DECLARR) {
1523                int i, ind, tp, num, *tpa = tree2->testparticles;
1524                double acc[TMPSIZE], cacc[TMPSIZE];
1525                double *bc = tree1->bc;
1526                double tmp, gm1 = tree1->tmass;
1527                num = tpa[0];
1528        for ( tp = 1; tp <= num; ++tp ) {
1529                ind = tpa[tp]; if ( ind < 1 ) continue;
1530/* rewrite this to avoid bad numerical attitude */
1531                        Accel(bc, &p[ind][0], acc);
1532                        Accel(cm, &p[ind][0], cacc);
1533                for ( i = 1; i <= DIM; ++i ) {
1534/* since test particles are in 2nd */
1535                        tmp = acc[i]-cacc[i];
1536                        p[ind][ACC+i] += -gm1*tmp; 
1537                                }
1538                        } 
1539        return 0; }
1540                       
1541
1542/* c1 and c2 are for possible future use */
1543int RecPertAccSubTest(double *cm,  struct tree *tree1, struct tree *tree2, double p DECLARR) {
1544        int ind1;
1545        ind1 = tree1->particleindex;
1546        if ( ind1 ) {
1547/* do nothing if tree1 has only test particles */
1548                if ( tree1->testparticles ) return 0;
1549                TestParticleAccSpec(cm, tree1, tree2, p); 
1550                return 0; }
1551/* if not, descend for both branches */
1552        RecPertAccSubTest(cm, tree1->p1, tree2, p);
1553        RecPertAccSubTest(cm, tree1->p2, tree2, p);
1554        return 0; }
1555
1556
1557int RecPertAcc(struct tree *tree, double stepsize, double p DECLARR) {
1558        double cacc[TMPSIZE];
1559/* do nothing if all in p2 are attached to p1: pure Kepler motion */
1560        if ( tree->p1->particleindex &&  tree->p2->particleindex  )  return 0; 
1561        if ( tree->p2->testparticles == 0 ) {
1562/* actual nonzero masses in p2 */
1563        Accel(tree->p1->bc, tree->p2->bc, cacc);
1564        RecPertAccSub(tree->p1, tree->p2, cacc, tree->p1->bc, tree->p2->bc, p); 
1565        return 0;}
1566/* test particles in p2 */
1567        RecPertAccSubTest(tree->p1->bc, tree->p1, tree->p2, p); 
1568        return 0; }
1569
1570
1571
1572
1573/*
1574// ===============================================================================
1575//
1576//   KEPLER MOTION USING STUMPFF FUCNTIONS
1577//      FOLLOWING DANBY, 1992
1578//      (quartic Newtonian solver, Stumpff function definitions,
1579//       f and g functions)
1580//      see also: Stiefel and Sheifele, Linear and Regular
1581//      Celestial Mechanics
1582//
1583// ===============================================================================
1584
1585// Stumpff functions (c0, c1, c2, c3)
1586*/
1587
1588/* CHANGE these numbers to tune accuracy */
1589/* this one can be changed runtime */
1590static double StumpffQSTol = KEPLERSOLVERTOLERANCE;
1591
1592/* the limit for which trig and hyp. trig functions are used */
1593#define STUMPFFPLAINLIM (5.0)
1594/* reduction limit */
1595#define STUMPFFRED (1.0/32.0)
1596#define STUMPFFSHORTTOL (1.0e-6)
1597
1598/* this one is not used now */
1599static double StumpffSerTol = 1.0e-25;
1600
1601
1602static int StumpffErrorCount = 0;
1603static int StumpffErrorCountLimit = 100;
1604
1605CheckStumpffErrorCount() {
1606        ++StumpffErrorCount;
1607        if ( StumpffErrorCount > StumpffErrorCountLimit ) {
1608                fprintf(FLOG, "\n Too many errors, aborted \n");
1609                fflush(FLOG);
1610                exit(1); }
1611        return 0; }
1612
1613double util_fact(int k) {
1614        int i; double fact = 1;
1615        for ( i = 1; i <= k; ++i ) {
1616                fact *= i; }
1617        return fact; }
1618
1619/*
1620In case we need the coefficients, produced by Mathematica
1621static double c2_0 = 0.5;
1622static double c2_1 = -1.0/24;
1623static double c2_2 = 1.0/720;
1624static double c2_3 = -1.0/40320;
1625static double c2_4 = 1.0/3628800;
1626static double c2_5 = -1.0/479001600;
1627static double c2_6 = 1.0/87178291200;
1628static double c2_7 = -1.0/20922789888000;
1629static double c2_8 = 1.0/6402373705728000;
1630static double c2_9 = -1.0/2432902008176640000;                         
1631static double c2_10 = 1.0/1124000727777607680000;
1632
1633static double c3_0 = 1.0/6;
1634static double c3_1 = -1.0/120;
1635static double c3_2 = 1.0/5040;
1636static double c3_3 = -1.0/362880;
1637static double c3_4 = 1.0/39916800;
1638static double c3_5 = -1.0/6227020800;
1639static double c3_6 = 1.0/1307674368000;
1640static double c3_7 = -1.0/355687428096000;
1641static double c3_8 = 1.0/121645100408832000;
1642static double c3_9 = -1.0/51090942171709440000;
1643static double c3_10 = 1.0/25852016738884976640000;
1644*/
1645
1646
1647/* This is what is actually used, generated by factoring the coefficients */
1648double HornerStumpff2(double x) {
1649        return
1650(1.0 - x*(1 - x*(1 - x*(1 - x*(1 - x *(1 - x
1651/* up to degree 10, if you need it
1652*(1 - x*(1 - x*(1 - (1 - x/462)*x/380)/306)/240)        */
1653/182)/132)/90)/56)/30)/12)/2;
1654        }
1655
1656double HornerStumpff3(double x) {
1657        return
1658 (1.0 - x*(1 - x*(1 - x*(1 - x*(1 - x*(1 - x
1659/* up to degree 10, if you need it
1660*(1 - x*(1 - x*(1 - (1 - x/506)*x/420)/342)/272)        */
1661/210)/156)/110)/72)/42)/20)/6;
1662        }
1663
1664
1665/* first 4 terms, ie., up to x^3 */
1666double Stumpff2Short(double x) {
1667        return (1.0 - x*(1 - x*(1 - x/56)/30)/12)/2;
1668        }
1669
1670double Stumpff3Short(double x) {
1671        return (1.0 - x*(1 - x*(1 - x/72)/42)/20)/6;
1672        }
1673
1674
1675
1676
1677/* Original for reference */
1678double StumpffSeriesCk(int k, double x) {
1679        int i;
1680        double term, sum = 0.0;
1681        double terms[100], num;
1682        term = 1.0;
1683        for ( i = 1; i <= 50; ++i ) {
1684                num = (k+2*i-1)*(k+2*i);
1685                term = term*(-x/num);
1686                terms[i] = term;
1687                if ( fabs(term) < StumpffSerTol ) break; 
1688                }
1689        for ( ; i > 0; --i )  sum += terms[i];
1690        sum = (1.0+sum)/util_fact(k);
1691        return sum; }
1692
1693int StumpffFunctions(double *c, double x) {
1694        double x2, tmp;
1695        int i;
1696/* What this is supposed to be
1697c[0] = StumpffSeriesCk(0, x);
1698c[1] = StumpffSeriesCk(1, x);
1699c[2] = StumpffSeriesCk(2, x);
1700c[3] = StumpffSeriesCk(3, x);
1701return 0;
1702*/
1703/* FOR REFERENCE, THIS COULD BE USED FOR LARGE VALUES, MAYBE
1704        if ( x > STUMPFFPLAINLIM ) {
1705                x2 = sqrt(x);
1706                c[0] = cos(x2);
1707                c[1] = sin(x2)/x2;
1708                c[2] = (1.0-cos(x2))/x;
1709                c[3] = (x2-sin(x2))/(x*x2);
1710                return 0; }
1711        if ( x < -STUMPFFPLAINLIM ) {
1712                x2 = sqrt(-x);
1713                c[0] = cosh(x2);
1714                c[1] = sinh(x2)/x2;
1715                c[2] = (cosh(x2)-1.0)/x;
1716                c[3] = (sinh(x2)-x2)/(x*x2);
1717                return 0; }
1718*/
1719        if ( fabs(x) < STUMPFFSHORTTOL ) {
1720                c[2] = Stumpff2Short(x); c[3] = Stumpff3Short(x);
1721                c[1] = 1.0 - x*c[3]; c[0] = 1.0 - x*c[2];
1722                return 0; }
1723        if ( fabs(x) < STUMPFFRED ) {
1724                c[2] = HornerStumpff2(x); c[3] = HornerStumpff3(x); 
1725                c[1] = 1.0 - x*c[3]; c[0] = 1.0 - x*c[2];
1726                return 0; }
1727/* reduce value to small */
1728        x2 = x;
1729        for ( i = 1; i <= 20; ++i ) {
1730                x2 = x2/4.0;
1731                if ( fabs(x2) < STUMPFFRED ) break; 
1732                }
1733        if ( i > 20 ) {
1734                fprintf(FLOG, "%s"," StumpfReduction No Good \n"); 
1735                CheckStumpffErrorCount();
1736                fflush(FLOG);
1737                return 0; } 
1738/* original
1739        c[2] = StumpffSeriesCk(2, x2);
1740        c[3] = StumpffSeriesCk(3, x2);
1741*/
1742        c[2] = HornerStumpff2(x2);
1743        c[3] = HornerStumpff3(x2);
1744        c[1] = 1.0 - x2*c[3];
1745        c[0] = 1.0 - x2*c[2];
1746        for ( ; i > 0; --i ) {
1747                tmp = c[2];
1748                c[2] = (c[1]*c[1])/2.0; 
1749                c[3] = (tmp+c[0]*c[3])/4.0; 
1750                c[1] = c[0]*c[1];
1751                c[0] = (2.0*c[0]*c[0])-1.0;
1752                }
1753        return i; }
1754
1755
1756/* ===============================================================================
1757*/
1758
1759/* NOTE:
1760It seems that there might be a systematic error here, not sure yet.
1761The Newton solvers tend to end up systematically on a preferred side
1762of the solution.  E.g., for x^2: if starting from above the solution,
1763it ends at a slightly larger value than the solution.  Starting from
1764 below, in the first step it will get to the previous case.  Since
1765the solution cannot be computed to machine precision, we could have
1766a systematic error.  Since the orbit spends more time near apocenter
1767than near pericenter, this seems to lead an increase in energy and
1768decrease in angular momentum.  Perhaps this could be randomized with
1769an additional bisection at the end?
1770*/
1771
1772double QuarticStumpffSolverIterator(double mu, double alpha, double s, 
1773                double r0, double u, double deltat) {
1774        double c[5];
1775        double d0, d1, d2, d3, s2, ssav = s, tmp1, tmp2, tmp3;
1776        double delta1, delta2, delta3;
1777        int i;
1778for ( i = 1; i <= 20; ++i ) {
1779        s2 = s*s;
1780        StumpffFunctions(c, alpha*s2);
1781        tmp1 = s*c[1];
1782        tmp2 = s2*c[2];
1783        tmp3 = -r0*alpha+mu;
1784        d0 = u*tmp2+mu*s2*s*c[3];
1785        d0 += r0*tmp1;
1786        d0 -= deltat; 
1787        d1 = u*tmp1+mu*tmp2; 
1788        d1 += r0*c[0];
1789        d2 = u*c[0]+tmp3*tmp1; 
1790        d3 = tmp3*c[0]-u*alpha*c[1]; 
1791        delta1 = -d0/d1;
1792        delta2 = -d0/(d1+delta1*d2/2);
1793        delta3 = -d0/(d1+delta1*d2/2+delta2*delta2*d3/6);
1794        ssav= s;
1795        s = s + delta3;
1796        if ( fabs(delta3/s) < StumpffQSTol ) return s; }
1797        fprintf(FLOG, "%s"," QUARTIC STUMPFF SOLVER NOT CONVERGENT? \n"); 
1798        fprintf(FLOG, "    s = %.16e \n", s); 
1799        fprintf(FLOG, "sprev = %.16e \n", ssav); 
1800        fprintf(FLOG, "delta3/s = %.16e \n", delta3/s); 
1801        fprintf(FLOG, "Stumpff x = %.16e \n", alpha*s2); 
1802        fprintf(FLOG, "deltas = %.16e ", delta1); fprintf(FLOG, "%.16e ", delta2); 
1803        fprintf(FLOG, "%.16e \n", delta3);
1804        CheckStumpffErrorCount();
1805        fflush(FLOG);
1806        return s; }
1807
1808
1809double QuarticStumpffSolver(double mu, double alpha, double s, 
1810                double r0, double u, double deltat) {
1811/* initial guess */
1812        double tmp = deltat/r0;
1813        double s0 = tmp;
1814        double r0dot = u/r0;
1815        tmp = tmp*tmp;
1816        s0 -= 0.5*r0dot*tmp;
1817/* third order in deltat */
1818        tmp = tmp*deltat/r0;
1819        s0 += (0.5*r0dot*r0dot-1/6*(mu-alpha*r0)/r0)*tmp;
1820        s = QuarticStumpffSolverIterator( mu, alpha, s0, r0, u, deltat);
1821        return s; }
1822
1823/* ===============================================================================
1824*/
1825
1826void IncfgAndDerivsStumpff(double *vals, double mu, double alpha, double s, 
1827                double r0, double u, double deltat) {
1828        double r, s2, tmp1, tmp2, tmp3;
1829        double c[5];
1830        s = QuarticStumpffSolver( mu, alpha, s, r0, u, deltat);
1831        s2 = s*s;
1832        StumpffFunctions(c, alpha*s2);
1833/* recompute c[1] from angular mom. */
1834        tmp1 = s*c[1];
1835        tmp2 = s2*c[2];
1836        tmp3 = u*tmp1+mu*tmp2;
1837/* r */
1838        r = tmp3+(r0*c[0]);
1839        vals[0] = r;
1840        tmp3 = -mu*tmp2;
1841/* f */
1842        vals[1] = tmp3/r0;
1843/* g, this is usually bad in terms of roundoff */
1844        vals[2] = r0*tmp1+u*tmp2;
1845/* derivs */
1846        vals[3] = -(mu*tmp1)/(r0*r);
1847        vals[4] =  tmp3/r;
1848        return; }
1849
1850void fgAndDerivsStumpff(double *vals, double mu, double alpha, double s, 
1851                double r0, double u, double deltat) {
1852        IncfgAndDerivsStumpff(vals, mu, alpha, s, r0, u, deltat);
1853        vals[1] += 1.0;
1854        vals[4] += 1.0;
1855        return; }
1856
1857double KeplerStumpffStep(int dim, double mu, double *xv, double deltat) {
1858        double fgvals[5], tmp;
1859        double *x = xv;
1860        double *v = xv+dim;
1861        double v02 = 0.0, u = 0.0, r0 = 0.0;
1862        double r0dot, alpha;
1863        int i;
1864        for ( i = 1; i <= dim; ++i ) {
1865                v02 += v[i]*v[i];
1866                u += x[i]*v[i];
1867                r0 += x[i]*x[i]; }
1868        r0 = sqrt(r0);
1869        r0dot = u/r0;
1870/* a = 1.0/(2.0/r0 - v02/mu); alpha = mu/a; */
1871        alpha = 2.0*mu/r0 - v02;
1872        fgAndDerivsStumpff(fgvals, mu, alpha, 0.0, r0, u, deltat);
1873/* now update */
1874        for ( i = 1; i <= dim; ++i ) {
1875                tmp = x[i];
1876                x[i] = fgvals[1]*x[i]+fgvals[2]*v[i];
1877                v[i] = fgvals[3]*tmp+fgvals[4]*v[i]; }
1878        return 0.0; }
1879
1880
1881/* return only the increment */
1882double IncrementKeplerStumpffStep(int dim, double mu, double *xv, double deltat) {
1883        double fg[5], tmp, tmp2;
1884        double *x = xv;
1885        double *v = xv+dim;
1886        double alpha, v02 = 0.0, u = 0.0, r02 = 0.0, r0;
1887        int i, ii;
1888        double tmp1, r, delta1;
1889        for ( i = 1; i <= dim; ++i ) {
1890                v02 += v[i]*v[i];
1891                u += x[i]*v[i];
1892                r02 += x[i]*x[i]; }
1893        r0 = sqrt(r02);
1894/* a = 1.0/(2.0/r0 - v02/mu); alpha = mu/a; */
1895        alpha = 2.0*mu/r0 - v02;
1896        IncfgAndDerivsStumpff(fg, mu, alpha, 0.0, r0, u, deltat);
1897
1898/* using fg'-gf' = 1 ( ' = dot) or total energy do not help much */
1899               
1900/* now update */
1901        for ( i = 1; i <= dim; ++i ) {
1902                tmp = x[i];
1903                x[i] = fg[1]*x[i]+fg[2]*v[i]; 
1904                v[i] = fg[3]*tmp+fg[4]*v[i]; 
1905                }
1906        return 0.0; }
1907
1908/* FOR TEST
1909double IncrementKeplerStumpffStep(int dim, double mu, double *xv, double deltat) {
1910        double tmp[7], tmp2[7];
1911        int i;
1912        for ( i = 1; i <= 2*dim; ++i ) { tmp[i] = tmp2[i] = xv[i]; }
1913IncrementKeplerStumpffStepSub(dim, mu, tmp, deltat);
1914        for ( i = 1; i <= 2*dim; ++i ) tmp2[i] += tmp[i];
1915IncrementKeplerStumpffStepSub(dim, mu, tmp2, -deltat);
1916        for ( i = 1; i <= 2*dim; ++i ) xv[i] = (tmp[i]-tmp2[i])/2.0;
1917        return 0.0; }
1918*/
1919
1920       
1921
1922
Note: See TracBrowser for help on using the repository browser.