void eraPvxpv(double a[2][3], double b[2][3], double axb[2][3]) /* ** - - - - - - - - - ** e r a P v x p v ** - - - - - - - - - ** ** Outer (=vector=cross) product of two pv-vectors. ** ** Given: ** a double[2][3] first pv-vector ** b double[2][3] second pv-vector ** ** Returned: ** axb double[2][3] a x b ** ** Notes: ** ** 1) If the position and velocity components of the two pv-vectors are ** ( ap, av ) and ( bp, bv ), the result, a x b, is the pair of ** vectors ( ap x bp, ap x bv + av x bp ). The two vectors are the ** cross-product of the two p-vectors and its derivative. ** ** 2) It is permissible to re-use the same array for any of the ** arguments. ** ** Called: ** eraCpv copy pv-vector ** eraPxp vector product of two p-vectors ** eraPpp p-vector plus p-vector ** ** Copyright (C) 2013-2016, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { double wa[2][3], wb[2][3], axbd[3], adxb[3]; /* Make copies of the inputs. */ eraCpv(a, wa); eraCpv(b, wb); /* a x b = position part of result. */ eraPxp(wa[0], wb[0], axb[0]); /* a x bdot + adot x b = velocity part of result. */ eraPxp(wa[0], wb[1], axbd); eraPxp(wa[1], wb[0], adxb); eraPpp(axbd, adxb, axb[1]); return; }
double eraSepp(double a[3], double b[3]) /* ** - - - - - - - - ** e r a S e p p ** - - - - - - - - ** ** Angular separation between two p-vectors. ** ** Given: ** a double[3] first p-vector (not necessarily unit length) ** b double[3] second p-vector (not necessarily unit length) ** ** Returned (function value): ** double angular separation (radians, always positive) ** ** Notes: ** ** 1) If either vector is null, a zero result is returned. ** ** 2) The angular separation is most simply formulated in terms of ** scalar product. However, this gives poor accuracy for angles ** near zero and pi. The present algorithm uses both cross product ** and dot product, to deliver full accuracy whatever the size of ** the angle. ** ** Called: ** eraPxp vector product of two p-vectors ** eraPm modulus of p-vector ** eraPdp scalar product of two p-vectors ** ** Copyright (C) 2013-2016, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { double axb[3], ss, cs, s; /* Sine of angle between the vectors, multiplied by the two moduli. */ eraPxp(a, b, axb); ss = eraPm(axb); /* Cosine of the angle, multiplied by the two moduli. */ cs = eraPdp(a, b); /* The angle. */ s = ((ss != 0.0) || (cs != 0.0)) ? atan2(ss, cs) : 0.0; return s; }
void eraH2fk5(double rh, double dh, double drh, double ddh, double pxh, double rvh, double *r5, double *d5, double *dr5, double *dd5, double *px5, double *rv5) /* ** - - - - - - - - - ** e r a H 2 f k 5 ** - - - - - - - - - ** ** Transform Hipparcos star data into the FK5 (J2000.0) system. ** ** Given (all Hipparcos, epoch J2000.0): ** rh double RA (radians) ** dh double Dec (radians) ** drh double proper motion in RA (dRA/dt, rad/Jyear) ** ddh double proper motion in Dec (dDec/dt, rad/Jyear) ** pxh double parallax (arcsec) ** rvh double radial velocity (km/s, positive = receding) ** ** Returned (all FK5, equinox J2000.0, epoch J2000.0): ** r5 double RA (radians) ** d5 double Dec (radians) ** dr5 double proper motion in RA (dRA/dt, rad/Jyear) ** dd5 double proper motion in Dec (dDec/dt, rad/Jyear) ** px5 double parallax (arcsec) ** rv5 double radial velocity (km/s, positive = receding) ** ** Notes: ** ** 1) This function transforms Hipparcos star positions and proper ** motions into FK5 J2000.0. ** ** 2) The proper motions in RA are dRA/dt rather than ** cos(Dec)*dRA/dt, and are per year rather than per century. ** ** 3) The FK5 to Hipparcos transformation is modeled as a pure ** rotation and spin; zonal errors in the FK5 catalog are not ** taken into account. ** ** 4) See also eraFk52h, eraFk5hz, eraHfk5z. ** ** Called: ** eraStarpv star catalog data to space motion pv-vector ** eraFk5hip FK5 to Hipparcos rotation and spin ** eraRv2m r-vector to r-matrix ** eraRxp product of r-matrix and p-vector ** eraTrxp product of transpose of r-matrix and p-vector ** eraPxp vector product of two p-vectors ** eraPmp p-vector minus p-vector ** eraPvstar space motion pv-vector to star catalog data ** ** Reference: ** ** F.Mignard & M.Froeschle, Astron. Astrophys. 354, 732-739 (2000). ** ** Copyright (C) 2013-2015, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { int i; double pvh[2][3], r5h[3][3], s5h[3], sh[3], wxp[3], vv[3], pv5[2][3]; /* Hipparcos barycentric position/velocity pv-vector (normalized). */ eraStarpv(rh, dh, drh, ddh, pxh, rvh, pvh); /* FK5 to Hipparcos orientation matrix and spin vector. */ eraFk5hip(r5h, s5h); /* Make spin units per day instead of per year. */ for ( i = 0; i < 3; s5h[i++] /= 365.25 ); /* Orient the spin into the Hipparcos system. */ eraRxp(r5h, s5h, sh); /* De-orient the Hipparcos position into the FK5 system. */ eraTrxp(r5h, pvh[0], pv5[0]); /* Apply spin to the position giving an extra space motion component. */ eraPxp(pvh[0], sh, wxp); /* Subtract this component from the Hipparcos space motion. */ eraPmp(pvh[1], wxp, vv); /* De-orient the Hipparcos space motion into the FK5 system. */ eraTrxp(r5h, vv, pv5[1]); /* FK5 pv-vector to spherical. */ eraPvstar(pv5, r5, d5, dr5, dd5, px5, rv5); return; }
double eraPap(double a[3], double b[3]) /* ** - - - - - - - ** e r a P a p ** - - - - - - - ** ** Position-angle from two p-vectors. ** ** Given: ** a double[3] direction of reference point ** b double[3] direction of point whose PA is required ** ** Returned (function value): ** double position angle of b with respect to a (radians) ** ** Notes: ** ** 1) The result is the position angle, in radians, of direction b with ** respect to direction a. It is in the range -pi to +pi. The ** sense is such that if b is a small distance "north" of a the ** position angle is approximately zero, and if b is a small ** distance "east" of a the position angle is approximately +pi/2. ** ** 2) The vectors a and b need not be of unit length. ** ** 3) Zero is returned if the two directions are the same or if either ** vector is null. ** ** 4) If vector a is at a pole, the result is ill-defined. ** ** Called: ** eraPn decompose p-vector into modulus and direction ** eraPm modulus of p-vector ** eraPxp vector product of two p-vectors ** eraPmp p-vector minus p-vector ** eraPdp scalar product of two p-vectors ** ** Copyright (C) 2013-2015, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { double am, au[3], bm, st, ct, xa, ya, za, eta[3], xi[3], a2b[3], pa; /* Modulus and direction of the a vector. */ eraPn(a, &am, au); /* Modulus of the b vector. */ bm = eraPm(b); /* Deal with the case of a null vector. */ if ((am == 0.0) || (bm == 0.0)) { st = 0.0; ct = 1.0; } else { /* The "north" axis tangential from a (arbitrary length). */ xa = a[0]; ya = a[1]; za = a[2]; eta[0] = -xa * za; eta[1] = -ya * za; eta[2] = xa*xa + ya*ya; /* The "east" axis tangential from a (same length). */ eraPxp(eta, au, xi); /* The vector from a to b. */ eraPmp(b, a, a2b); /* Resolve into components along the north and east axes. */ st = eraPdp(a2b, xi); ct = eraPdp(a2b, eta); /* Deal with degenerate cases. */ if ((st == 0.0) && (ct == 0.0)) ct = 1.0; } /* Position angle. */ pa = atan2(st, ct); return pa; }
void eraHfk5z(double rh, double dh, double date1, double date2, double *r5, double *d5, double *dr5, double *dd5) /* ** - - - - - - - - - ** e r a H f k 5 z ** - - - - - - - - - ** ** Transform a Hipparcos star position into FK5 J2000.0, assuming ** zero Hipparcos proper motion. ** ** Given: ** rh double Hipparcos RA (radians) ** dh double Hipparcos Dec (radians) ** date1,date2 double TDB date (Note 1) ** ** Returned (all FK5, equinox J2000.0, date date1+date2): ** r5 double RA (radians) ** d5 double Dec (radians) ** dr5 double FK5 RA proper motion (rad/year, Note 4) ** dd5 double Dec proper motion (rad/year, Note 4) ** ** Notes: ** ** 1) The TT date date1+date2 is a Julian Date, apportioned in any ** convenient way between the two arguments. For example, ** JD(TT)=2450123.7 could be expressed in any of these ways, ** among others: ** ** date1 date2 ** ** 2450123.7 0.0 (JD method) ** 2451545.0 -1421.3 (J2000 method) ** 2400000.5 50123.2 (MJD method) ** 2450123.5 0.2 (date & time method) ** ** The JD method is the most natural and convenient to use in ** cases where the loss of several decimal digits of resolution ** is acceptable. The J2000 method is best matched to the way ** the argument is handled internally and will deliver the ** optimum resolution. The MJD method and the date & time methods ** are both good compromises between resolution and convenience. ** ** 2) The proper motion in RA is dRA/dt rather than cos(Dec)*dRA/dt. ** ** 3) The FK5 to Hipparcos transformation is modeled as a pure rotation ** and spin; zonal errors in the FK5 catalogue are not taken into ** account. ** ** 4) It was the intention that Hipparcos should be a close ** approximation to an inertial frame, so that distant objects have ** zero proper motion; such objects have (in general) non-zero ** proper motion in FK5, and this function returns those fictitious ** proper motions. ** ** 5) The position returned by this function is in the FK5 J2000.0 ** reference system but at date date1+date2. ** ** 6) See also eraFk52h, eraH2fk5, eraFk5zhz. ** ** Called: ** eraS2c spherical coordinates to unit vector ** eraFk5hip FK5 to Hipparcos rotation and spin ** eraRxp product of r-matrix and p-vector ** eraSxp multiply p-vector by scalar ** eraRxr product of two r-matrices ** eraTrxp product of transpose of r-matrix and p-vector ** eraPxp vector product of two p-vectors ** eraPv2s pv-vector to spherical ** eraAnp normalize angle into range 0 to 2pi ** ** Reference: ** ** F.Mignard & M.Froeschle, 2000, Astron.Astrophys. 354, 732-739. ** ** Copyright (C) 2013, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { double t, ph[3], r5h[3][3], s5h[3], sh[3], vst[3], rst[3][3], r5ht[3][3], pv5e[2][3], vv[3], w, r, v; /* Time interval from fundamental epoch J2000.0 to given date (JY). */ t = ((date1 - DJ00) + date2) / DJY; /* Hipparcos barycentric position vector (normalized). */ eraS2c(rh, dh, ph); /* FK5 to Hipparcos orientation matrix and spin vector. */ eraFk5hip(r5h, s5h); /* Rotate the spin into the Hipparcos system. */ eraRxp(r5h, s5h, sh); /* Accumulated Hipparcos wrt FK5 spin over that interval. */ eraSxp(t, s5h, vst); /* Express the accumulated spin as a rotation matrix. */ eraRv2m(vst, rst); /* Rotation matrix: accumulated spin, then FK5 to Hipparcos. */ eraRxr(r5h, rst, r5ht); /* De-orient & de-spin the Hipparcos position into FK5 J2000.0. */ eraTrxp(r5ht, ph, pv5e[0]); /* Apply spin to the position giving a space motion. */ eraPxp(sh, ph, vv); /* De-orient & de-spin the Hipparcos space motion into FK5 J2000.0. */ eraTrxp(r5ht, vv, pv5e[1]); /* FK5 position/velocity pv-vector to spherical. */ eraPv2s(pv5e, &w, d5, &r, dr5, dd5, &v); *r5 = eraAnp(w); return; }
void eraLd(double bm, double p[3], double q[3], double e[3], double em, double dlim, double p1[3]) /* ** - - - - - - ** e r a L d ** - - - - - - ** ** Apply light deflection by a solar-system body, as part of ** transforming coordinate direction into natural direction. ** ** Given: ** bm double mass of the gravitating body (solar masses) ** p double[3] direction from observer to source (unit vector) ** q double[3] direction from body to source (unit vector) ** e double[3] direction from body to observer (unit vector) ** em double distance from body to observer (au) ** dlim double deflection limiter (Note 4) ** ** Returned: ** p1 double[3] observer to deflected source (unit vector) ** ** Notes: ** ** 1) The algorithm is based on Expr. (70) in Klioner (2003) and ** Expr. (7.63) in the Explanatory Supplement (Urban & Seidelmann ** 2013), with some rearrangement to minimize the effects of machine ** precision. ** ** 2) The mass parameter bm can, as required, be adjusted in order to ** allow for such effects as quadrupole field. ** ** 3) The barycentric position of the deflecting body should ideally ** correspond to the time of closest approach of the light ray to ** the body. ** ** 4) The deflection limiter parameter dlim is phi^2/2, where phi is ** the angular separation (in radians) between source and body at ** which limiting is applied. As phi shrinks below the chosen ** threshold, the deflection is artificially reduced, reaching zero ** for phi = 0. ** ** 5) The returned vector p1 is not normalized, but the consequential ** departure from unit magnitude is always negligible. ** ** 6) The arguments p and p1 can be the same array. ** ** 7) To accumulate total light deflection taking into account the ** contributions from several bodies, call the present function for ** each body in succession, in decreasing order of distance from the ** observer. ** ** 8) For efficiency, validation is omitted. The supplied vectors must ** be of unit magnitude, and the deflection limiter non-zero and ** positive. ** ** References: ** ** Urban, S. & Seidelmann, P. K. (eds), Explanatory Supplement to ** the Astronomical Almanac, 3rd ed., University Science Books ** (2013). ** ** Klioner, Sergei A., "A practical relativistic model for micro- ** arcsecond astrometry in space", Astr. J. 125, 1580-1597 (2003). ** ** Called: ** eraPdp scalar product of two p-vectors ** eraPxp vector product of two p-vectors ** ** Copyright (C) 2013-2014, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { int i; double qpe[3], qdqpe, w, eq[3], peq[3]; /* q . (q + e). */ for (i = 0; i < 3; i++) { qpe[i] = q[i] + e[i]; } qdqpe = eraPdp(q, qpe); /* 2 x G x bm / ( em x c^2 x ( q . (q + e) ) ). */ w = bm * ERFA_SRS / em / ERFA_GMAX(qdqpe,dlim); /* p x (e x q). */ eraPxp(e, q, eq); eraPxp(p, eq, peq); /* Apply the deflection. */ for (i = 0; i < 3; i++) { p1[i] = p[i] + w*peq[i]; } /* Finished. */ }
void eraLtecm(double epj, double rm[3][3]) /* ** - - - - - - - - - ** e r a L t e c m ** - - - - - - - - - ** ** ICRS equatorial to ecliptic rotation matrix, long-term. ** ** Given: ** epj double Julian epoch (TT) ** ** Returned: ** rm double[3][3] ICRS to ecliptic rotation matrix ** ** Notes: ** ** 1) The matrix is in the sense ** ** E_ep = rm x P_ICRS, ** ** where P_ICRS is a vector with respect to ICRS right ascension ** and declination axes and E_ep is the same vector with respect to ** the (inertial) ecliptic and equinox of epoch epj. ** ** 2) P_ICRS is a free vector, merely a direction, typically of unit ** magnitude, and not bound to any particular spatial origin, such ** as the Earth, Sun or SSB. No assumptions are made about whether ** it represents starlight and embodies astrometric effects such as ** parallax or aberration. The transformation is approximately that ** between mean J2000.0 right ascension and declination and ecliptic ** longitude and latitude, with only frame bias (always less than ** 25 mas) to disturb this classical picture. ** ** 3) The Vondrak et al. (2011, 2012) 400 millennia precession model ** agrees with the IAU 2006 precession at J2000.0 and stays within ** 100 microarcseconds during the 20th and 21st centuries. It is ** accurate to a few arcseconds throughout the historical period, ** worsening to a few tenths of a degree at the end of the ** +/- 200,000 year time span. ** ** Called: ** eraLtpequ equator pole, long term ** eraLtpecl ecliptic pole, long term ** eraPxp vector product ** eraPn normalize vector ** ** References: ** ** Vondrak, J., Capitaine, N. and Wallace, P., 2011, New precession ** expressions, valid for long time intervals, Astron.Astrophys. 534, ** A22 ** ** Vondrak, J., Capitaine, N. and Wallace, P., 2012, New precession ** expressions, valid for long time intervals (Corrigendum), ** Astron.Astrophys. 541, C1 ** ** Copyright (C) 2013-2017, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { /* Frame bias (IERS Conventions 2010, Eqs. 5.21 and 5.33) */ const double dx = -0.016617 * ERFA_DAS2R, de = -0.0068192 * ERFA_DAS2R, dr = -0.0146 * ERFA_DAS2R; double p[3], z[3], w[3], s, x[3], y[3]; /* Equator pole. */ eraLtpequ(epj, p); /* Ecliptic pole (bottom row of equatorial to ecliptic matrix). */ eraLtpecl(epj, z); /* Equinox (top row of matrix). */ eraPxp(p, z, w); eraPn(w, &s, x); /* Middle row of matrix. */ eraPxp(z, x, y); /* Combine with frame bias. */ rm[0][0] = x[0] - x[1]*dr + x[2]*dx; rm[0][1] = x[0]*dr + x[1] + x[2]*de; rm[0][2] = - x[0]*dx - x[1]*de + x[2]; rm[1][0] = y[0] - y[1]*dr + y[2]*dx; rm[1][1] = y[0]*dr + y[1] + y[2]*de; rm[1][2] = - y[0]*dx - y[1]*de + y[2]; rm[2][0] = z[0] - z[1]*dr + z[2]*dx; rm[2][1] = z[0]*dr + z[1] + z[2]*de; rm[2][2] = - z[0]*dx - z[1]*de + z[2]; }
void palDvxv ( double va[3], double vb[3], double vc[3] ) { eraPxp( va, vb, vc ); }