void eraAtciqz(double rc, double dc, eraASTROM *astrom, double *ri, double *di) /* ** - - - - - - - - - - ** e r a A t c i q z ** - - - - - - - - - - ** ** Quick ICRS to CIRS transformation, given precomputed star- ** independent astrometry parameters, and assuming zero parallax and ** proper motion. ** ** Use of this function is appropriate when efficiency is important and ** where many star positions are to be transformed for one date. The ** star-independent parameters can be obtained by calling one of the ** functions eraApci[13], eraApcg[13], eraApco[13] or eraApcs[13]. ** ** The corresponding function for the case of non-zero parallax and ** proper motion is eraAtciq. ** ** Given: ** rc,dc double ICRS astrometric RA,Dec (radians) ** astrom eraASTROM* star-independent astrometry parameters: ** pmt double PM time interval (SSB, Julian years) ** eb double[3] SSB to observer (vector, au) ** eh double[3] Sun to observer (unit vector) ** em double distance from Sun to observer (au) ** v double[3] barycentric observer velocity (vector, c) ** bm1 double sqrt(1-|v|^2): reciprocal of Lorenz factor ** bpn double[3][3] bias-precession-nutation matrix ** along double longitude + s' (radians) ** xpl double polar motion xp wrt local meridian (radians) ** ypl double polar motion yp wrt local meridian (radians) ** sphi double sine of geodetic latitude ** cphi double cosine of geodetic latitude ** diurab double magnitude of diurnal aberration vector ** eral double "local" Earth rotation angle (radians) ** refa double refraction constant A (radians) ** refb double refraction constant B (radians) ** ** Returned: ** ri,di double CIRS RA,Dec (radians) ** ** Note: ** ** All the vectors are with respect to BCRS axes. ** ** 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: ** eraS2c spherical coordinates to unit vector ** eraLdsun light deflection due to Sun ** eraAb stellar aberration ** eraRxp product of r-matrix and p-vector ** eraC2s p-vector to spherical ** eraAnp normalize angle into range +/- pi ** ** Copyright (C) 2013-2015, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { double pco[3], pnat[3], ppr[3], pi[3], w; /* BCRS coordinate direction (unit vector). */ eraS2c(rc, dc, pco); /* Light deflection by the Sun, giving BCRS natural direction. */ eraLdsun(pco, astrom->eh, astrom->em, pnat); /* Aberration, giving GCRS proper direction. */ eraAb(pnat, astrom->v, astrom->em, astrom->bm1, ppr); /* Bias-precession-nutation, giving CIRS proper direction. */ eraRxp(astrom->bpn, ppr, pi); /* CIRS RA,Dec. */ eraC2s(pi, &w, di); *ri = eraAnp(w); /* Finished. */ }
void eraAticq(double ri, double di, eraASTROM *astrom, double *rc, double *dc) /* ** - - - - - - - - - ** e r a A t i c q ** - - - - - - - - - ** ** Quick CIRS RA,Dec to ICRS astrometric place, given the star- ** independent astrometry parameters. ** ** Use of this function is appropriate when efficiency is important and ** where many star positions are all to be transformed for one date. ** The star-independent astrometry parameters can be obtained by ** calling one of the functions eraApci[13], eraApcg[13], eraApco[13] ** or eraApcs[13]. ** ** Given: ** ri,di double CIRS RA,Dec (radians) ** astrom eraASTROM* star-independent astrometry parameters: ** pmt double PM time interval (SSB, Julian years) ** eb double[3] SSB to observer (vector, au) ** eh double[3] Sun to observer (unit vector) ** em double distance from Sun to observer (au) ** v double[3] barycentric observer velocity (vector, c) ** bm1 double sqrt(1-|v|^2): reciprocal of Lorenz factor ** bpn double[3][3] bias-precession-nutation matrix ** along double longitude + s' (radians) ** xpl double polar motion xp wrt local meridian (radians) ** ypl double polar motion yp wrt local meridian (radians) ** sphi double sine of geodetic latitude ** cphi double cosine of geodetic latitude ** diurab double magnitude of diurnal aberration vector ** eral double "local" Earth rotation angle (radians) ** refa double refraction constant A (radians) ** refb double refraction constant B (radians) ** ** Returned: ** rc,dc double ICRS astrometric RA,Dec (radians) ** ** Notes: ** ** 1) Only the Sun is taken into account in the light deflection ** correction. ** ** 2) Iterative techniques are used for the aberration and light ** deflection corrections so that the functions eraAtic13 (or ** eraAticq) and eraAtci13 (or eraAtciq) are accurate inverses; ** even at the edge of the Sun's disk the discrepancy is only about ** 1 nanoarcsecond. ** ** Called: ** eraS2c spherical coordinates to unit vector ** eraTrxp product of transpose of r-matrix and p-vector ** eraZp zero p-vector ** eraAb stellar aberration ** eraLdsun light deflection by the Sun ** eraC2s p-vector to spherical ** eraAnp normalize angle into range +/- pi ** ** Copyright (C) 2013-2015, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { int j, i; double pi[3], ppr[3], pnat[3], pco[3], w, d[3], before[3], r2, r, after[3]; /* CIRS RA,Dec to Cartesian. */ eraS2c(ri, di, pi); /* Bias-precession-nutation, giving GCRS proper direction. */ eraTrxp(astrom->bpn, pi, ppr); /* Aberration, giving GCRS natural direction. */ eraZp(d); for (j = 0; j < 2; j++) { r2 = 0.0; for (i = 0; i < 3; i++) { w = ppr[i] - d[i]; before[i] = w; r2 += w*w; } r = sqrt(r2); for (i = 0; i < 3; i++) { before[i] /= r; } eraAb(before, astrom->v, astrom->em, astrom->bm1, after); r2 = 0.0; for (i = 0; i < 3; i++) { d[i] = after[i] - before[i]; w = ppr[i] - d[i]; pnat[i] = w; r2 += w*w; } r = sqrt(r2); for (i = 0; i < 3; i++) { pnat[i] /= r; } } /* Light deflection by the Sun, giving BCRS coordinate direction. */ eraZp(d); for (j = 0; j < 5; j++) { r2 = 0.0; for (i = 0; i < 3; i++) { w = pnat[i] - d[i]; before[i] = w; r2 += w*w; } r = sqrt(r2); for (i = 0; i < 3; i++) { before[i] /= r; } eraLdsun(before, astrom->eh, astrom->em, after); r2 = 0.0; for (i = 0; i < 3; i++) { d[i] = after[i] - before[i]; w = pnat[i] - d[i]; pco[i] = w; r2 += w*w; } r = sqrt(r2); for (i = 0; i < 3; i++) { pco[i] /= r; } } /* ICRS astrometric RA,Dec. */ eraC2s(pco, &w, dc); *rc = eraAnp(w); /* Finished. */ }
void eraAtciq(double rc, double dc, double pr, double pd, double px, double rv, eraASTROM *astrom, double *ri, double *di) /* ** - - - - - - - - - ** e r a A t c i q ** - - - - - - - - - ** ** Quick ICRS, epoch J2000.0, to CIRS transformation, given precomputed ** star-independent astrometry parameters. ** ** Use of this function is appropriate when efficiency is important and ** where many star positions are to be transformed for one date. The ** star-independent parameters can be obtained by calling one of the ** functions eraApci[13], eraApcg[13], eraApco[13] or eraApcs[13]. ** ** If the parallax and proper motions are zero the eraAtciqz function ** can be used instead. ** ** Given: ** rc,dc double ICRS RA,Dec at J2000.0 (radians) ** pr double RA proper motion (radians/year; Note 3) ** pd double Dec proper motion (radians/year) ** px double parallax (arcsec) ** rv double radial velocity (km/s, +ve if receding) ** astrom eraASTROM* star-independent astrometry parameters: ** pmt double PM time interval (SSB, Julian years) ** eb double[3] SSB to observer (vector, au) ** eh double[3] Sun to observer (unit vector) ** em double distance from Sun to observer (au) ** v double[3] barycentric observer velocity (vector, c) ** bm1 double sqrt(1-|v|^2): reciprocal of Lorenz factor ** bpn double[3][3] bias-precession-nutation matrix ** along double longitude + s' (radians) ** xpl double polar motion xp wrt local meridian (radians) ** ypl double polar motion yp wrt local meridian (radians) ** sphi double sine of geodetic latitude ** cphi double cosine of geodetic latitude ** diurab double magnitude of diurnal aberration vector ** eral double "local" Earth rotation angle (radians) ** refa double refraction constant A (radians) ** refb double refraction constant B (radians) ** ** Returned: ** ri,di double CIRS RA,Dec (radians) ** ** Notes: ** ** 1) All the vectors are with respect to BCRS axes. ** ** 2) Star data for an epoch other than J2000.0 (for example from the ** Hipparcos catalog, which has an epoch of J1991.25) will require a ** preliminary call to eraPmsafe before use. ** ** 3) The proper motion in RA is dRA/dt rather than cos(Dec)*dRA/dt. ** ** Called: ** eraPmpx proper motion and parallax ** eraLdsun light deflection by the Sun ** eraAb stellar aberration ** eraRxp product of r-matrix and pv-vector ** eraC2s p-vector to spherical ** eraAnp normalize angle into range 0 to 2pi ** ** Copyright (C) 2013-2016, NumFOCUS Foundation. ** Derived, with permission, from the SOFA library. See notes at end of file. */ { double pco[3], pnat[3], ppr[3], pi[3], w; /* Proper motion and parallax, giving BCRS coordinate direction. */ eraPmpx(rc, dc, pr, pd, px, rv, astrom->pmt, astrom->eb, pco); /* Light deflection by the Sun, giving BCRS natural direction. */ eraLdsun(pco, astrom->eh, astrom->em, pnat); /* Aberration, giving GCRS proper direction. */ eraAb(pnat, astrom->v, astrom->em, astrom->bm1, ppr); /* Bias-precession-nutation, giving CIRS proper direction. */ eraRxp(astrom->bpn, ppr, pi); /* CIRS RA,Dec. */ eraC2s(pi, &w, di); *ri = eraAnp(w); /* Finished. */ }