int wcs_pv2sip(const char* wcsinfn, int ext, const char* wcsoutfn, anbool scamp_head_file, double* xy, int Nxy, int imageW, int imageH, anbool forcetan) { qfits_header* hdr = NULL; double* radec = NULL; int rtn = -1; tan_t tanwcs; double x,y, px,py; double xyz[3]; double* xorig = NULL; double* yorig = NULL; double* rddist = NULL; int i, j; // 1 x y r x2 xy y2 x3 x2y xy2 y3 r3 x4 x3y x2y2 xy3 y4 // x5 x4y x3y2 x2y3 xy4 y5 r5 x6 x5y x4y2, x3y3 x2y4 xy5 y6 // x7 x6y x5y2 x4y3 x3y4 x2y5 xy6 y7 r7 int xp[] = { 0, 1, 0, 0, 2, 1, 0, 3, 2, 1, 0, 0, 4, 3, 2, 1, 0, 5, 4, 3, 2, 1, 5, 0, 6, 5, 4, 3, 2, 1, 0, 7, 6, 5, 4, 3, 2, 1, 0, 0}; int yp[] = { 0, 0, 1, 0, 0, 1, 2, 0, 1, 2, 3, 0, 0, 1, 2, 3, 4, 0, 1, 2, 3, 4, 0, 0, 0, 1, 2, 3, 4, 5, 6, 0, 1, 2, 3, 4, 5, 6, 7, 0}; int rp[] = { 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 7}; double xpows[8]; double ypows[8]; double rpows[8]; double pv1[40]; double pv2[40]; double r; if (scamp_head_file) { size_t sz = 0; char* txt; char* prefix; int np; int nt; unsigned char* txthdr; sl* lines; int i; txt = file_get_contents(wcsinfn, &sz, TRUE); if (!txt) { ERROR("Failed to read file %s", wcsinfn); goto bailout; } lines = sl_split(NULL, txt, "\n"); prefix = "SIMPLE = T / Standard FITS file " "BITPIX = 8 / ASCII or bytes array " "NAXIS = 0 / Minimal header " "EXTEND = T / There may be FITS ext " "WCSAXES = 2 / "; np = strlen(prefix); nt = np + FITS_LINESZ * sl_size(lines); txthdr = malloc(nt); memset(txthdr, ' ', np + FITS_LINESZ * sl_size(lines)); memcpy(txthdr, prefix, np); for (i=0; i<sl_size(lines); i++) memcpy(txthdr + np + i*FITS_LINESZ, sl_get(lines, i), strlen(sl_get(lines, i))); sl_free2(lines); hdr = qfits_header_read_hdr_string(txthdr, nt); free(txthdr); free(txt); } else { char* ct; hdr = anqfits_get_header2(wcsinfn, ext); ct = fits_get_dupstring(hdr, "CTYPE1"); if ((ct && streq(ct, "RA---TPV")) || forcetan) { // http://iraf.noao.edu/projects/ccdmosaic/tpv.html logmsg("Replacing CTYPE1 = %s header with RA---TAN\n", ct); fits_update_value(hdr, "CTYPE1", "RA---TAN"); } ct = fits_get_dupstring(hdr, "CTYPE2"); if ((ct && streq(ct, "DEC--TPV")) || forcetan) { logmsg("Replacing CTYPE2 = %s header with DEC--TAN\n", ct); fits_update_value(hdr, "CTYPE2", "DEC--TAN"); } } if (!hdr) { ERROR("Failed to read header: file %s, ext %i\n", wcsinfn, ext); goto bailout; } tan_read_header(hdr, &tanwcs); for (i=0; i<sizeof(pv1)/sizeof(double); i++) { char key[10]; sprintf(key, "PV1_%i", i); pv1[i] = qfits_header_getdouble(hdr, key, 0.0); sprintf(key, "PV2_%i", i); pv2[i] = qfits_header_getdouble(hdr, key, 0.0); } xorig = malloc(Nxy * sizeof(double)); yorig = malloc(Nxy * sizeof(double)); rddist = malloc(2 * Nxy * sizeof(double)); for (j=0; j<Nxy; j++) { xorig[j] = xy[2*j+0]; yorig[j] = xy[2*j+1]; tan_pixelxy2iwc(&tanwcs, xorig[j], yorig[j], &x, &y); r = sqrt(x*x + y*y); xpows[0] = ypows[0] = rpows[0] = 1.0; for (i=1; i<sizeof(xpows)/sizeof(double); i++) { xpows[i] = xpows[i-1]*x; ypows[i] = ypows[i-1]*y; rpows[i] = rpows[i-1]*r; } px = py = 0; for (i=0; i<sizeof(xp)/sizeof(int); i++) { px += pv1[i] * xpows[xp[i]] * ypows[yp[i]] * rpows[rp[i]]; py += pv2[i] * ypows[xp[i]] * xpows[yp[i]] * rpows[rp[i]]; } tan_iwc2xyzarr(&tanwcs, px, py, xyz); xyzarr2radecdeg(xyz, rddist+2*j, rddist+2*j+1); } // { starxy_t sxy; tweak_t* t; il* imgi; il* refi; int sip_order = 5; int sip_inv_order = 5; sxy.N = Nxy; sxy.x = xorig; sxy.y = yorig; imgi = il_new(256); refi = il_new(256); for (i=0; i<Nxy; i++) { il_append(imgi, i); il_append(refi, i); } t = tweak_new(); t->sip->a_order = t->sip->b_order = sip_order; t->sip->ap_order = t->sip->bp_order = sip_inv_order; tweak_push_wcs_tan(t, &tanwcs); tweak_push_ref_ad_array(t, rddist, Nxy); tweak_push_image_xy(t, &sxy); tweak_push_correspondence_indices(t, imgi, refi, NULL, NULL); tweak_go_to(t, TWEAK_HAS_LINEAR_CD); if (imageW) t->sip->wcstan.imagew = imageW; if (imageH) t->sip->wcstan.imageh = imageH; sip_write_to_file(t->sip, wcsoutfn); tweak_free(t); } rtn = 0; bailout: free(xorig); free(yorig); free(rddist); qfits_header_destroy(hdr); free(radec); return rtn; }
sip_t* wcs_pv2sip_header(qfits_header* hdr, double* xy, int Nxy, double stepsize, double xlo, double xhi, double ylo, double yhi, int imageW, int imageH, int order, anbool forcetan, int doshift) { double* radec = NULL; int rtn = -1; tan_t tanwcs; double x,y, px,py; double* rddist = NULL; int i, j; int nx, ny; double xstep, ystep; sip_t* sip = NULL; /** From http://iraf.noao.edu/projects/mosaic/tpv.html p = PV1_ xi' = p0 + p1 * xi + p2 * eta + p3 * r + p4 * xi^2 + p5 * xi * eta + p6 * eta^2 + p7 * xi^3 + p8 * xi^2 * eta + p9 * xi * eta^2 + p10 * eta^3 + p11 * r^3 + p12 * xi^4 + p13 * xi^3 * eta + p14 * xi^2 * eta^2 + p15 * xi * eta^3 + p16 * eta^4 + p17 * xi^5 + p18 * xi^4 * eta + p19 * xi^3 * eta^2 + p20 * xi^2 * eta^3 + p21 * xi * eta^4 + p22 * eta^5 + p23 * r^5 + p24 * xi^6 + p25 * xi^5 * eta + p26 * xi^4 * eta^2 + p27 * xi^3 * eta^3 + p28 * xi^2 * eta^4 + p29 * xi * eta^5 + p30 * eta^6 p31 * xi^7 + p32 * xi^6 * eta + p33 * xi^5 * eta^2 + p34 * xi^4 * eta^3 + p35 * xi^3 * eta^4 + p36 * xi^2 * eta^5 + p37 * xi * eta^6 + p38 * eta^7 + p39 * r^7 p = PV2_ eta' = p0 + p1 * eta + p2 * xi + p3 * r + p4 * eta^2 + p5 * eta * xi + p6 * xi^2 + p7 * eta^3 + p8 * eta^2 * xi + p9 * eta * xi^2 + p10 * xi^3 + p11 * r^3 + p12 * eta^4 + p13 * eta^3 * xi + p14 * eta^2 * xi^2 + p15 * eta * xi^3 + p16 * xi^4 + p17 * eta^5 + p18 * eta^4 * xi + p19 * eta^3 * xi^2 + p20 * eta^2 * xi^3 + p21 * eta * xi^4 + p22 * xi^5 + p23 * r^5 + p24 * eta^6 + p25 * eta^5 * xi + p26 * eta^4 * xi^2 + p27 * eta^3 * xi^3 + p28 * eta^2 * xi^4 + p29 * eta * xi^5 + p30 * xi^6 p31 * eta^7 + p32 * eta^6 * xi + p33 * eta^5 * xi^2 + p34 * eta^4 * xi^3 + p35 * eta^3 * xi^4 + p36 * eta^2 * xi^5 + p37 * eta * xi^6 + p38 * xi^7 + p39 * r^7 Note the "cross-over" -- the xi' powers are in terms of xi,eta while the eta' powers are in terms of eta,xi. */ // 1 x y r x2 xy y2 x3 x2y xy2 y3 r3 x4 x3y x2y2 xy3 y4 // x5 x4y x3y2 x2y3 xy4 y5 r5 x6 x5y x4y2, x3y3 x2y4 xy5 y6 // x7 x6y x5y2 x4y3 x3y4 x2y5 xy6 y7 r7 int xp[] = { 0, 1, 0, 0, 2, 1, 0, 3, 2, 1, 0, 0, 4, 3, 2, 1, 0, 5, 4, 3, 2, 1, 0, 0, 6, 5, 4, 3, 2, 1, 0, 7, 6, 5, 4, 3, 2, 1, 0, 0}; int yp[] = { 0, 0, 1, 0, 0, 1, 2, 0, 1, 2, 3, 0, 0, 1, 2, 3, 4, 0, 1, 2, 3, 4, 5, 0, 0, 1, 2, 3, 4, 5, 6, 0, 1, 2, 3, 4, 5, 6, 7, 0}; int rp[] = { 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 3, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 5, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 7}; double xpows[8]; double ypows[8]; double rpows[8]; double pv1[40]; double pv2[40]; double r; char* ct; ct = fits_get_dupstring(hdr, "CTYPE1"); if ((ct && streq(ct, "RA---TPV")) || forcetan) { // http://iraf.noao.edu/projects/ccdmosaic/tpv.html logmsg("Replacing CTYPE1 = %s header with RA---TAN\n", ct); fits_update_value(hdr, "CTYPE1", "RA---TAN"); } ct = fits_get_dupstring(hdr, "CTYPE2"); if ((ct && streq(ct, "DEC--TPV")) || forcetan) { logmsg("Replacing CTYPE2 = %s header with DEC--TAN\n", ct); fits_update_value(hdr, "CTYPE2", "DEC--TAN"); } tan_read_header(hdr, &tanwcs); if (log_get_level() >= LOG_VERB) { printf("Read TAN header:\n"); tan_print(&tanwcs); } if (imageW && (imageW != tanwcs.imagew)) { logmsg("Overriding image width %f with user-specified %i\n", tanwcs.imagew, imageW); tanwcs.imagew = imageW; } if (imageH && (imageH != tanwcs.imageh)) { logmsg("Overriding image height %f with user-specified %i\n", tanwcs.imageh, imageH); tanwcs.imageh = imageH; } for (i=0; i<sizeof(pv1)/sizeof(double); i++) { char key[10]; double defaultval; if (i == 1) { defaultval = 1.0; } else { defaultval = 0.0; } sprintf(key, "PV1_%i", i); pv1[i] = qfits_header_getdouble(hdr, key, defaultval); sprintf(key, "PV2_%i", i); pv2[i] = qfits_header_getdouble(hdr, key, defaultval); } // choose grid for evaluating TAN-PV WCS if (xlo == 0 && xhi == 0) { xlo = 1.; xhi = tanwcs.imagew; } if (ylo == 0 && yhi == 0) { ylo = 1.; yhi = tanwcs.imageh; } assert(xhi >= xlo); assert(yhi >= ylo); if (stepsize == 0) stepsize = 100.; nx = MAX(2, round((xhi - xlo)/stepsize)); ny = MAX(2, round((yhi - ylo)/stepsize)); xstep = (xhi - xlo) / (double)(nx - 1); ystep = (yhi - ylo) / (double)(ny - 1); logverb("Stepping from x = %g to %g, steps of %g\n", xlo, xhi, xstep); logverb("Stepping from y = %g to %g, steps of %g\n", ylo, yhi, ystep); Nxy = nx * ny; if (xy == NULL) { int k = 0; xy = malloc(Nxy * 2 * sizeof(double)); for (i=0; i<ny; i++) { y = ylo + i*ystep; for (j=0; j<nx; j++) { x = xlo + j*xstep; //if (i == 0) //printf("x=%f\n", x); xy[k] = x; k++; xy[k] = y; k++; } //printf("y=%f\n", y); } assert(k == (Nxy*2)); } // distorted RA,Dec rddist = malloc(2 * Nxy * sizeof(double)); for (j=0; j<Nxy; j++) { double ix = xy[2*j+0]; double iy = xy[2*j+1]; tan_pixelxy2iwc(&tanwcs, ix, iy, &x, &y); // "x,y" here are most commonly known as "xi, eta". r = sqrt(x*x + y*y); // compute powers of x,y xpows[0] = ypows[0] = rpows[0] = 1.0; for (i=1; i<sizeof(xpows)/sizeof(double); i++) { xpows[i] = xpows[i-1]*x; ypows[i] = ypows[i-1]*y; rpows[i] = rpows[i-1]*r; } px = py = 0; for (i=0; i<sizeof(xp)/sizeof(int); i++) { px += pv1[i] * xpows[xp[i]] * ypows[yp[i]] * rpows[rp[i]]; // here's the "cross-over" mentioned above py += pv2[i] * ypows[xp[i]] * xpows[yp[i]] * rpows[rp[i]]; } // Note that the PV terms *include* a linear term, so no need // to re-add x,y to px,py. tan_iwc2radec(&tanwcs, px, py, rddist + 2*j, rddist + 2*j + 1); } sip = malloc(sizeof(sip_t)); assert(sip); { double* starxyz; starxyz = malloc(3 * Nxy * sizeof(double)); for (i=0; i<Nxy; i++) radecdegarr2xyzarr(rddist + i*2, starxyz + i*3); memset(sip, 0, sizeof(sip_t)); rtn = fit_sip_coefficients(starxyz, xy, NULL, Nxy, &tanwcs, order, order, sip); assert(rtn == 0); if (log_get_level() >= LOG_VERB) { printf("Fit SIP:\n"); sip_print(sip); } // FIXME? -- use xlo,xhi,ylo,yhi here?? Not clear. sip_compute_inverse_polynomials(sip, 0, 0, 0, 0, 0, 0); if (log_get_level() >= LOG_VERB) { printf("Fit SIP inverse polynomials:\n"); sip_print(sip); } free(starxyz); } free(rddist); free(radec); return sip; }