source: proiecte/NBody/NBI/NBIparalel.c @ 162

Last change on this file since 162 was 162, checked in by (none), 14 years ago
File size: 46.4 KB
Line 
1
2
3#define CODENAME "NBI"
4#define VERSION    "1"
5
6/* dimension of physical space */
7#define DIM                     3
8/* max number of bodies, obviously */
9#define MAXNBODIES              2100
10
11/* TOLERANCES */
12
13/* the extrapolation scheme stops when
14        Tj,k - T(j-1),k becomes smaller than this */
15#define GBSTOL                  5.0e-32
16
17/* The solution of Kepler's equation is accepted when
18   when iterates change by an amount smaller than this  */
19#define KEPLERSOLVERTOLERANCE   5.0e-16
20
21/* what method and order to use to initialize the CSN integrator, GBS or HWH  */
22#define CSNINITGBS              1
23#define CSNINITGBSORDER         4       
24#define CSNINITHWH              2
25#define CSNINITHWHORDER         4
26int csninit = CSNINITGBS;
27
28
29/* which extrapolation sequence to use in the GBS method */
30#define EXSEQHARM int extrapolseq[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10};
31#define EXSEQBUL  int extrapolseq[] = { 0, 1, 2, 3, 4, 6, 8, 12, 16, 24, 32};
32
33/* cange this line accordingly */
34EXSEQBUL
35
36
37/* MACROS TO DISTINGUISH BETWEEN TEST PARTICLES AND
38   PARTICLES WITH FINITE MASS
39*/
40
41#define NOTTESTPARTICLE         1
42#define TESTPARTICLE            0
43
44/* the header files we need, only standard stuff */
45#include <stdio.h>
46#include <stdlib.h>
47#include <string.h>
48#include <math.h>
49#include "omp.h"
50
51/* Macro definitions */
52#define NUMMETHODS      4
53#define SKP 1
54#define HWH 2
55#define CSN 3
56#define GBS 4
57
58/* it is really the number of stages for GBS */
59char *methodnames[] = { "   ", "SKP", "HWH", "CSN", "GBS" };
60#define MAXCSN  15
61#define MAXGBS  10
62int MAXORDERS[] = { 0, 4, 4, MAXCSN, MAXGBS };
63
64/*
65   MACROS WHICH DEFINE WHAT IS STORED FOR THE PARTICLES
66        from 1 to 6: x y z vx vy vz
67                  7: G*mass
68        from 8 to 10:  accelerations in x, y and z
69   e.g., particle i's y velocity is
70         parts[i][5]
71*/
72
73#define TMPSIZE      (2*DIM+1)
74/* Note: we need only the G*mass for the actual integrations */ 
75#define GMASSCO      (7)
76#define ACC          (7)
77/* you can add more data for particles here
78   but keep NPARTDATA to specify the total number of
79   coordinates and other components  */
80#define NPARTDATA    (ACC+3)
81/* this is used from index 1 in both indeces */
82#define DECLARR [MAXNBODIES+1][NPARTDATA+1]
83
84/* Cowell-Stormer-Newman array */
85#define DECLSTORMER [MAXNBODIES+1][DIM+1][MAXCSN+2]
86/* Gragg-Bulirsch-Stoer array */
87#define DECLBS [MAXNBODIES+1][2*DIM+1][2][MAXGBS+2]
88
89#define DIM2         (2*DIM)
90
91/* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
92                GLOBAL DATA
93        ARRAYS THAT HOLD THE PARTICLES AND AUXILIARY DATA
94*/
95
96/* ARRAYS THAT HOLD THE PARTICLES' DATA.
97        This can be moved into main() without
98        any adverse effect.  */
99double parts DECLARR;
100
101/* This holds test particle indices.
102   It seems easier to keep this global. It is used, e.g.,
103   by PlainAccel to avoid the computation of force between
104   test particles.  */
105int tpstore[MAXNBODIES+1];
106
107/* this is needed for error messages from deep within the integrator */
108FILE *FLOG;
109
110/* Local error monitoring for the Cowell-Stormer integrator */
111double CSNerror = 0;
112
113
114/* FUNCTION AND STRUCT DECLARATIONS */
115
116/* the tree to store Jacobi coordinates */
117struct tree {
118        double tmass; /* total G*mass */
119        double rm1;   /* total mass / mass1 */
120        double rm2;   /* total mass / mass2 */ 
121        double bc[DIM2+1]; /* barycenter: position and velocity
122                                for actual particle these  are pos. and vel. */
123        double acc[DIM+1];
124/* nonzero particle index indicates that this is one with finite mass
125 or it is an array of test particles which can be processed in a simple loop,
126   zero means that this is a node in the tree
127*/
128        int particleindex;  /* nonzero indicates actual particle(s) */
129        struct tree *p1;
130        struct tree *p2;
131        int *testparticles; /* nonzero indicates array of test particles */
132        };
133
134
135int CStep(int ord, double h, int nump, double parts DECLARR, 
136                struct tree *treein);
137
138int GBSStep(int ord, double h, int nump, double parts DECLARR);
139
140/* makes a simple chain, usual planetary motion case */
141struct tree *MakeSimpleTree(int nump);
142/* makes an arbitrary tree from a textual description, i.e.,
143  parses the text and builds the tree */
144struct tree *MakeTree(FILE *f);
145
146/* Fill tree with barycenters, used to initialize the tree */
147int FillCMinTree(struct tree *tree, double p DECLARR);
148/* copy particle coordinates from tree;
149   since the tree is the primary storage for the coordinates
150   of particles with finite mass, one needs to call this
151   before accessing data in the "parts" array
152 */
153int ParticlesFromTree(struct tree *tree, double p DECLARR);
154
155/* The actual HWH integrator */
156int HWHStep(int order, struct tree *tree, double stepsize, double p DECLARR);
157
158/* Kinetic-potential energy split */
159int SKPStep(int order, double stepsize, int nump, double p DECLARR);
160
161
162/* access integer arrays in trees that hold the indeces
163        of test particles attached to that point,
164        returns pointer to testparticles,
165        included for the sake of completeness, not actually used */
166int **GetTestParticles(struct tree *tree, int ind);
167
168
169/*    VERY SIMPLE IO, you can write your own if you need something more */
170int readdata(FILE *f, int num, double p DECLARR) {
171        int i, j; double tmp;
172        for ( i = 1; i <= num; ++i ) {
173                        fscanf(f, "%le", &tmp); p[i][GMASSCO] = tmp; }
174        for ( i = 1; i <= num; ++i ) {
175                for ( j = 1; j <= 6; ++j ) {
176                        fscanf(f, "%le", &tmp); p[i][j] = tmp; } }
177        return num; }
178
179// X
180int printdataA(int num, FILE *f, double p DECLARR) {
181        int i, j;
182        for ( i = 1; i <= num; ++i ) {
183                for ( j = 1; j <= DIM; ++j ) { 
184                        fprintf(f, "%.14le  ", p[i][j]); 
185                } 
186                fprintf(f, "%c", '\n');
187                for ( j = DIM+1; j <= 2*DIM; ++j ) { 
188                        fprintf(f, "%.14le  ", p[i][j]); 
189                } 
190                fprintf(f, "%c", '\n');
191        }
192        return num; 
193}
194
195// X
196int printdata(int num, FILE *f, double p DECLARR, int c1, int c2) {
197        int i, j;
198        for ( i = 1; i <= num; ++i ) {
199                for ( j = c1; j <= c2; ++j ) {
200                        fprintf(f, "%.14le  ", p[i][j]); } 
201                    fprintf(f, "%c", '\n');
202                }
203        return num; 
204}
205
206//X total energy for particles with finite mass
207double energy(int nump, int dim, double p DECLARR, int *tps) {
208        double h = 0, tmp, tmp2;
209        int i,j,k;
210       
211        for ( i = 1; i <= nump; ++i ) {
212                if ( tps[i] == 0 ) 
213                        continue;
214                tmp = 0.0;
215               
216                #pragma omp parallel for reduction(+:tmp) schedule(runtime)
217                for ( k = 1; k <= dim; ++k ) {
218                        tmp = tmp + p[i][k+dim]*p[i][k+dim]; 
219                }
220                tmp = tmp*p[i][GMASSCO]/2.0;
221               
222                h = h + tmp;
223                for ( j = i+1; j <= nump; ++j ) {
224                        tmp = 0.0;
225                        #pragma omp parallel for reduction(+:tmp) schedule(runtime)
226                        for ( k = 1; k <= dim; ++k ) {
227                                tmp2 = p[i][k]-p[j][k];
228                                tmp = tmp+tmp2*tmp2; 
229                        }
230                        tmp = 1/sqrt(tmp);
231                        h = h - tmp*p[i][GMASSCO]*p[j][GMASSCO]; 
232                }
233        }
234       
235        return h; 
236}
237
238// X z component of angular momentum
239double angmomz(int nump, int dim, double p DECLARR, int *tps) {
240        double az = 0.0;
241        int i;
242       
243        for ( i = 1; i <= nump; ++i ) {
244                if ( tps[i] == 0 ) 
245                        continue;
246                // x*vy - y*vx 
247                az = az + (p[i][1]*p[i][2+dim] - p[i][2]*p[i][1+dim]) * p[i][GMASSCO]; 
248        }
249        return az; 
250}
251
252// X transforms all particles to coordinates relative to barycenter
253int tobarycenter(int nump, int dim, double p DECLARR) {
254        double bc[10];
255        int i, j;
256       
257        #pragma omp parallel private(j)
258        #pragma omp for schedule(runtime)
259        for (j = 0; j <= 2*dim; ++j ) 
260                bc[j] = 0.0;
261       
262        for ( i = nump; i >= 1; --i ) {
263                bc[0] += p[i][GMASSCO];
264                #pragma omp parallel private(j)
265                #pragma omp for schedule(runtime)
266                for (j = 1; j <= 2*dim; ++j ) 
267                        bc[j] += p[i][j]*p[i][GMASSCO]; 
268        }
269       
270        #pragma omp parallel private(j)
271        #pragma omp for schedule(runtime)
272        for (j = 1; j <= 2*dim; ++j ) 
273                bc[j] = bc[j]/bc[0];
274               
275   
276        for ( i = nump; i >= 1; --i ) {
277                #pragma omp parallel private(j)
278                #pragma omp for schedule(runtime)
279                for (j = 1; j <= 2*dim; ++j ) 
280                        p[i][j] -= bc[j]; 
281        }
282        return 0; 
283}
284
285
286
287
288/* """"""""""""""""""""""""""""""""""""""""""""""""""""""""""
289        MAIN
290*/
291
292int main(int argc, char **argv) {
293        FILE *f, *fin, *flog, *flocerr;
294        char sicode[100];
295        int icode; /* SKP = 1 HWH = 2 CSN = 3 GBS = 4 */
296        struct tree *tree; int num, num2;
297        int ii, i, order = 1, nin, nps;
298        double stepsize;
299        double time = 0.0, ftime, ptime;
300        char fnamin[1000]; char id[1000];
301        double az, az0, toten, toten0;
302
303        // get command-line argument
304        if ( argc < 2 ) { exit(1); }
305        sscanf(argv[1], "%s", fnamin); 
306
307        // read input file
308        fin = fopen(fnamin, "r");
309        fscanf(fin, "%s", id); strcpy(fnamin, id);
310       
311        // create output and log file
312        // id - name of output file     
313        f = fopen(id, "w"); strcat(fnamin, ".log"); flog = fopen(fnamin, "w"); 
314        FLOG = flog;
315
316        #define FAILEDIN { fclose(fin); fclose(f); fclose(flog); exit(1); }
317
318        fprintf(flog, "ID = %s\n", id);
319        fscanf(fin, "%s", sicode); 
320        fscanf(fin, "%d", &order); fprintf(flog, "ORDER = %2d\n", order);
321        icode = 0;
322        for ( icode = 1; icode <= NUMMETHODS; ++icode ) {
323                if ( strcmp(sicode, methodnames[icode] ) == 0 ) { 
324                        if ( order < 1 || order > MAXORDERS[icode] ) {
325                                fprintf(flog, " Order %2d for method %s is not implemented \n", order, sicode);
326                                FAILEDIN
327                        }
328                        break; 
329                } 
330        }
331       
332        if ( icode > NUMMETHODS ) {
333                fprintf(flog, " Invalid method code %s\n", sicode);
334                FAILEDIN
335        }
336
337
338        fscanf(fin, "%le", &stepsize); fprintf(flog, "STEPSIZE = %.12f\n", stepsize);
339        fscanf(fin, "%le", &ftime); fprintf(flog, "FINAL TIME = %.12f\n", ftime);
340        fscanf(fin, "%le", &ptime); fprintf(flog, "PRINTING TIME = %.12f\n", ptime);
341        fscanf(fin, "%d", &num); fprintf(flog, "NUMBER OF PARTS = %2d\n", num);
342       
343        if ( num > MAXNBODIES ) {
344                fprintf(flog, " Too many particles, change MAXNBODIES macro in code ");
345                FAILEDIN
346        }
347
348        /* initialize tpstore */
349       
350        for ( i = 1; i <= num+2; ++i ) { 
351                tpstore[i] = NOTTESTPARTICLE; 
352        }
353
354        /* this will TESTPARTICLE for test particles in tpstore */
355        tree = MakeTree(fin);
356
357        /* For CSN create error file */
358        if ( icode == CSN ) {
359                strcpy(fnamin, id); strcat(fnamin, ".lerr"); flocerr = fopen(fnamin, "w"); 
360        }
361
362        num2 = readdata(fin, num, parts);
363        /* DONE READING ALL INPUTS */
364        fclose(fin);
365
366        /* zero out test particles' tpstore, again, just in case */
367        for ( i = 1; i <= num; ++i ) { 
368                if ( parts[i][GMASSCO] == 0.0 ) { 
369                        tpstore[i] = TESTPARTICLE; 
370                } 
371        }
372
373        if ( num != num2 ) exit(1);
374        fprintf(flog, "GMs\n");
375        printdata(num, flog, parts, GMASSCO, GMASSCO); 
376       
377        fprintf(flog, "Inits \n");
378        printdataA(num, flog, parts);
379       
380        fflush(flog);
381
382        /* WATCH OUT FOR THIS
383           The transformation to barycenter
384           tends to improve things.
385        */
386        tobarycenter(num, DIM, parts);
387        fprintf(flog, "Actual Barycentric Inits \n"); 
388        printdataA(num, flog, parts);
389
390        fprintf(flog, "Actual integrator is %s \n", sicode); 
391        fprintf(flog, "CODE: %s \n", CODENAME); 
392        fprintf(flog, "VERSION: %s \n", VERSION); 
393        fflush(flog);
394
395
396        // z component of total angular momentum, it is used to monitor accumulation of roundoff
397        az0 = angmomz(num, DIM, parts, tpstore);
398        toten0 = energy(num, DIM, parts, tpstore);
399
400
401
402        /* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
403                ACTUAL INTEGRATION
404        */
405       
406        // the length of the inner loop, i.e.,  between printing
407        nin = (int) (fabs(ptime)/fabs(stepsize)); 
408       
409        if ( nin == 0 ) {
410                fprintf(flog, "%s\n", " Stepsize is larger than printing time, aborted");
411                exit(1); 
412        }
413       
414        // the length of the outer loop, i.e., printings from beginning to end
415        nps = (int) (fabs(ftime)/fabs(ptime)); 
416
417        fprintf(flog, "nin = %d nps = %d", nin, nps);
418       
419        /* WARNING: VERY LONG INTEGRATIONS CAN PRODUCE VERY LARGE
420                INTEGERS AND ONE HAS TO MAKE SURE THAT THEY ARE STILL
421                SMALLER THAN THE LARGEST REPRESENTABLE INTEGER
422          ANOTHER WARNING: ROUNDOFF EFFECTS MAY CAUSE nin AND/OR nps TO BE
423          DIFFERENT FROM WHAT THE USER EXPECTS
424        */
425
426
427        /* print the first data point, i.e., the initial data */ 
428        #define PRINTDATA fprintf(f, "%.12f   ", time);\
429                az = angmomz(num, DIM, parts, tpstore); \
430                toten = energy(num, DIM, parts, tpstore); \
431                fprintf(f, "%.16e  ", (toten-toten0)/toten0); \
432                fprintf(f, "%.16e\n", (az-az0)/az0); \
433                printdataA(num, f, parts);
434
435        #define PRINTERROR \
436                if ( icode == CSN ) { \
437                        fprintf(flocerr, "%.12f    %.12e\n", time, CSNerror); }
438
439        PRINTDATA
440        PRINTERROR
441
442
443        /* MAIN LOOP */
444        for ( ii = 1; ii <= nps; ++ii ) {
445                for ( i = 1; i <= nin; ++i ) {
446                        switch(icode) {
447                        case SKP:
448                          SKPStep(order, stepsize, num, parts); break;
449                        case HWH:
450                          HWHStep(order, tree, stepsize, parts); break;
451                        case CSN:
452                          CStep(order, stepsize, num, parts, tree); break;
453                        case GBS:
454                          GBSStep(order, stepsize, num, parts); break;
455                        default: 
456                                fprintf(flog, "Not implemented\n");
457                                fclose(flog); fclose(f); exit(1); break;
458                        } 
459               
460                /* DONE ONE TIME STEP
461                ADD MORE CODE HERE IF NECESSARY
462                variables you can use:
463                double time;
464                double parts DECLARR;
465                Note that actual time should be computed as
466                actualtime = time+stepsize*i
467                */
468
469                }
470
471                time += stepsize*nin;
472                fprintf(f, "end of step, time = %f\n", time);
473               
474                if ( icode == HWH ) { ParticlesFromTree(tree, parts); }
475
476                /* CHANGE THIS PART OR ADD YOUR CODE IF NECESSARY */
477       
478                PRINTDATA
479                PRINTERROR
480
481                if ( stepsize > 0.0 && time > ftime ) break;
482                if ( stepsize < 0.0 && time < ftime ) break;
483
484        }
485
486
487        fclose(f);
488        fprintf(flog, "\n Finished \n");
489        fclose(flog);
490        if ( icode == CSN ) fclose(flocerr);
491       
492        return 0;
493}
494
495
496/* >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
497        FOR TEST: THE SIMPLE SPLIT, I.E., Kinetic - Potential
498        It does not provide special treatment for test particles.
499*/
500int PlainSplitKineticStep(int num, double stepsize, double p DECLARR) {
501        int i, k;
502        for ( i = 1; i <= num; ++i ) {
503                for ( k = 1; k <= DIM; ++k ) {
504                        p[i][k] += stepsize*p[i][k+DIM]; } }
505        return 0; }
506
507int PlainSplitPotentialStep(int num, double stepsize, double p DECLARR) {
508        int i, j, k;
509        double tmp, tmp2, tmp3;
510        for ( i = 1; i <= num; ++i ) {
511                tmp = 0.0;
512                for ( j = i+1; j <= num; ++j ) {
513                        tmp2 = 0.0;
514                        for ( k = 1; k <= DIM; ++k ) {
515                                tmp3 = p[i][k] - p[j][k];
516                                tmp2 = tmp2 + tmp3*tmp3; }
517                        tmp3 = 1.0/sqrt(tmp2);
518                        tmp2 = tmp3/tmp2;
519                        for ( k = 1; k <= DIM; ++k ) {
520                                tmp3 = p[i][k] - p[j][k];
521                                tmp = -tmp3*tmp2*p[j][GMASSCO]*stepsize;
522                                p[i][k+DIM] += tmp;
523                                tmp = tmp3*tmp2*p[i][GMASSCO]*stepsize;
524                                p[j][k+DIM] += tmp; }
525                        }
526        }
527        return 0; }
528
529
530
531/* >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
532        COWELL-STORMER INTEGRATOR */
533
534
535/* this computes accelerations, used both for CSN and GBS,
536        not to be confused with Accel used for HWH */
537
538int PlainAccel(int num, double p DECLARR) {
539        int i, j, k; int tpi, tpj;
540        double tmp, tmp2, tmp3;
541        double massi, massj;
542        for ( i = 1; i <= num; ++i ) {
543                for ( k = 1; k <= DIM; ++k ) { p[i][ACC+k] = 0.0; } }
544/* for testing things
545for ( i = 1; i <= num; ++i ) {
546        for ( k = 1; k <= DIM; ++k ) { p[i][ACC+k] = 1.0; } }
547return 0;
548*/
549        for ( i = num; i >= 1; --i ) { 
550                tpi = tpstore[i];
551                massi = p[i][GMASSCO];
552                for ( j = i+1; j <= num; ++j ) {
553                        tpj = tpstore[j];
554                        if ( tpi == TESTPARTICLE && tpj == TESTPARTICLE ) continue;
555                        tmp2 = 0.0;
556                        for ( k = 1; k <= DIM; ++k ) {
557                                tmp3 = p[i][k] - p[j][k];
558                                tmp2 = tmp2 + tmp3*tmp3; }
559                        tmp3 = 1.0/sqrt(tmp2);
560                        tmp2 = tmp3/tmp2;
561                        massj = p[j][GMASSCO];
562                        for ( k = 1; k <= DIM; ++k ) {
563                                tmp3 = p[i][k] - p[j][k];
564                                tmp = tmp3*tmp2;
565                        if ( tpj == NOTTESTPARTICLE ) {
566                                p[i][ACC+k] -= (tmp*massj); }
567                        if ( tpi == NOTTESTPARTICLE ) {
568                                p[j][ACC+k] += (tmp*massi); }
569                                }
570                        }
571        }
572        return 0; }
573
574/* Cowell-Stormer-Newman 14(6)th order */
575
576
577/* Newman's 14(6)th order Stormer's coeffs */
578
579double CSalpha[] =  {   0.0,
580                        1.000000000000000000000,
581                        0.000000000000000000000,
582                        0.083333333333333333333,
583                        0.083333333333333333333,
584                        0.079166666666666666666,
585                        0.075000000000000000000,
586                        0.071345899470899470899,
587                        0.068204365079365079365,
588                        0.065495756172839506172,
589                        0.063140432098765432098,
590                        0.061072649861712361712,
591                        0.059240564123376623376,
592                        0.057603625837453135733,
593                        0.056129980884507174190,
594                        0.054794379107071147139
595                                                 };
596 
597 
598double CSbeta[] =  {    0.0,
599                        1.50000000000000000000,
600                        0.33333333333333333333,
601                        0.37500000000000000000,
602                        0.35277777777777777777,
603                        0.33402777777777777777,
604                        0.31924603174603174603,
605                        0.30736607142857142857,
606                        0.29757660934744268077,
607                        0.28933077050264550264,
608                        0.28225737868098979210,
609                        0.27609762576993479771,
610                        0.27066578505957226195,
611                        0.26582499331955247133,
612                        0.26147199790503706399,
613                        0.25752728193649391773
614                                               };
615 
616       
617
618double  sumfrombottom(double *a, int nr, double *coef) {
619        int i;
620        double sum = 0.0;
621        for ( i = nr; i >= 1;  --i ) {
622                sum += a[i]*coef[i]; }
623        return sum; }
624
625int differencing(double *a, int nr, double in) {
626        int i, j, p;
627        double tmp1, tmp2;
628        tmp1 = a[1]; a[1] = in;
629        for ( i = 2; i <= nr; ++i ) {
630                tmp2 = a[i]; a[i] = a[i-1]-tmp1; tmp1 = tmp2; } 
631        return 0; }
632
633
634
635static int CStepCount = 0;
636static int CStepInit = 0;
637static struct tree *initCStree;
638
639int CStep(int ord, double h, int nump, 
640                double parts DECLARR, struct tree *treein) {
641        double *a, asav, vx, sum, tmp; int p,j, i;
642static double df DECLSTORMER;
643/* actual local order to determine how terms to sum */
644if ( CStepInit == 0 ) {
645        initCStree = treein;
646        CStepInit = 1;
647        for ( i = 1; i <= ord+1; ++i ) {
648        for ( p = 1; p <= nump; ++p ) {
649        for ( j = 1; j <= DIM; ++j ) {
650                df[p][j][i] = 0.0; } } }        }
651if ( CStepInit == 1 ) {
652/* initialization stage */
653                for ( p = 1; p <= nump; ++p ) {
654                for ( j = 1; j <= DIM; ++j ) {
655                        df[p][j][0] = -parts[p][j]; } }
656
657                if ( csninit == CSNINITHWH ) {
658                        HWHStep(CSNINITHWHORDER, initCStree, h, parts); 
659                        ParticlesFromTree(initCStree, parts);  }
660                else { GBSStep(CSNINITGBSORDER, h, nump, parts); }
661                if ( CStepCount > ord ) CStepInit = 2;
662                for ( p = 1; p <= nump; ++p ) {
663                for ( j = 1; j <= DIM; ++j ) {
664                        df[p][j][0] += parts[p][j]; } }
665                }
666else {
667/* actual run */
668CSNerror = 0;
669for ( p = 1; p <= nump; ++p ) {
670for ( j = 1; j <= DIM; ++j ) {
671/* formula: x(n+1) - x(n) = (x(n) - x(n-1)) + h^2*sum
672        with X(n) = x(n)-x(n-1)
673        X(n+1) = X(n) + h^2*sum
674        x(n+1) = X(N+1) + x(n)
675X(n) is in a[0]
676*/
677        a = df[p][j];
678        asav = a[0];
679        sum = sumfrombottom(a, ord, CSalpha);
680        sum *= h*h;
681        a[0] += sum;
682        parts[p][j] += a[0];
683/* vx = (1.0/h)*(x(n)-x(n-1)) + h*sum */
684        sum = h*sumfrombottom(a, ord, CSbeta);
685        vx = asav/h + sum;
686        parts[p][DIM+j] = vx;
687/* record error */
688        tmp = fabs(CSalpha[ord]*a[ord]*h*h);
689        if ( CSNerror < tmp ) CSNerror = tmp;
690        } } }
691PlainAccel(nump, parts);
692for ( p = 1; p <= nump; ++p ) {
693for ( j = 1; j <= DIM; ++j ) {
694        differencing(df[p][j], ord, parts[p][ACC+j]);
695        parts[p][ACC+j] = 0.0;
696        } }
697++CStepCount;
698return 0; }
699
700
701
702/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
703        BULIRSCH-STOER
704        actually, GRAGG-BULIRSCH-STOER, following Hairer et al., 1993
705*/
706
707/* Gragg's symmetric integrator, this is used to generate Tj,1  */
708int GStep(int nump, double parts DECLARR, double s DECLARR, double h, int len) {
709int k, j, p, d, knext, kprev, m;
710double dif, difsum, tmp1, astep;
711static double y0 DECLARR;
712static double y1 DECLARR;
713/* save originial in y0 */
714for ( p = 1; p <= nump; ++p ) {
715                y0[p][GMASSCO] = parts[p][GMASSCO];
716for ( d = 1; d <= DIM; ++d ) { 
717                y0[p][d] = parts[p][d]; 
718                y0[p][d+DIM] = parts[p][d+DIM]; } }
719/* consult Hairer et al. why we need this
720        in short: we need a method which has only even powers
721        of h in the expansion of its error, i.e., a symmetric method
722 */
723h = h/2.0;
724PlainAccel(nump, y0);
725for ( p = 1; p <= nump; ++p ) {
726/* load masses into y1 */
727                y1[p][GMASSCO] = y0[p][GMASSCO];
728for ( d = 1; d <= DIM; ++d ) { 
729                y1[p][d] = y0[p][d]+h*y0[p][DIM+d]; 
730                y1[p][d+DIM] = y0[p][d+DIM]+h*y0[p][ACC+d]; } }
731/* compute y2 (m=1), y3' (m=2), y4 (m=3), y5' (m=4) y(2n+1)' (m=2*n)
732        n = len
733        */
734for ( m = 1; m <= 2*len; ++m ) {
735        if ( (m%2) ) {
736        for ( p = 1; p <= nump; ++p ) {
737        for ( d = 1; d <= DIM; ++d ) { 
738                y1[p][d] = y0[p][d]+2.0*h*y1[p][d+DIM]; 
739                s[p][d] = y1[p][d]; y0[p][d] = y1[p][d]; 
740                }  } }
741        else {
742        PlainAccel(nump, y1);
743        for ( p = 1; p <= nump; ++p ) {
744        for ( d = 1; d <= DIM; ++d ) { 
745                y0[p][d+DIM] = y1[p][d+DIM]; 
746                y1[p][d+DIM] = y0[p][d+DIM]+2.0*h*y1[p][ACC+d]; 
747                s[p][d+DIM] = (y0[p][d+DIM]+y1[p][d+DIM])/2.0; 
748                } } 
749                }
750        }
751return 0; }
752
753
754int GBSStep(int bsord, double h, int nump, double parts DECLARR) {
755int k, j, jfinal, p, d, jnext, jprev, len;
756static double gbsarr DECLBS;
757double dif, difsum, tmp1, astep;
758int *n = extrapolseq;
759static double s DECLARR;
760/* for testing GStep
761        astep = h; len = 1;
762        GStep(nump, parts, s, astep, len);
763for ( p = 1; p <= nump; ++p ) {
764for ( d = 1; d <= 2*DIM; ++d ) { parts[p][d] = s[p][d]; } }
765return 0;
766*/
767jfinal = bsord;
768for ( j = 1; j <= bsord; ++j ) {
769        difsum = 0.0;
770/* map j and j-1 to storage index */
771        jnext = ((j)%2); 
772        jprev = ((j-1)%2);
773/* compute Tj,1 */
774        len = n[j];
775        astep = h/len;
776        GStep(nump, parts, s, astep, len);
777        for ( p = 1; p <= nump; ++p ) {
778        for ( d = 1; d <= 2*DIM; ++d ) { gbsarr[p][d][jnext][1] = s[p][d]; } }
779/* compute Tj,k */
780for ( k = 1; k <= j-1; ++k ) {
781/* Aitken-Neville for nonsymmetric underlying method */
782        tmp1 = ((double) n[j])/((double) n[j-k]); 
783/* Aitken-Neville for symmetric underlying method (square of previous) */
784        tmp1 = tmp1*tmp1;
785for ( p = 1; p <= nump; ++p ) {
786for ( d = 1; d <= 2*DIM; ++d ) {
787        dif = (gbsarr[p][d][jnext][k] - gbsarr[p][d][jprev][k]);
788        if ( k == j-1 ) difsum += fabs(dif);
789        gbsarr[p][d][jnext][k+1] = gbsarr[p][d][jnext][k] + dif/(tmp1-1.0);
790        } }
791        } 
792/* stop computing more is reached roundoff level */
793        if ( j > 1 && difsum <= GBSTOL ) { jfinal = j; break; }
794        }
795j = jfinal; k = j; jnext = ((j)%2); 
796for ( p = 1; p <= nump; ++p ) {
797for ( d = 1; d <= 2*DIM; ++d ) { parts[p][d] = gbsarr[p][d][jnext][k]; } } 
798return 0; }
799       
800
801
802
803
804
805/* >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
806        Actual  HWH integration */
807int MoveCMandKepler(struct tree *tree, double *cmove, double stepsize, 
808                double p DECLARR);
809int KeplerOnTree(struct tree *tree, double stepsize, double p DECLARR);
810int RecursivePerturbations(struct tree *tree, double stepsize, double p DECLARR);
811int RecPertAcc(struct tree *tree, double stepsize, double p DECLARR);
812
813/* fourth order coeffs */
814static double Coeff4x0 = 0.0;
815static double Coeff4x1 = 0.0;
816
817static int SympCalled = 0;
818
819void initSympCoefs() {
820        double cbrt2 = cbrt(2.0);
821        Coeff4x0 = -cbrt2/(2.0-cbrt2);
822        Coeff4x1 = 1.0/(2.0-cbrt2);
823        SympCalled = 1; }
824       
825
826static int HWHCalled = 0;
827
828// X
829int HWHStep(int order, struct tree *tree, double stepsize, double p DECLARR) {
830        if ( SympCalled == 0 ) { 
831                initSympCoefs(); 
832        }
833        if ( HWHCalled == 0) {
834                // initialize tree
835                FillCMinTree(tree, p); 
836                HWHCalled = 1; 
837        }
838        if ( order == 1 ) {
839                KeplerOnTree(tree, stepsize, p);
840                RecursivePerturbations(tree, stepsize, p);
841                return 0; 
842        }
843        if ( order == 2 ) {
844                KeplerOnTree(tree, stepsize/2.0, p);
845                RecursivePerturbations(tree, stepsize, p);
846                KeplerOnTree(tree, stepsize/2.0, p); 
847                return 0; 
848        }
849        if ( order == 4 ) {
850                HWHStep(2, tree, stepsize*Coeff4x1, p);
851                HWHStep(2, tree, stepsize*Coeff4x0, p);
852                HWHStep(2, tree, stepsize*Coeff4x1, p);
853                return 0; 
854        }
855        printf(" Order %2d  is not implemented, Aborted", order);
856        exit(1); 
857}
858
859
860int SKPStep(int order, double stepsize, int nump, double p DECLARR) {
861        if ( SympCalled == 0 ) { initSympCoefs(); }
862        if ( order == 1 ) {
863                PlainSplitKineticStep(nump, stepsize, p);
864                PlainSplitPotentialStep(nump, stepsize, p);
865                return 0; }
866        if ( order == 2 ) {
867                PlainSplitKineticStep(nump, stepsize/2.0, p);
868                PlainSplitPotentialStep(nump, stepsize, p);
869                PlainSplitKineticStep(nump, stepsize/2.0, p);
870                return 0; }
871        if ( order == 4 ) {
872                SKPStep(2, stepsize*Coeff4x1, nump, p);
873                SKPStep(2, stepsize*Coeff4x0, nump, p);
874                SKPStep(2, stepsize*Coeff4x1, nump, p);
875                return 0; }
876        printf(" Order %2d  is not implemented, Aborted", order);
877        exit(1); }
878
879
880
881/* +++++++++++++++++++++++++++++++++++++++++++++++++++++
882/       TREE
883/
884/
885*/
886
887
888int inittree(struct tree *tmp) {
889        int i; tmp->tmass = 0.0; tmp->particleindex = 0;
890        tmp->p1 = 0; tmp->p2 = 0;
891        tmp->testparticles = 0; tmp->rm1 = 0; tmp->rm2 = 0;
892        for ( i = 1; i <= DIM2; ++i ) { tmp->bc[i] = 0; }
893        for ( i = 1; i <= DIM; ++i ) { tmp->acc[i] = 0; }
894        return 0; }
895
896
897struct tree *newtree() {
898        struct tree *tmp = (struct tree*) malloc(sizeof(struct tree)); 
899        inittree(tmp); return tmp; }
900
901void deletetree(struct tree *p) { 
902        if ( p->p1 ) deletetree(p->p1);
903        if ( p->p2 ) deletetree(p->p2);
904        free((void *) p); }
905
906
907struct tree *MakeSimpleTree(int nump) {
908        int i;
909        struct tree *root;
910        struct tree *t1, *p1, *p2;
911        root = newtree();
912        t1 = root;
913        for ( i = nump; i > 1; --i ) {
914                p1 = newtree();
915                p2 = newtree();
916                p1->particleindex = 0;
917                p2->particleindex = i;
918                t1->p1 = p1; t1->p2 = p2;
919                t1 = p1; }
920        t1->particleindex = 1;
921        return root; }
922
923int eatwhite(FILE *f) {
924        int c; 
925        for ( ; ; ) {
926                c = getc(f);
927                if ( !isspace(c) ) break; 
928        }
929        ungetc(c, f); 
930        return c;
931}
932
933
934/* just an array to hold indeces read in,
935        portions of this are assigned to
936        tpa's below */
937int tmpReadList[MAXNBODIES+1000];
938/* bookkeeping */
939int freelist = 0;
940
941int *ReadList(FILE *f) {
942        int i1, i2; int c; int i; int *arr;
943        c = eatwhite(f);
944        if ( c != '[' ) {
945                printf(" Error: bad list, aborted \n"); exit(1); 
946        }
947        c = getc(f);
948        fscanf(f, "%d - %d", &i1, &i2);
949        c = eatwhite(f);
950        if ( c != ']' ) {
951                printf(" Error: bad list, aborted \n"); exit(1); 
952        }
953        c = getc(f); 
954        arr = tmpReadList+freelist;
955       
956        /* put zero in tpstore */
957        for ( i = i1; i <= i2; ++i ) tpstore[i] = TESTPARTICLE;
958        for ( i = i1; i <= i2; ++i ) arr[i-i1+1] = i;
959       
960        arr[0] = i2-i1+1;
961        freelist += i2-i1+2;
962        return arr; 
963}
964       
965               
966
967struct tree *MakeTree(FILE *f) {
968        int i; int c; int pind; int *tpa;
969        struct tree *root;
970        struct tree *p1, *p2;
971        c = eatwhite(f);
972        if ( c != '(' ) {
973                if ( isdigit(c) ) {
974                        // single particle
975                        fscanf(f,"%d", &pind);
976                        root = newtree();
977                        root->particleindex = pind;
978                        return root;
979                } 
980                // set of test particles
981                tpa = ReadList(f);
982                root = newtree();
983                root->particleindex = tpa[1];
984                root->testparticles = tpa;
985                return root;
986        }
987        c = getc(f);
988        root = newtree();
989        p1 = MakeTree(f);
990        p2 = MakeTree(f);
991        root->p1 = p1; root->p2 = p2;
992        root->particleindex = 0;
993        c = eatwhite(f);
994       
995        if ( c != ')' ) {
996                printf(" Error: bad tree, aborted \n"); exit(1); 
997        }
998        c = getc(f);
999        return root; 
1000}
1001
1002
1003
1004/* not used now but might come handy */
1005struct tree *getparticleintree(struct tree *tree, int ind) {
1006        struct tree *tmp;
1007        if ( tree->particleindex ) { 
1008                if ( tree->particleindex == ind ) { return tree; }
1009                return 0; }
1010        tmp = getparticleintree(tree->p1, ind);
1011        if ( tmp != 0 ) return tmp;
1012        tmp = getparticleintree(tree->p2, ind);
1013        if ( tmp != 0 ) return tmp;
1014        return 0; }
1015       
1016
1017/* not used now but might come handy */
1018int **GetTestParticles(struct tree *tree, int ind) {
1019        struct tree *tmp = getparticleintree(tree, ind);
1020        if ( tmp == 0 ) return 0;
1021        return &(tmp->testparticles); }
1022       
1023
1024/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1025        ACTUAL NUMERICAL STUFF ON TREE
1026*/
1027
1028// X fills in initial data: masses, mass ratios, barycenters
1029int FillCMinTree(struct tree *tree, double p DECLARR) {
1030        int j; double *bc = tree->bc;
1031        int i, ind, tp, num, *tpa = tree->testparticles;
1032       
1033        #pragma omp parallel for private(j) schedule(runtime)
1034        for ( j = 1; j <= DIM; ++j ) 
1035                tree->acc[j] = 0.0;
1036       
1037        tree->tmass = 0.0; tree->rm1 = 0.0; tree->rm2 = 0.0;
1038        #pragma omp parallel for private(j) schedule(runtime)
1039        for ( j = 1; j <= DIM2; ++j ) 
1040                bc[j] = 0.0;
1041               
1042        ind = tree->particleindex;
1043        if ( ind ) {
1044            if ( tpa != 0 ) { 
1045                        num = tpa[0];
1046                        for ( tp = 1; tp <= num; ++tp ) { 
1047                                ind = tpa[tp]; 
1048                                if ( ind < 1 ) 
1049                                        continue;
1050                                #pragma omp parallel for private(i) schedule(runtime)
1051                                for ( i = 1; i <= DIM; ++i ) 
1052                                        p[ind][ACC+i] = 0.0;
1053                                } 
1054                        } 
1055            else { 
1056                        tree->tmass = p[ind][GMASSCO];
1057                        #pragma omp parallel for private(j) schedule(runtime)
1058                        for ( j = 1; j <= 2*DIM; ++j ) 
1059                                bc[j] = p[ind][j]; 
1060                }
1061            return 0; 
1062        } 
1063        FillCMinTree(tree->p1, p);
1064        FillCMinTree(tree->p2, p);
1065        tree->tmass = tree->p1->tmass + tree->p2->tmass;
1066        tree->rm1 = (tree->p1->tmass)/(tree->tmass); 
1067        tree->rm2 = (tree->p2->tmass)/(tree->tmass);
1068       
1069        #pragma omp parallel for private(j) schedule(runtime)
1070        for ( j = 1; j <= DIM2; ++j ) {
1071                bc[j] = tree->p1->bc[j] * tree->rm1 + tree->p2->bc[j] * tree->rm2; 
1072        } 
1073        return 0; 
1074}
1075
1076/* copy from tree into p */
1077int ParticlesFromTree(struct tree *tree, double p DECLARR) {
1078        int j; double *bc = tree->bc;
1079        int i, ind, *tpa = tree->testparticles;
1080        ind = tree->particleindex;
1081        if ( ind ) {
1082           if ( tpa == 0 ) {
1083                tree->tmass = p[ind][GMASSCO];
1084                for ( j = 1; j <= DIM2; ++j ) p[ind][j] = bc[j]; }
1085           return 0; } 
1086        ParticlesFromTree(tree->p1, p);
1087        ParticlesFromTree(tree->p2, p);
1088        return 0; }
1089
1090
1091/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1092        KEPLER STEP ON TREE
1093*/
1094
1095/* computes Stumpff functions several ways, depending on x */
1096int StumpffFunctions(double *vals, double x);
1097/* advances the Kepler motion 1/2*v^2 - mu/|x| for time deltat */
1098double KeplerStumpffStep(int dim, double mu, double *xv, double deltat);
1099/* return only the increment in xv */
1100double IncrementKeplerStumpffStep(int dim, double mu, double *xv, double deltat);
1101
1102// X NOTE: this does not CHANGE centers of masses
1103int KeplerOnTree(struct tree *tree, double stepsize, double p DECLARR) { 
1104        // move the center of mass
1105        double cmove[10]; 
1106        int i;
1107        #pragma omp parallel for private(i) schedule(runtime)
1108        for ( i = 1; i <= 2*DIM; ++i ) 
1109                cmove[i] = 0.0;
1110
1111        // move by 1/2*(Tmass)*v^2 */
1112        #pragma omp parallel for private(i) schedule(runtime)
1113        for ( i = 1; i <= DIM; ++i ) 
1114                cmove[i] = (tree->bc[DIM+i])*stepsize; 
1115
1116        MoveCMandKepler(tree, cmove, stepsize, p);
1117        return 0; 
1118}
1119
1120
1121// X Assume that test particles are attached to 2nd
1122int DoTestPartKepler(struct tree *tree, double *cmove, double stepsize, double p DECLARR) {
1123        int i, tp, num, *tpa, ind;
1124        double mu, tmpkep[TMPSIZE], *bc;
1125        tpa = tree->p2->testparticles; 
1126        num = tpa[0];
1127        if ( num > 0 ) {
1128                bc = tree->p1->bc; 
1129                mu = tree->tmass;
1130                for ( tp = 1; tp <= num; ++tp ) { 
1131                        ind = tpa[tp]; 
1132                        if ( ind < 1 ) continue;
1133                        #pragma omp parallel for private(i) schedule(runtime)
1134                        for ( i = 1; i <= DIM2; ++i )   
1135                                tmpkep[i] = p[ind][i] - bc[i];
1136                        IncrementKeplerStumpffStep(DIM, mu, tmpkep, stepsize); 
1137                        for ( i = 1; i <= DIM2; ++i ) 
1138                                p[ind][i] += (cmove[i]+tmpkep[i]);
1139                } 
1140        }
1141        return 0; 
1142}
1143
1144// X move the center of mass and do a Kepler step on relative coordinates, transmit the new "moves" down
1145int MoveCMandKepler(struct tree *tree, double *cmove, double stepsize, double p DECLARR) {
1146        int i, *tpa, ind;
1147        double tmp[TMPSIZE], tmpkep[TMPSIZE], *xv1, *xv2, mu;
1148        double *bc = tree->bc;
1149       
1150        for ( i = 1; i <= DIM2; ++i ) { 
1151                bc[i] += cmove[i]; 
1152        }
1153        // do nothing else for actual particles
1154        if ( tree->particleindex ) { 
1155                return 0; 
1156        }
1157        if ( tree->p2->testparticles ) {       
1158                // this is a pair where p2 is an array of test particles
1159                DoTestPartKepler(tree, cmove, stepsize, p); 
1160                MoveCMandKepler(tree->p1, cmove, stepsize, p); 
1161                return 0; 
1162        }
1163        mu = tree->tmass;
1164        xv1 = tree->p1->bc; 
1165        xv2 = tree->p2->bc;
1166        #pragma omp parallel for private(i) schedule(runtime)
1167        for ( i = 1; i <= DIM2; ++i ) { 
1168                tmpkep[i] = xv2[i]-xv1[i]; 
1169        }
1170        // do Kepler step on tmpkep, we get back the increment
1171        IncrementKeplerStumpffStep(DIM, mu, tmpkep, stepsize); 
1172        #pragma omp parallel for private(i) schedule(runtime)
1173        // use only the relative change
1174        for ( i = 1; i <= DIM2; ++i ) { 
1175                tmp[i] = cmove[i] - tmpkep[i]*(tree->rm2); 
1176        }
1177        MoveCMandKepler(tree->p1, tmp, stepsize, p);
1178        #pragma omp parallel for private(i) schedule(runtime)
1179        for ( i = 1; i <= DIM2; ++i ) { 
1180                tmp[i] = cmove[i] + tmpkep[i]*(tree->rm1); 
1181        }
1182        MoveCMandKepler(tree->p2, tmp, stepsize, p);
1183        return 0; 
1184}
1185       
1186/*  ++++++++++++++++++++++++++++++++++++++++++++++++++++++
1187        PERTURBATION ACCELERATIONS
1188*/
1189
1190/* X simple acceleration */
1191int Accel(double *x1, double *x2, double *acc) {
1192        double tmp[10], r = 0.0; int i; 
1193        for ( i = 1; i <= DIM; ++i ) {
1194                tmp[i] = x1[i]-x2[i];
1195                r = r + tmp[i]*tmp[i]; 
1196        }
1197        r = sqrt(r);
1198        r = r*r*r;
1199        for ( i = 1; i <= DIM; ++i ) { 
1200                acc[i] = -tmp[i]/r; 
1201        }
1202        return 0; 
1203}
1204
1205int AdvanceTestParticleVel(struct tree *tree, double stepsize, double p DECLARR) {
1206                int i, ind, tp, num, *tpa = tree->testparticles;
1207                if ( tpa == 0 ) return 0; num = tpa[0];
1208        for ( tp = 1; tp <= num; ++tp ) {
1209                ind = tpa[tp]; if ( ind < 1 ) continue;
1210                for ( i = 1; i <= DIM; ++i ) {
1211                        p[ind][DIM+i] += stepsize*p[ind][ACC+i]; 
1212                        p[ind][ACC+i] = 0.0; }
1213                        } 
1214        return 0; }
1215                       
1216       
1217/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1218        PERTURBATIONS ON TREE
1219*/
1220// X
1221int RecursivePerturbations(struct tree *tree, double stepsize, double p DECLARR) {
1222        int j, ind; double *bc = tree->bc;
1223        /* advance attached test particles
1224        if ( tree->particleindex ) {
1225        /* if particle we are ready to advance v */
1226        if ( tree->testparticles ) {
1227                AdvanceTestParticleVel(tree, stepsize, p); 
1228                return 0; 
1229        }
1230        // this is a nonzero mass actual particle, advance v
1231                for ( j = 1; j <= DIM; ++j ) {
1232                        bc[DIM+j] += stepsize*tree->acc[j]; 
1233                        tree->acc[j] = 0.0; 
1234                }
1235                return 0; 
1236       
1237        /* if not, first we take care of accs. for pairs with the two parts 
1238                IMPORTANT TO DO THIS FIRST
1239        */
1240        RecPertAcc(tree, stepsize, p);
1241        // now descend and do the same for both parts
1242        RecursivePerturbations(tree->p2, stepsize, p);
1243        RecursivePerturbations(tree->p1, stepsize, p);
1244        /* now regenerate center of mass etc. */
1245        if ( tree->p2->testparticles ) {
1246                #pragma omp parallel private(j)
1247                #pragma omp for schedule(runtime)
1248                for ( j = 1; j <= DIM2; ++j ) { 
1249                        bc[j] = tree->p1->bc[j]; 
1250                }
1251                return 0; 
1252        }
1253       
1254        #pragma omp parallel private(j)
1255        #pragma omp for schedule(runtime)
1256        for ( j = 1; j <= DIM2; ++j ) {
1257                bc[j] = tree->p1->bc[j] * tree->rm1 + tree->p2->bc[j] * tree->rm2; 
1258        }
1259        return 0; 
1260}
1261
1262/* assume test particles in tree2 for dir > 1, reverse signs for dir < 0 */
1263int TestParticleAcc(int dir, double *cacc, struct tree *tree1, struct tree *tree2, 
1264                                double p DECLARR) {
1265                int i, ind, tp, num, *tpa = tree2->testparticles;
1266                double acc[TMPSIZE];
1267                double *bc = tree1->bc;
1268                double tmp, gm1 = tree1->tmass;
1269                num = tpa[0];
1270        for ( tp = 1; tp <= num; ++tp ) {
1271                ind = tpa[tp]; if ( ind < 1 ) continue;
1272/* rewrite this to avoid bad numerical attitude */
1273                        Accel(bc, &p[ind][0], acc);
1274                for ( i = 1; i <= DIM; ++i ) {
1275/* since test particles are in 2nd */
1276                        tmp = acc[i]-cacc[i]*dir;
1277                        p[ind][ACC+i] += -gm1*tmp; }
1278                        } 
1279        return 0; }
1280                       
1281
1282/* X c1 and c2 are for possible future use */
1283int RecPertAccSub(struct tree *tree1, 
1284        struct tree *tree2, double *cacc, double *c1, double *c2, double p DECLARR) {
1285        /* do the work if both are particles */
1286        int i, ind1, ind2; double acc[TMPSIZE];
1287        double tmp, gm1 = tree1->tmass, gm2 = tree2->tmass;
1288        ind1 = tree1->particleindex;
1289        ind2 = tree2->particleindex;
1290        if ( ind1 && ind2 ) {
1291                /* do nothing if both are test particle arrays */
1292                if ( tree1->testparticles  && tree2->testparticles  ) 
1293                        return 0; 
1294                /* test particles first */
1295                if ( tree2->testparticles ) {
1296                        TestParticleAcc(1, cacc, tree1, tree2, p); 
1297                        return 0; 
1298                }
1299                if ( tree1->testparticles ) {
1300                        TestParticleAcc( -1, cacc, tree2, tree1, p); 
1301                        return 0; 
1302                }
1303                /* if actual particles with nonzero mass */
1304                Accel(tree1->bc, tree2->bc, acc);
1305                for ( i = 1; i <= DIM; ++i ) {
1306                        /* rewrite this to avoid bad numerical attitude */
1307                        tmp = acc[i]-cacc[i];
1308                        tree1->acc[i] += gm2*tmp; 
1309                        tree2->acc[i] += -gm1*tmp; 
1310                }
1311                return 0; 
1312        }
1313        /* if not, descend but only for one branch, the other is taken care in turn */
1314        if ( !ind1 ) {
1315                RecPertAccSub(tree1->p2, tree2, cacc, c1, c2, p);
1316                RecPertAccSub(tree1->p1, tree2, cacc, c1, c2, p);
1317                return 0; 
1318        }
1319        RecPertAccSub(tree1, tree2->p1, cacc, c1, c2, p);
1320        RecPertAccSub(tree1, tree2->p2, cacc, c1, c2, p);
1321        return 0; 
1322}
1323
1324/* assume test particles in tree2 */
1325int TestParticleAccSpec(double *cm, struct tree *tree1, struct tree *tree2, double p DECLARR) {
1326        int i, ind, tp, num, *tpa = tree2->testparticles;
1327        double acc[TMPSIZE], cacc[TMPSIZE];
1328        double *bc = tree1->bc;
1329        double tmp, gm1 = tree1->tmass;
1330        num = tpa[0];
1331        for ( tp = 1; tp <= num; ++tp ) {
1332                ind = tpa[tp]; if ( ind < 1 ) continue;
1333                /* rewrite this to avoid bad numerical attitude */
1334                Accel(bc, &p[ind][0], acc);
1335                Accel(cm, &p[ind][0], cacc);
1336                for ( i = 1; i <= DIM; ++i ) {
1337                        /* since test particles are in 2nd */
1338                        tmp = acc[i]-cacc[i];
1339                        p[ind][ACC+i] += -gm1*tmp; 
1340                }
1341        } 
1342        return 0; 
1343}
1344                       
1345
1346/* c1 and c2 are for possible future use */
1347int RecPertAccSubTest(double *cm,  struct tree *tree1, struct tree *tree2, double p DECLARR) {
1348        int ind1;
1349        ind1 = tree1->particleindex;
1350        if ( ind1 ) {
1351                /* do nothing if tree1 has only test particles */
1352                if ( tree1->testparticles ) 
1353                        return 0;
1354                TestParticleAccSpec(cm, tree1, tree2, p); 
1355                return 0; 
1356        }
1357        /* if not, descend for both branches */
1358        RecPertAccSubTest(cm, tree1->p1, tree2, p);
1359        RecPertAccSubTest(cm, tree1->p2, tree2, p);
1360        return 0; 
1361}
1362
1363
1364int RecPertAcc(struct tree *tree, double stepsize, double p DECLARR) {
1365        double cacc[TMPSIZE];
1366        /* do nothing if all in p2 are attached to p1: pure Kepler motion */
1367        if ( tree->p1->particleindex &&  tree->p2->particleindex  ) 
1368                return 0; 
1369        if ( tree->p2->testparticles == 0 ) {
1370                /* actual nonzero masses in p2 */
1371                Accel(tree->p1->bc, tree->p2->bc, cacc);
1372               
1373                RecPertAccSub(tree->p1, tree->p2, cacc, tree->p1->bc, tree->p2->bc, p); 
1374                return 0;
1375        }
1376        /* test particles in p2 */
1377        RecPertAccSubTest(tree->p1->bc, tree->p1, tree->p2, p); 
1378        return 0; 
1379}
1380
1381
1382
1383
1384/*
1385// ===============================================================================
1386//
1387//   KEPLER MOTION USING STUMPFF FUCNTIONS
1388//      FOLLOWING DANBY, 1992
1389//      (quartic Newtonian solver, Stumpff function definitions,
1390//       f and g functions)
1391//      see also: Stiefel and Sheifele, Linear and Regular
1392//      Celestial Mechanics
1393//
1394// ===============================================================================
1395
1396// Stumpff functions (c0, c1, c2, c3)
1397*/
1398
1399/* CHANGE these numbers to tune accuracy */
1400/* this one can be changed runtime */
1401static double StumpffQSTol = KEPLERSOLVERTOLERANCE;
1402
1403/* the limit for which trig and hyp. trig functions are used */
1404#define STUMPFFPLAINLIM (5.0)
1405/* reduction limit */
1406#define STUMPFFRED (1.0/32.0)
1407#define STUMPFFSHORTTOL (1.0e-6)
1408
1409/* this one is not used now */
1410static double StumpffSerTol = 1.0e-25;
1411
1412static int StumpffErrorCount = 0;
1413static int StumpffErrorCountLimit = 100;
1414
1415CheckStumpffErrorCount() {
1416        ++StumpffErrorCount;
1417        if ( StumpffErrorCount > StumpffErrorCountLimit ) {
1418                fprintf(FLOG, "\n Too many errors, aborted \n");
1419                fflush(FLOG);
1420                exit(1); 
1421        }
1422        return 0; 
1423}
1424
1425double util_fact(int k) {
1426        int i; double fact = 1;
1427        for ( i = 1; i <= k; ++i ) {
1428                fact *= i; 
1429        }
1430        return fact; 
1431}
1432
1433/*
1434In case we need the coefficients, produced by Mathematica
1435static double c2_0 = 0.5;
1436static double c2_1 = -1.0/24;
1437static double c2_2 = 1.0/720;
1438static double c2_3 = -1.0/40320;
1439static double c2_4 = 1.0/3628800;
1440static double c2_5 = -1.0/479001600;
1441static double c2_6 = 1.0/87178291200;
1442static double c2_7 = -1.0/20922789888000;
1443static double c2_8 = 1.0/6402373705728000;
1444static double c2_9 = -1.0/2432902008176640000;                         
1445static double c2_10 = 1.0/1124000727777607680000;
1446
1447static double c3_0 = 1.0/6;
1448static double c3_1 = -1.0/120;
1449static double c3_2 = 1.0/5040;
1450static double c3_3 = -1.0/362880;
1451static double c3_4 = 1.0/39916800;
1452static double c3_5 = -1.0/6227020800;
1453static double c3_6 = 1.0/1307674368000;
1454static double c3_7 = -1.0/355687428096000;
1455static double c3_8 = 1.0/121645100408832000;
1456static double c3_9 = -1.0/51090942171709440000;
1457static double c3_10 = 1.0/25852016738884976640000;
1458*/
1459
1460
1461/* This is what is actually used, generated by factoring the coefficients */
1462double HornerStumpff2(double x) {
1463        return
1464(1.0 - x*(1 - x*(1 - x*(1 - x*(1 - x *(1 - x
1465/* up to degree 10, if you need it
1466*(1 - x*(1 - x*(1 - (1 - x/462)*x/380)/306)/240)        */
1467/182)/132)/90)/56)/30)/12)/2;
1468        }
1469
1470double HornerStumpff3(double x) {
1471        return
1472 (1.0 - x*(1 - x*(1 - x*(1 - x*(1 - x*(1 - x
1473/* up to degree 10, if you need it
1474*(1 - x*(1 - x*(1 - (1 - x/506)*x/420)/342)/272)        */
1475/210)/156)/110)/72)/42)/20)/6;
1476        }
1477
1478
1479/* first 4 terms, ie., up to x^3 */
1480double Stumpff2Short(double x) {
1481        return (1.0 - x*(1 - x*(1 - x/56)/30)/12)/2;
1482        }
1483
1484double Stumpff3Short(double x) {
1485        return (1.0 - x*(1 - x*(1 - x/72)/42)/20)/6;
1486        }
1487
1488
1489
1490
1491// Original for reference
1492double StumpffSeriesCk(int k, double x) {
1493        int i;
1494        double term, sum = 0.0;
1495        double terms[100], num;
1496        term = 1.0;
1497        for ( i = 1; i <= 50; ++i ) {
1498                num = (k+2*i-1)*(k+2*i);
1499                term = term*(-x/num);
1500                terms[i] = term;
1501                if ( fabs(term) < StumpffSerTol ) break; 
1502                }
1503        for ( ; i > 0; --i )  sum += terms[i];
1504        sum = (1.0+sum)/util_fact(k);
1505        return sum; }
1506
1507// X
1508int StumpffFunctions(double *c, double x) {
1509        double x2, tmp;
1510        int i;
1511/* What this is supposed to be
1512c[0] = StumpffSeriesCk(0, x);
1513c[1] = StumpffSeriesCk(1, x);
1514c[2] = StumpffSeriesCk(2, x);
1515c[3] = StumpffSeriesCk(3, x);
1516return 0;
1517*/
1518/* FOR REFERENCE, THIS COULD BE USED FOR LARGE VALUES, MAYBE
1519        if ( x > STUMPFFPLAINLIM ) {
1520                x2 = sqrt(x);
1521                c[0] = cos(x2);
1522                c[1] = sin(x2)/x2;
1523                c[2] = (1.0-cos(x2))/x;
1524                c[3] = (x2-sin(x2))/(x*x2);
1525                return 0; }
1526        if ( x < -STUMPFFPLAINLIM ) {
1527                x2 = sqrt(-x);
1528                c[0] = cosh(x2);
1529                c[1] = sinh(x2)/x2;
1530                c[2] = (cosh(x2)-1.0)/x;
1531                c[3] = (sinh(x2)-x2)/(x*x2);
1532                return 0; }
1533*/
1534        if ( fabs(x) < STUMPFFSHORTTOL ) {
1535                c[2] = Stumpff2Short(x); 
1536                c[3] = Stumpff3Short(x);
1537                c[1] = 1.0 - x*c[3]; 
1538                c[0] = 1.0 - x*c[2];
1539                return 0;
1540        }
1541        if ( fabs(x) < STUMPFFRED ) {
1542                c[2] = HornerStumpff2(x); 
1543                c[3] = HornerStumpff3(x); 
1544                c[1] = 1.0 - x*c[3]; 
1545                c[0] = 1.0 - x*c[2];
1546                return 0;
1547        }
1548        // reduce value to small
1549        x2 = x;
1550        for ( i = 1; i <= 20; ++i ) {
1551                x2 = x2/4.0;
1552                if ( fabs(x2) < STUMPFFRED ) 
1553                        break; 
1554        }
1555        if ( i > 20 ) {
1556                fprintf(FLOG, "%s"," StumpfReduction No Good \n"); 
1557                CheckStumpffErrorCount();
1558                fflush(FLOG);
1559                return 0; 
1560        } 
1561        /* original
1562                c[2] = StumpffSeriesCk(2, x2);
1563                c[3] = StumpffSeriesCk(3, x2);
1564        */
1565        c[2] = HornerStumpff2(x2);
1566        c[3] = HornerStumpff3(x2);
1567        c[1] = 1.0 - x2*c[3];
1568        c[0] = 1.0 - x2*c[2];
1569        for ( ; i > 0; --i ) {
1570                tmp = c[2];
1571                c[2] = (c[1]*c[1])/2.0; 
1572                c[3] = (tmp+c[0]*c[3])/4.0; 
1573                c[1] = c[0]*c[1];
1574                c[0] = (2.0*c[0]*c[0])-1.0;
1575        }
1576        return i; 
1577}
1578
1579
1580/* ===============================================================================
1581*/
1582
1583/* NOTE:
1584It seems that there might be a systematic error here, not sure yet.
1585The Newton solvers tend to end up systematically on a preferred side
1586of the solution.  E.g., for x^2: if starting from above the solution,
1587it ends at a slightly larger value than the solution.  Starting from
1588 below, in the first step it will get to the previous case.  Since
1589the solution cannot be computed to machine precision, we could have
1590a systematic error.  Since the orbit spends more time near apocenter
1591than near pericenter, this seems to lead an increase in energy and
1592decrease in angular momentum.  Perhaps this could be randomized with
1593an additional bisection at the end?
1594*/
1595// X
1596double QuarticStumpffSolverIterator(double mu, double alpha, double s, 
1597                double r0, double u, double deltat) {
1598        double c[5];
1599        double d0, d1, d2, d3, s2, ssav = s, tmp1, tmp2, tmp3;
1600        double delta1, delta2, delta3;
1601        int i;
1602        for ( i = 1; i <= 20; ++i ) {
1603                s2 = s*s;
1604                StumpffFunctions(c, alpha*s2);
1605                tmp1 = s*c[1];
1606                tmp2 = s2*c[2];
1607                tmp3 = -r0*alpha+mu;
1608                d0 = u*tmp2+mu*s2*s*c[3];
1609                d0 += r0*tmp1;
1610                d0 -= deltat; 
1611                d1 = u*tmp1+mu*tmp2; 
1612                d1 += r0*c[0];
1613                d2 = u*c[0]+tmp3*tmp1; 
1614                d3 = tmp3*c[0]-u*alpha*c[1]; 
1615                delta1 = -d0/d1;
1616                delta2 = -d0/(d1+delta1*d2/2);
1617                delta3 = -d0/(d1+delta1*d2/2+delta2*delta2*d3/6);
1618                ssav= s;
1619                s = s + delta3;
1620                if ( fabs(delta3/s) < StumpffQSTol ) return s; 
1621        }
1622        fprintf(FLOG, "%s"," QUARTIC STUMPFF SOLVER NOT CONVERGENT? \n"); 
1623        fprintf(FLOG, "    s = %.16e \n", s); 
1624        fprintf(FLOG, "sprev = %.16e \n", ssav); 
1625        fprintf(FLOG, "delta3/s = %.16e \n", delta3/s); 
1626        fprintf(FLOG, "Stumpff x = %.16e \n", alpha*s2); 
1627        fprintf(FLOG, "deltas = %.16e ", delta1); fprintf(FLOG, "%.16e ", delta2); 
1628        fprintf(FLOG, "%.16e \n", delta3);
1629        CheckStumpffErrorCount();
1630        fflush(FLOG);
1631        return s; 
1632}
1633
1634// X
1635double QuarticStumpffSolver(double mu, double alpha, double s, double r0, double u, double deltat) {
1636        /* initial guess */
1637        double tmp = deltat/r0;
1638        double s0 = tmp;
1639        double r0dot = u/r0;
1640        tmp = tmp*tmp;
1641        s0 -= 0.5*r0dot*tmp;
1642        /* third order in deltat */
1643        tmp = tmp*deltat/r0;
1644        s0 += (0.5*r0dot*r0dot-1/6*(mu-alpha*r0)/r0)*tmp;
1645        s = QuarticStumpffSolverIterator( mu, alpha, s0, r0, u, deltat);
1646        return s; 
1647}
1648
1649/* ===============================================================================
1650*/
1651// X
1652void IncfgAndDerivsStumpff(double *vals, double mu, double alpha, double s, 
1653                double r0, double u, double deltat) {
1654        double r, s2, tmp1, tmp2, tmp3;
1655        double c[5];
1656        s = QuarticStumpffSolver( mu, alpha, s, r0, u, deltat);
1657        s2 = s*s;
1658        StumpffFunctions(c, alpha*s2);
1659        //recompute c[1] from angular mom.
1660        tmp1 = s*c[1];
1661        tmp2 = s2*c[2];
1662        tmp3 = u*tmp1+mu*tmp2;
1663        /* r */
1664        r = tmp3+(r0*c[0]);
1665        vals[0] = r;
1666        tmp3 = -mu*tmp2;
1667        /* f */
1668        vals[1] = tmp3/r0;
1669        /* g, this is usually bad in terms of roundoff */
1670        vals[2] = r0*tmp1+u*tmp2;
1671        /* derivs */
1672        vals[3] = -(mu*tmp1)/(r0*r);
1673        vals[4] =  tmp3/r;
1674        return; 
1675}
1676
1677void fgAndDerivsStumpff(double *vals, double mu, double alpha, double s, 
1678                double r0, double u, double deltat) {
1679        IncfgAndDerivsStumpff(vals, mu, alpha, s, r0, u, deltat);
1680        vals[1] += 1.0;
1681        vals[4] += 1.0;
1682        return; }
1683
1684double KeplerStumpffStep(int dim, double mu, double *xv, double deltat) {
1685        double fgvals[5], tmp;
1686        double *x = xv;
1687        double *v = xv+dim;
1688        double v02 = 0.0, u = 0.0, r0 = 0.0;
1689        double r0dot, alpha;
1690        int i;
1691        for ( i = 1; i <= dim; ++i ) {
1692                v02 += v[i]*v[i];
1693                u += x[i]*v[i];
1694                r0 += x[i]*x[i]; }
1695        r0 = sqrt(r0);
1696        r0dot = u/r0;
1697/* a = 1.0/(2.0/r0 - v02/mu); alpha = mu/a; */
1698        alpha = 2.0*mu/r0 - v02;
1699        fgAndDerivsStumpff(fgvals, mu, alpha, 0.0, r0, u, deltat);
1700/* now update */
1701        for ( i = 1; i <= dim; ++i ) {
1702                tmp = x[i];
1703                x[i] = fgvals[1]*x[i]+fgvals[2]*v[i];
1704                v[i] = fgvals[3]*tmp+fgvals[4]*v[i]; }
1705        return 0.0; }
1706
1707
1708// X return only the increment
1709double IncrementKeplerStumpffStep(int dim, double mu, double *xv, double deltat) {
1710        double fg[5], tmp, tmp2;
1711        double *x = xv;
1712        double *v = xv+dim;
1713        double alpha, v02 = 0.0, u = 0.0, r02 = 0.0, r0;
1714        int i, ii;
1715        double tmp1, r, delta1;
1716        for ( i = 1; i <= dim; ++i ) {
1717                v02 += v[i]*v[i];
1718                u += x[i]*v[i];
1719                r02 += x[i]*x[i]; 
1720        }
1721        r0 = sqrt(r02);
1722        /* a = 1.0/(2.0/r0 - v02/mu); alpha = mu/a; */
1723        alpha = 2.0*mu/r0 - v02;
1724        IncfgAndDerivsStumpff(fg, mu, alpha, 0.0, r0, u, deltat);
1725
1726        // using fg'-gf' = 1 ( ' = dot) or total energy do not help much
1727               
1728        // now update
1729        for ( i = 1; i <= dim; ++i ) {
1730                tmp = x[i];
1731                x[i] = fg[1]*x[i]+fg[2]*v[i]; 
1732                v[i] = fg[3]*tmp+fg[4]*v[i]; 
1733        }
1734        return 0.0; 
1735}
1736
1737/* FOR TEST
1738double IncrementKeplerStumpffStep(int dim, double mu, double *xv, double deltat) {
1739        double tmp[7], tmp2[7];
1740        int i;
1741        for ( i = 1; i <= 2*dim; ++i ) { tmp[i] = tmp2[i] = xv[i]; }
1742IncrementKeplerStumpffStepSub(dim, mu, tmp, deltat);
1743        for ( i = 1; i <= 2*dim; ++i ) tmp2[i] += tmp[i];
1744IncrementKeplerStumpffStepSub(dim, mu, tmp2, -deltat);
1745        for ( i = 1; i <= 2*dim; ++i ) xv[i] = (tmp[i]-tmp2[i])/2.0;
1746        return 0.0; }
1747*/
1748
1749       
1750
1751
Note: See TracBrowser for help on using the repository browser.