// Pixels to RA,Dec in degrees. void sip_pixelxy2radec(const sip_t* sip, double px, double py, double *ra, double *dec) { if (has_distortions(sip)) { double U, V; sip_distortion(sip, px, py, &U, &V); // Run a normal TAN conversion on the distorted pixel coords. tan_pixelxy2radec(&(sip->wcstan), U, V, ra, dec); } else // Run a normal TAN conversion tan_pixelxy2radec(&(sip->wcstan), px, py, ra, dec); }
// Pixels to RA,Dec in degrees. void tpv_pixelxy2radec(const tpv_t* tpv, double px, double py, double *ra, double *dec) { double U, V; tpv_distortion(tpv, px, py, &U, &V); // Run a normal TAN conversion on the distorted pixel coords. tan_pixelxy2radec(&(tpv->wcstan), U, V, ra, dec); }
int fit_sip_wcs_2(const double* starxyz, const double* fieldxy, const double* weights, int M, int sip_order, int inv_order, int W, int H, int crpix_center, double* crpix, int doshift, sip_t* sipout) { tan_t wcs; memset(&wcs, 0, sizeof(tan_t)); // Start with a TAN if (fit_tan_wcs(starxyz, fieldxy, M, &wcs, NULL)) { ERROR("Failed to fit for TAN WCS"); return -1; } if (crpix || crpix_center) { double cx,cy; double cr,cd; if (crpix) { cx = crpix[0]; cy = crpix[1]; } else { int i; if (W == 0) { for (i=0; i<M; i++) { W = MAX(W, (int)ceil(fieldxy[2*i + 0])); } } if (H == 0) { for (i=0; i<M; i++) { H = MAX(H, (int)ceil(fieldxy[2*i + 1])); } } cx = 1. + 0.5 * W; cy = 1. + 0.5 * H; } tan_pixelxy2radec(&wcs, cx, cy, &cr, &cd); wcs.crpix[0] = cx; wcs.crpix[1] = cy; wcs.crval[0] = cr; wcs.crval[1] = cd; } wcs.imagew = W; wcs.imageh = H; return fit_sip_wcs(starxyz, fieldxy, weights, M, &wcs, sip_order, inv_order, doshift, sipout); }
// Given a pixel offset (shift in image plane), adjust the WCS // CRVAL to the position given by CRPIX + offset. // Why not just // sip_pixelxy2radec(wcs, crpix0 +- xs, crpix1 +- ys, // wcs->wcstan.crval+0, wcs->wcstan.crval+1); // The answer is that "North" changes when you move the reference point. void wcs_shift(tan_t* wcs, double xs, double ys) { // UNITS: xs/ys in pixels // crvals in degrees, nx/nyref and theta in degrees double crpix0, crpix1, crval0; double theta, sintheta, costheta; double newcrval0, newcrval1; double newcd00,newcd01,newcd10,newcd11; // Save old vals crpix0 = wcs->crpix[0]; crpix1 = wcs->crpix[1]; crval0 = wcs->crval[0]; // compute the desired projection of the new tangent point by // shifting the projection of the current tangent point wcs->crpix[0] += xs; wcs->crpix[1] += ys; // now reproject the old crpix[xy] into shifted wcs tan_pixelxy2radec(wcs, crpix0, crpix1, &newcrval0, &newcrval1); // Restore crpix wcs->crpix[0] = crpix0; wcs->crpix[1] = crpix1; // RA,DEC coords of new tangent point wcs->crval[0] = newcrval0; wcs->crval[1] = newcrval1; theta = -deg2rad(newcrval0 - crval0); // deltaRA = new minus old RA; theta *= sin(deg2rad(newcrval1)); // multiply by the sin of the NEW Dec; at equator this correctly // evals to zero sintheta = sin(theta); costheta = cos(theta); // Fix the CD matrix since "northwards" has changed due to moving RA newcd00 = costheta * wcs->cd[0][0] - sintheta * wcs->cd[0][1]; newcd01 = sintheta * wcs->cd[0][0] + costheta * wcs->cd[0][1]; newcd10 = costheta * wcs->cd[1][0] - sintheta * wcs->cd[1][1]; newcd11 = sintheta * wcs->cd[1][0] + costheta * wcs->cd[1][1]; wcs->cd[0][0] = newcd00; wcs->cd[0][1] = newcd01; wcs->cd[1][0] = newcd10; wcs->cd[1][1] = newcd11; }
static int fit_tan_wcs_solve(const double* starxyz, const double* fieldxy, const double* weights, int N, const double* crpix, const tan_t* tanin, tan_t* tanout, double* p_scale) { int i, j, k; double field_cm[2] = {0, 0}; double cov[4] = {0, 0, 0, 0}; double R[4] = {0, 0, 0, 0}; double scale; // projected star coordinates double* p; // relative field coordinates double* f; double pcm[2] = {0, 0}; double w = 0; double totalw; gsl_matrix* A; gsl_matrix* U; gsl_matrix* V; gsl_vector* S; gsl_vector* work; gsl_matrix_view vcov; gsl_matrix_view vR; double crxyz[3]; double star_cm[3] = {0, 0, 0}; assert(((tanin != NULL) && (crpix != NULL)) || ((tanin == NULL) && (crpix == NULL))); if (tanin) { // default vals... memcpy(tanout, tanin, sizeof(tan_t)); } else { memset(tanout, 0, sizeof(tan_t)); } // -allocate and fill "p" and "f" arrays. ("projected" and "field") p = malloc(N * 2 * sizeof(double)); f = malloc(N * 2 * sizeof(double)); // -get field center-of-mass totalw = 0.0; for (i=0; i<N; i++) { w = (weights ? weights[i] : 1.0); field_cm[0] += w * fieldxy[i*2 + 0]; field_cm[1] += w * fieldxy[i*2 + 1]; totalw += w; } field_cm[0] /= totalw; field_cm[1] /= totalw; // Subtract it out. for (i=0; i<N; i++) { f[2*i+0] = fieldxy[2*i+0] - field_cm[0]; f[2*i+1] = fieldxy[2*i+1] - field_cm[1]; } if (tanin) { // Use original WCS to set the center of projection to the new crpix. tan_pixelxy2xyzarr(tanin, crpix[0], crpix[1], crxyz); for (i=0; i<N; i++) { Unused anbool ok; // -project the stars around crval ok = star_coords(starxyz + i*3, crxyz, TRUE, p + 2*i, p + 2*i + 1); assert(ok); } } else { // -get the star center-of-mass (this will become the tangent point CRVAL) for (i=0; i<N; i++) { w = (weights ? weights[i] : 1.0); star_cm[0] += w * starxyz[i*3 + 0]; star_cm[1] += w * starxyz[i*3 + 1]; star_cm[2] += w * starxyz[i*3 + 2]; } normalize_3(star_cm); // -project the stars around their center of mass for (i=0; i<N; i++) { Unused anbool ok; ok = star_coords(starxyz + i*3, star_cm, TRUE, p + 2*i, p + 2*i + 1); assert(ok); } } // -compute the center of mass of the projected stars and subtract it out. for (i=0; i<N; i++) { w = (weights ? weights[i] : 1.0); pcm[0] += w * p[2*i + 0]; pcm[1] += w * p[2*i + 1]; } pcm[0] /= totalw; pcm[1] /= totalw; for (i=0; i<N; i++) { p[2*i + 0] -= pcm[0]; p[2*i + 1] -= pcm[1]; } // -compute the covariance between field positions and projected // positions of the corresponding stars. for (i=0; i<N; i++) for (j=0; j<2; j++) for (k=0; k<2; k++) cov[j*2 + k] += p[i*2 + k] * f[i*2 + j]; for (i=0; i<4; i++) assert(isfinite(cov[i])); // -run SVD V = gsl_matrix_alloc(2, 2); S = gsl_vector_alloc(2); work = gsl_vector_alloc(2); vcov = gsl_matrix_view_array(cov, 2, 2); vR = gsl_matrix_view_array(R, 2, 2); A = &(vcov.matrix); // The Jacobi version doesn't always compute an orthonormal U if S has zeros. //gsl_linalg_SV_decomp_jacobi(A, V, S); gsl_linalg_SV_decomp(A, V, S, work); // the U result is written to A. U = A; gsl_vector_free(S); gsl_vector_free(work); // R = V U' gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, V, U, 0.0, &(vR.matrix)); gsl_matrix_free(V); for (i=0; i<4; i++) assert(isfinite(R[i])); // -compute scale: make the variances equal. { double pvar, fvar; pvar = fvar = 0.0; for (i=0; i<N; i++) { w = (weights ? weights[i] : 1.0); for (j=0; j<2; j++) { pvar += w * square(p[i*2 + j]); fvar += w * square(f[i*2 + j]); } } scale = sqrt(pvar / fvar); } // -compute WCS parameters. scale = rad2deg(scale); tanout->cd[0][0] = R[0] * scale; // CD1_1 tanout->cd[0][1] = R[1] * scale; // CD1_2 tanout->cd[1][0] = R[2] * scale; // CD2_1 tanout->cd[1][1] = R[3] * scale; // CD2_2 assert(isfinite(tanout->cd[0][0])); assert(isfinite(tanout->cd[0][1])); assert(isfinite(tanout->cd[1][0])); assert(isfinite(tanout->cd[1][1])); if (tanin) { // CRPIX is fixed. tanout->crpix[0] = crpix[0]; tanout->crpix[1] = crpix[1]; // Set CRVAL temporarily... tan_pixelxy2radec(tanin, crpix[0], crpix[1], tanout->crval+0, tanout->crval+1); // Shift CRVAL so that the center of the quad is in the right place. { double ix,iy; double dx,dy; double dxyz[3]; tan_pixelxy2iwc(tanout, field_cm[0], field_cm[1], &ix, &iy); dx = rad2deg(pcm[0]) - ix; dy = rad2deg(pcm[1]) - iy; tan_iwc2xyzarr(tanout, dx, dy, dxyz); xyzarr2radecdeg(dxyz, tanout->crval + 0, tanout->crval + 1); } } else { tanout->crpix[0] = field_cm[0]; tanout->crpix[1] = field_cm[1]; xyzarr2radecdegarr(star_cm, tanout->crval); // FIXME -- we ignore pcm. It should get added back in (after // multiplication by CD in the appropriate units) to either crval or // crpix. It's a very small correction probably of the same size // as the other approximations we're making. } if (p_scale) *p_scale = scale; free(p); free(f); return 0; }
void tan_pixelxy2radecarr(const tan_t* wcs_tan, double px, double py, double *radec) { tan_pixelxy2radec(wcs_tan, px, py, radec+0, radec+1); }
void test_wcslib_equals_tan(CuTest* tc) { #ifndef WCSLIB_EXISTS printf("\n\nWarning, WCSLIB support was not compiled in, so no WCSLIB functionality will be tested.\n\n\n"); return; #endif anwcs_t* anwcs = NULL; tan_t* tan; char* tmpfile = create_temp_file("test-anwcs-wcs", "/tmp"); double x, y, x2, y2, ra, dec, ra2, dec2; int ok, ok2; int i; // from http://live.astrometry.net/status.php?job=alpha-201004-29242410 // test-anwcs-1.wcs const char* wcsliteral = "SIMPLE = T / Standard FITS file BITPIX = 8 / ASCII or bytes array NAXIS = 0 / Minimal header EXTEND = T / There may be FITS ext CTYPE1 = 'RA---TAN' / TAN (gnomic) projection CTYPE2 = 'DEC--TAN' / TAN (gnomic) projection WCSAXES = 2 / no comment EQUINOX = 2000.0 / Equatorial coordinates definition (yr) LONPOLE = 180.0 / no comment LATPOLE = 0.0 / no comment CRVAL1 = 83.7131676182 / RA of reference point CRVAL2 = -5.10104333945 / DEC of reference point CRPIX1 = 221.593284607 / X reference pixel CRPIX2 = 169.655508041 / Y reference pixel CUNIT1 = 'deg ' / X pixel scale units CUNIT2 = 'deg ' / Y pixel scale units CD1_1 = 1.55258090814E-06 / Transformation matrix CD1_2 = 0.00081692280013 / no comment CD2_1 = -0.00081692280013 / no comment CD2_2 = 1.55258090814E-06 / no comment IMAGEW = 900 / Image width, in pixels. IMAGEH = 600 / Image height, in pixels. DATE = '2010-04-14T12:12:18' / Date this file was created. HISTORY Created by the Astrometry.net suite. HISTORY For more details, see http://astrometry.net . HISTORY Subversion URL HISTORY http://astrometry.net/svn/branches/astrometry/alpha/quads/ HISTORY Subversion revision 5409 HISTORY Subversion date 2007-10-09 13:49:13 -0400 (Tue, 09 Oct HISTORY 2007) HISTORY This WCS header was created by the program \"blind\". COMMENT -- blind solver parameters: -- COMMENT Index(0): /data1/INDEXES/200/index-219 COMMENT Index(1): /data1/INDEXES/200/index-218 COMMENT Index(2): /data1/INDEXES/200/index-217 COMMENT Index(3): /data1/INDEXES/200/index-216 COMMENT Index(4): /data1/INDEXES/200/index-215 COMMENT Index(5): /data1/INDEXES/200/index-214 COMMENT Index(6): /data1/INDEXES/200/index-213 COMMENT Index(7): /data1/INDEXES/200/index-212 COMMENT Index(8): /data1/INDEXES/200/index-211 COMMENT Index(9): /data1/INDEXES/200/index-210 COMMENT Index(10): /data1/INDEXES/200/index-209 COMMENT Index(11): /data1/INDEXES/200/index-208 COMMENT Index(12): /data1/INDEXES/200/index-207 COMMENT Index(13): /data1/INDEXES/200/index-206 COMMENT Index(14): /data1/INDEXES/200/index-205 COMMENT Index(15): /data1/INDEXES/200/index-204-00 COMMENT Index(16): /data1/INDEXES/200/index-204-01 COMMENT Index(17): /data1/INDEXES/200/index-204-02 COMMENT Index(18): /data1/INDEXES/200/index-204-03 COMMENT Index(19): /data1/INDEXES/200/index-204-04 COMMENT Index(20): /data1/INDEXES/200/index-204-05 COMMENT Index(21): /data1/INDEXES/200/index-204-06 COMMENT Index(22): /data1/INDEXES/200/index-204-07 COMMENT Index(23): /data1/INDEXES/200/index-204-08 COMMENT Index(24): /data1/INDEXES/200/index-204-09 COMMENT Index(25): /data1/INDEXES/200/index-204-10 COMMENT Index(26): /data1/INDEXES/200/index-204-11 COMMENT Index(27): /data1/INDEXES/200/index-203-00 COMMENT Index(28): /data1/INDEXES/200/index-203-01 COMMENT Index(29): /data1/INDEXES/200/index-203-02 COMMENT Index(30): /data1/INDEXES/200/index-203-03 COMMENT Index(31): /data1/INDEXES/200/index-203-04 COMMENT Index(32): /data1/INDEXES/200/index-203-05 COMMENT Index(33): /data1/INDEXES/200/index-203-06 COMMENT Index(34): /data1/INDEXES/200/index-203-07 COMMENT Index(35): /data1/INDEXES/200/index-203-08 COMMENT Index(36): /data1/INDEXES/200/index-203-09 COMMENT Index(37): /data1/INDEXES/200/index-203-10 COMMENT Index(38): /data1/INDEXES/200/index-203-11 COMMENT Index(39): /data1/INDEXES/200/index-202-00 COMMENT Index(40): /data1/INDEXES/200/index-202-01 COMMENT Index(41): /data1/INDEXES/200/index-202-02 COMMENT Index(42): /data1/INDEXES/200/index-202-03 COMMENT Index(43): /data1/INDEXES/200/index-202-04 COMMENT Index(44): /data1/INDEXES/200/index-202-05 COMMENT Index(45): /data1/INDEXES/200/index-202-06 COMMENT Index(46): /data1/INDEXES/200/index-202-07 COMMENT Index(47): /data1/INDEXES/200/index-202-08 COMMENT Index(48): /data1/INDEXES/200/index-202-09 COMMENT Index(49): /data1/INDEXES/200/index-202-10 COMMENT Index(50): /data1/INDEXES/200/index-202-11 COMMENT Index(51): /data1/INDEXES/200/index-201-00 COMMENT Index(52): /data1/INDEXES/200/index-201-01 COMMENT Index(53): /data1/INDEXES/200/index-201-02 COMMENT Index(54): /data1/INDEXES/200/index-201-03 COMMENT Index(55): /data1/INDEXES/200/index-201-04 COMMENT Index(56): /data1/INDEXES/200/index-201-05 COMMENT Index(57): /data1/INDEXES/200/index-201-06 COMMENT Index(58): /data1/INDEXES/200/index-201-07 COMMENT Index(59): /data1/INDEXES/200/index-201-08 COMMENT Index(60): /data1/INDEXES/200/index-201-09 COMMENT Index(61): /data1/INDEXES/200/index-201-10 COMMENT Index(62): /data1/INDEXES/200/index-201-11 COMMENT Index(63): /data1/INDEXES/200/index-200-00 COMMENT Index(64): /data1/INDEXES/200/index-200-01 COMMENT Index(65): /data1/INDEXES/200/index-200-02 COMMENT Index(66): /data1/INDEXES/200/index-200-03 COMMENT Index(67): /data1/INDEXES/200/index-200-04 COMMENT Index(68): /data1/INDEXES/200/index-200-05 COMMENT Index(69): /data1/INDEXES/200/index-200-06 COMMENT Index(70): /data1/INDEXES/200/index-200-07 COMMENT Index(71): /data1/INDEXES/200/index-200-08 COMMENT Index(72): /data1/INDEXES/200/index-200-09 COMMENT Index(73): /data1/INDEXES/200/index-200-10 COMMENT Index(74): /data1/INDEXES/200/index-200-11 COMMENT Field name: field.xy.fits COMMENT Field scale lower: 0.4 arcsec/pixel COMMENT Field scale upper: 720 arcsec/pixel COMMENT X col name: X COMMENT Y col name: Y COMMENT Start obj: 0 COMMENT End obj: 200 COMMENT Solved_in: solved COMMENT Solved_out: solved COMMENT Solvedserver: (null) COMMENT Parity: 2 COMMENT Codetol: 0.01 COMMENT Verify distance: 0 arcsec COMMENT Verify pixels: 1 pix COMMENT Maxquads: 0 COMMENT Maxmatches: 0 COMMENT Cpu limit: 0 s COMMENT Time limit: 0 s COMMENT Total time limit: 0 s COMMENT Total CPU limit: 600 s COMMENT Tweak: no COMMENT -- COMMENT -- properties of the matching quad: -- COMMENT quadno: 686636 COMMENT stars: 1095617,1095660,1095623,1095618 COMMENT field: 6,5,24,35 COMMENT code error: 0.00868071 COMMENT noverlap: 42 COMMENT nconflict: 1 COMMENT nfield: 88 COMMENT nindex: 139 COMMENT scale: 2.94093 arcsec/pix COMMENT parity: 1 COMMENT quads tried: 2166080 COMMENT quads matched: 2079562 COMMENT quads verified: 1747182 COMMENT objs tried: 0 COMMENT cpu time: 117.82 COMMENT -- AN_JOBID= 'alpha-201004-29242410' / Astrometry.net job ID END "; if (write_file(tmpfile, wcsliteral, strlen(wcsliteral))) { ERROR("Failed to write WCS to temp file %s", tmpfile); CuFail(tc, "failed to write WCS to temp file"); } tan = tan_read_header_file(tmpfile, NULL); CuAssertPtrNotNull(tc, tan); for (i=0; i<2; i++) { if (i == 0) { anwcs = anwcs_open_wcslib(tmpfile, 0); } else if (i == 1) { anwcs = anwcs_open_sip(tmpfile, 0); } CuAssertPtrNotNull(tc, anwcs); /* printf("ANWCS:\n"); anwcs_print(anwcs, stdout); printf("TAN:\n"); tan_print_to(tan, stdout); */ /* this wcs has: crval=(83.7132, -5.10104) crpix=(221.593, 169.656) CD = ( 1.5526e-06 0.00081692 ) ( -0.00081692 1.5526e-06 ) */ // check crval <-> crpix. x = tan->crpix[0]; y = tan->crpix[1]; ra = 1; ra2 = 2; dec = 3; dec2 = 4; tan_pixelxy2radec(tan, x, y, &ra, &dec); CuAssertDblEquals(tc, tan->crval[0], ra, 1e-6); CuAssertDblEquals(tc, tan->crval[1], dec, 1e-6); ok = anwcs_pixelxy2radec(anwcs, x, y, &ra2, &dec2); CuAssertIntEquals(tc, 0, ok); CuAssertDblEquals(tc, tan->crval[0], ra2, 1e-6); CuAssertDblEquals(tc, tan->crval[1], dec2, 1e-6); ra = tan->crval[0]; dec = tan->crval[1]; x = 1; x2 = 2; y = 3; y2 = 4; ok = tan_radec2pixelxy(tan, ra, dec, &x, &y); CuAssertIntEquals(tc, TRUE, ok); CuAssertDblEquals(tc, tan->crpix[0], x, 1e-6); CuAssertDblEquals(tc, tan->crpix[1], y, 1e-6); ok2 = anwcs_radec2pixelxy(anwcs, ra, dec, &x2, &y2); CuAssertIntEquals(tc, 0, ok2); CuAssertDblEquals(tc, tan->crpix[0], x2, 1e-6); CuAssertDblEquals(tc, tan->crpix[1], y2, 1e-6); // check pixel (0,0). x = y = 0.0; ra = 1; ra2 = 2; dec = 3; dec2 = 4; tan_pixelxy2radec(tan, x, y, &ra, &dec); ok = anwcs_pixelxy2radec(anwcs, x, y, &ra2, &dec2); CuAssertIntEquals(tc, 0, ok); CuAssertDblEquals(tc, ra, ra2, 1e-6); CuAssertDblEquals(tc, dec, dec2, 1e-6); // check RA,Dec (85, -4) ra = 85; dec = -4; x = 1; x2 = 2; y = 3; y2 = 4; ok = tan_radec2pixelxy(tan, ra, dec, &x, &y); CuAssertIntEquals(tc, TRUE, ok); ok2 = anwcs_radec2pixelxy(anwcs, ra, dec, &x2, &y2); CuAssertIntEquals(tc, 0, ok2); printf("x,y (%g,%g) vs (%g,%g)\n", x, y, x2, y2); CuAssertDblEquals(tc, x, x2, 1e-6); CuAssertDblEquals(tc, y, y2, 1e-6); anwcs_free(anwcs); } free(tan); }