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