PaGMO  1.1.5
propagateKEP.cpp
1 /*****************************************************************************
2  * Copyright (C) 2004-2015 The PaGMO development team, *
3  * Advanced Concepts Team (ACT), European Space Agency (ESA) *
4  * *
5  * https://github.com/esa/pagmo *
6  * *
7  * act@esa.int *
8  * *
9  * This program is free software; you can redistribute it and/or modify *
10  * it under the terms of the GNU General Public License as published by *
11  * the Free Software Foundation; either version 2 of the License, or *
12  * (at your option) any later version. *
13  * *
14  * This program is distributed in the hope that it will be useful, *
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
17  * GNU General Public License for more details. *
18  * *
19  * You should have received a copy of the GNU General Public License *
20  * along with this program; if not, write to the *
21  * Free Software Foundation, Inc., *
22  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
23  *****************************************************************************/
24 
25 #include <cmath>
26 
27 #include "Astro_Functions.h"
28 #include "propagateKEP.h"
29 
30 /*
31  Origin: MATLAB code programmed by Dario Izzo (ESA/ACT)
32 
33  C++ version by Tamas Vinko (ESA/ACT)
34 
35  Inputs:
36  r0: column vector for the non dimensional position
37  v0: column vector for the non dimensional velocity
38  t: non dimensional time
39 
40  Outputs:
41  r: column vector for the non dimensional position
42  v: column vector for the non dimensional velocity
43 
44  Comments: The function works in non dimensional units, it takes an
45  initial condition and it propagates it as in a kepler motion analytically.
46 */
47 
48 using namespace std;
49 
50 void propagateKEP(const double *r0_in, const double *v0_in, const double &t, const double &mu,
51  double *r, double *v)
52 {
53 
54 /*
55  The matrix DD will be almost always the unit matrix, except for orbits
56  with little inclination in which cases a rotation is performed so that
57  par2IC is always defined
58 */
59 
60  double DD[9] = {1, 0, 0,
61  0, 1, 0,
62  0, 0, 1};
63  double h[3];
64  double ih[3] = {0,0,0};
65  double temp1[3] = {0,0,0}, temp2[3] = {0,0,0};
66  double E[6];
67  double normh, M, M0;
68  double r0[3], v0[3];
69 
70  int i;
71 
72  for (i=0; i<3; i++)
73  {
74  r0[i] = r0_in[i];
75  v0[i] = v0_in[i];
76  }
77 
78  vett(r0, v0, h);
79 
80  normh=norm2(h);
81 
82  for (i=0; i<3; i++)
83  ih[i] = h[i]/normh;
84 
85  if (fabs(fabs(ih[2])-1.0) < 1e-3) // the abs is needed in cases in which the orbit is retrograde,
86  { // that would held ih=[0,0,-1]!!
87  DD[0] = 1; DD[1] = 0; DD[2] = 0;
88  DD[3] = 0; DD[4] = 0; DD[5] = 1;
89  DD[6] = 0; DD[7] = -1; DD[8] = 0;
90 
91  // Random rotation matrix that make the Euler angles well defined for the case
92  // For orbits with little inclination another ref. frame is used.
93 
94  for (int i=0; i<3; i++)
95  {
96  temp1[0] += DD[i]*r0[i];
97  temp1[1] += DD[i+3]*r0[i];
98  temp1[2] += DD[i+6]*r0[i];
99  temp2[0] += DD[i]*v0[i];
100  temp2[1] += DD[i+3]*v0[i];
101  temp2[2] += DD[i+6]*v0[i];
102  }
103  for (int i=0; i<3; i++)
104  {
105  r0[i] = temp1[i];
106  temp1[i] = 0.0;
107  v0[i] = temp2[i];
108  temp2[i] = 0.0;
109  }
110  // for practical reason we take the transpose of the matrix DD here (will be used at the end of the function)
111  DD[0] = 1; DD[1] = 0; DD[2] = 0;
112  DD[3] = 0; DD[4] = 0; DD[5] = -1;
113  DD[6] = 0; DD[7] = 1; DD[8] = 0;
114  }
115 
116  IC2par(r0, v0, mu, E);
117  if (E[1] < 1.0)
118  {
119  M0 = E[5] - E[1]*sin(E[5]);
120  M=M0+sqrt(mu/pow(E[0],3))*t;
121  }
122  else
123  {
124  M0 = E[1]*tan(E[5]) - log(tan(0.5*E[5] + M_PI_4));
125  M=M0+sqrt(mu/pow(-E[0],3))*t;
126  }
127 
128  E[5]=Mean2Eccentric(M, E[1]);
129  par2IC(E, mu, r, v);
130 
131 
132  for (int j=0; j<3; j++)
133  {
134  temp1[0] += DD[j]*r[j];
135  temp1[1] += DD[j+3]*r[j];
136  temp1[2] += DD[j+6]*r[j];
137  temp2[0] += DD[j]*v[j];
138  temp2[1] += DD[j+3]*v[j];
139  temp2[2] += DD[j+6]*v[j];
140  }
141  for (int i=0; i<3; i++)
142  {
143  r[i] = temp1[i];
144  v[i] = temp2[i];
145  }
146 
147  return;
148 }
149 
150 
151 
152 
153 /*
154  Origin: MATLAB code programmed by Dario Izzo (ESA/ACT)
155 
156  C++ version by Tamas Vinko (ESA/ACT) 12/09/2006
157 
158  Inputs:
159  r0: column vector for the position
160  v0: column vector for the velocity
161 
162  Outputs:
163  E: Column Vectors containing the six keplerian parameters,
164  (a,e,i,OM,om,Eccentric Anomaly (or Gudermannian whenever e>1))
165 
166  Comments: The parameters returned are, of course, referred to the same
167  ref. frame in which r0,v0 are given. Units have to be consistent, and
168  output angles are in radians
169  The algorithm used is quite common and can be found as an example in Bate,
170  Mueller and White. It goes singular for zero inclination
171 */
172 
173 void IC2par(const double *r0, const double *v0, const double &mu, double *E)
174 {
175  double k[3];
176  double h[3];
177  double Dum_Vec[3];
178  double n[3];
179  double evett[3];
180 
181  double p = 0.0;
182  double temp =0.0;
183  double R0, ni;
184  int i;
185 
186  vett(r0, v0, h);
187 
188  for (i=0; i<3; i++)
189  p += h[i]*h[i];
190 
191  p/=mu;
192 
193 
194  k[0] = 0; k[1] = 0; k[2] = 1;
195  vett(k, h, n);
196 
197 
198  for (i=0; i<3; i++)
199  temp += pow(n[i], 2);
200 
201  temp = sqrt(temp);
202 
203  for (i=0; i<3; i++)
204  n[i] /= temp;
205 
206  R0 = norm2(r0);
207 
208  vett(v0, h, Dum_Vec);
209 
210  for (i=0; i<3; i++)
211  evett[i] = Dum_Vec[i]/mu - r0[i]/R0;
212 
213  double e = 0.0;
214  for (i=0; i<3; i++)
215  e += pow(evett[i], 2);
216 
217  E[0] = p/(1-e);
218  E[1] = sqrt(e);
219  e = E[1];
220 
221  E[2] = acos(h[2]/norm2(h));
222 
223  temp = 0.0;
224  for (i=0; i<3; i++)
225  temp+=n[i]*evett[i];
226 
227  E[4] = acos(temp/e);
228 
229  if (evett[2] < 0) E[4] = 2*M_PI - E[4];
230 
231  E[3] = acos(n[0]);
232  if (n[1] < 0) E[3] = 2*M_PI-E[3];
233 
234  temp = 0.0;
235  for (i=0; i<3; i++)
236  temp+=evett[i]*r0[i];
237 
238  ni = acos(temp/e/R0); // danger, the argument could be improper.
239 
240  temp = 0.0;
241  for (i=0; i<3; i++)
242  temp+=r0[i]*v0[i];
243 
244  if (temp<0.0) ni = 2*M_PI - ni;
245 
246  if (e<1.0)
247  E[5] = 2.0*atan(sqrt((1-e)/(1+e))*tan(ni/2.0)); // algebraic kepler's equation
248  else
249  E[5] =2.0*atan(sqrt((e-1)/(e+1))*tan(ni/2.0)); // algebraic equivalent of kepler's equation in terms of the Gudermannian
250 }
251 
252 /*
253  Origin: MATLAB code programmed by Dario Izzo (ESA/ACT)
254 
255  C++ version by Tamas Vinko (ESA/ACT)
256 
257  Usage: [r0,v0] = IC2par(E,mu)
258 
259  Outputs:
260  r0: column vector for the position
261  v0: column vector for the velocity
262 
263  Inputs:
264  E: Column Vectors containing the six keplerian parameters,
265  (a (negative for hyperbolas),e,i,OM,om,Eccentric Anomaly)
266  mu: gravitational constant
267 
268  Comments: The parameters returned are, of course, referred to the same
269  ref. frame in which r0,v0 are given. a can be given either in kms or AUs,
270  but has to be consistent with mu.All the angles must be given in radians.
271  This function does work for hyperbolas as well.
272 */
273 
274 void par2IC(const double *E, const double &mu, double *r0, double *v0)
275 {
276  double a = E[0];
277  double e = E[1];
278  double i = E[2];
279  double omg = E[3];
280  double omp = E[4];
281  double EA = E[5];
282  double b, n, xper, yper, xdotper, ydotper;
283  double R[3][3];
284  double cosomg, cosomp, sinomg, sinomp, cosi, sini;
285  double dNdZeta;
286 
287  // Grandezze definite nel piano dell'orbita
288 
289  if (e<1.0)
290  {
291  b = a*sqrt(1-e*e);
292  n = sqrt(mu/(a*a*a));
293  xper=a*(cos(EA)-e);
294  yper=b*sin(EA);
295 
296  xdotper = -(a*n*sin(EA))/(1-e*cos(EA));
297  ydotper=(b*n*cos(EA))/(1-e*cos(EA));
298  }
299  else
300  {
301  b = -a*sqrt(e*e-1);
302  n = sqrt(-mu/(a*a*a));
303 
304  dNdZeta = e * (1+tan(EA)*tan(EA))-(0.5+0.5*pow(tan(0.5*EA + M_PI_4),2))/tan(0.5*EA+ M_PI_4);
305 
306  xper = a/cos(EA) - a*e;
307  yper = b*tan(EA);
308 
309  xdotper = a*tan(EA)/cos(EA)*n/dNdZeta;
310  ydotper = b/pow(cos(EA), 2)*n/dNdZeta;
311  }
312 
313  // Matrice di trasformazione da perifocale a ECI
314 
315  cosomg = cos(omg);
316  cosomp = cos(omp);
317  sinomg = sin(omg);
318  sinomp = sin(omp);
319  cosi = cos(i);
320  sini = sin(i);
321 
322 
323  R[0][0]=cosomg*cosomp-sinomg*sinomp*cosi;
324  R[0][1]=-cosomg*sinomp-sinomg*cosomp*cosi;
325  R[0][2]=sinomg*sini;
326  R[1][0]=sinomg*cosomp+cosomg*sinomp*cosi;
327  R[1][1]=-sinomg*sinomp+cosomg*cosomp*cosi;
328  R[1][2]=-cosomg*sini;
329  R[2][0]=sinomp*sini;
330  R[2][1]=cosomp*sini;
331  R[2][2]=cosi;
332 
333  // Posizione nel sistema inerziale
334 
335 
336  double temp[3] = {xper, yper, 0.0};
337  double temp2[3] = {xdotper, ydotper, 0};
338 
339  for (int j = 0; j<3; j++)
340  {
341  r0[j] = 0.0; v0[j] = 0.0;
342  for (int k = 0; k<3; k++)
343  {
344  r0[j]+=R[j][k]*temp[k];
345  v0[j]+=R[j][k]*temp2[k];
346  }
347  }
348  return;
349 }
STL namespace.