예제 #1
0
int coadd_add_image(coadd_t* ca, const number* img,
					const number* weightimg,
					number weight, const anwcs_t* wcs) {
	int W, H;
	int i, j;
	int xlo,xhi,ylo,yhi;
	check_bounds_t cb;

	W = anwcs_imagew(wcs);
	H = anwcs_imageh(wcs);

	// if check_bounds:
	cb.xlo = W;
	cb.xhi = 0;
	cb.ylo = H;
	cb.yhi = 0;
	cb.wcs = ca->wcs;
	anwcs_walk_image_boundary(wcs, 50, check_bounds, &cb);
	xlo = MAX(0,     floor(cb.xlo));
	xhi = MIN(ca->W,  ceil(cb.xhi)+1);
	ylo = MAX(0,     floor(cb.ylo));
	yhi = MIN(ca->H,  ceil(cb.yhi)+1);
	logmsg("Image projects to output image region: [%i,%i), [%i,%i)\n", xlo, xhi, ylo, yhi);

	for (i=ylo; i<yhi; i++) {
		for (j=xlo; j<xhi; j++) {
			double ra, dec;
			double px, py;
			double wt;
			double val;

			// +1 for FITS
			if (anwcs_pixelxy2radec(ca->wcs, j+1, i+1, &ra, &dec)) {
				ERROR("Failed to project pixel (%i,%i) through output WCS\n", j, i);
				continue;
			}
			if (anwcs_radec2pixelxy(wcs, ra, dec, &px, &py)) {
				ERROR("Failed to project pixel (%i,%i) through input WCS\n", j, i);
				continue;
			}
			// -1 for FITS
			px -= 1;
			py -= 1;

			if (px < 0 || px >= W)
				continue;
			if (py < 0 || py >= H)
				continue;

			val = ca->resample_func(px, py, img, weightimg, W, H, &wt,
									ca->resample_token);
			ca->img[i*ca->W + j] += val * weight;
			ca->weight[i*ca->W + j] += wt * weight;
		}
		logverb("Row %i of %i\n", i+1, ca->H);
	}
	return 0;
}
예제 #2
0
void test_mercator_2(CuTest* tc) {
    double ra,dec;
    double ra2,dec2;
    anwcs_t* wcs = anwcs_create_mercator_2(180., 0., 128.5, 128.5,
                                           1., 256, 256, TRUE);
    printf("RA,Dec in corners:\n");
    anwcs_pixelxy2radec(wcs, 1., 1., &ra, &dec);
    CuAssertDblEquals(tc,  85.0, dec, 0.1);
    CuAssertDblEquals(tc, 360.0, ra,  1.0);
    printf("  1,1 -> %.1f, %.1f\n", ra, dec);
    anwcs_pixelxy2radec(wcs, 256., 256., &ra, &dec);
    printf("  W,H -> %.1f, %.1f\n", ra, dec);
    CuAssertDblEquals(tc, -85.0, dec, 0.1);
    CuAssertDblEquals(tc, 0.0, ra,  1.0);

    printf("Delta near the center:\n");
    anwcs_pixelxy2radec(wcs, 128., 128., &ra, &dec);
    anwcs_pixelxy2radec(wcs, 129., 128., &ra2, &dec2);
    printf("dx = 1 -> dRA, dDec %.1f, %.1f\n", ra2-ra, dec2-dec);
    anwcs_pixelxy2radec(wcs, 128., 129., &ra2, &dec2);
    printf("dy = 1 -> dRA, dDec %.1f, %.1f\n", ra2-ra, dec2-dec);
}
예제 #3
0
int fvec(const gsl_vector *x, void *params, gsl_vector *f)
{
  double ra,dec,xp,yp;
  double xi = gsl_vector_get(x,0);
  double yi = gsl_vector_get(x,1);
  anwcs_pixelxy2radec(wcs1, xi, yi, &ra, &dec);
  anwcs_radec2pixelxy(wcs2, ra, dec, &xp, &yp);
  xp = xp - xi;
  yp = yp - yi;
  gsl_vector_set(f,0,xp);
  gsl_vector_set(f,1,yp);
  return GSL_SUCCESS;
}
예제 #4
0
int wcs_xy2rd(const char* wcsfn, int ext,
			  const char* xylsfn, const char* rdlsfn,
              const char* xcol, const char* ycol,
			  int forcetan,
			  int forcewcslib,
              il* fields) {
	rdlist_t* rdls = NULL;
	xylist_t* xyls = NULL;
	anwcs_t* wcs = NULL;
	int i;
    int rtn = -1;
    anbool alloced_fields = FALSE;

	// read WCS.
	if (forcewcslib) {
		wcs = anwcs_open_wcslib(wcsfn, ext);
	} else if (forcetan) {
		wcs = anwcs_open_tan(wcsfn, ext);
	} else {
		wcs = anwcs_open(wcsfn, ext);
	}
	if (!wcs) {
		ERROR("Failed to read WCS file \"%s\", extension %i", wcsfn, ext);
		return -1;
	}

	// read XYLS.
	xyls = xylist_open(xylsfn);
	if (!xyls) {
		ERROR("Failed to read an xylist from file %s", xylsfn);
		goto bailout;
	}
    xylist_set_include_flux(xyls, FALSE);
    xylist_set_include_background(xyls, FALSE);
	if (xcol)
		xylist_set_xname(xyls, xcol);
	if (ycol)
		xylist_set_yname(xyls, ycol);

	// write RDLS.
	rdls = rdlist_open_for_writing(rdlsfn);
	if (!rdls) {
		ERROR("Failed to open file %s to write RDLS.\n", rdlsfn);
		goto bailout;
	}
	if (rdlist_write_primary_header(rdls)) {
		ERROR("Failed to write header to RDLS file %s.\n", rdlsfn);
		goto bailout;
	}

    if (!fields) {
        alloced_fields = TRUE;
        fields = il_new(16);
    }
	if (!il_size(fields)) {
		// add all fields.
		int NF = xylist_n_fields(xyls);
		for (i=1; i<=NF; i++)
			il_append(fields, i);
	}

	logverb("Processing %zu extensions...\n", il_size(fields));
	for (i=0; i<il_size(fields); i++) {
		int fieldind = il_get(fields, i);
        starxy_t xy;
        rd_t rd;
		int j;

        if (!xylist_read_field_num(xyls, fieldind, &xy)) {
			ERROR("Failed to read xyls file %s, field %i", xylsfn, fieldind);
			goto bailout;
        }

		if (rdlist_write_header(rdls)) {
			ERROR("Failed to write rdls field header to %s", rdlsfn);
			goto bailout;
		}

        rd_alloc_data(&rd, starxy_n(&xy));

		for (j=0; j<starxy_n(&xy); j++) {
            double x, y, ra, dec;
            x = starxy_getx(&xy, j);
            y = starxy_gety(&xy, j);
			anwcs_pixelxy2radec(wcs, x, y, &ra, &dec);
            rd_setra (&rd, j, ra);
            rd_setdec(&rd, j, dec);
		}

        if (rdlist_write_field(rdls, &rd)) {
            ERROR("Failed to write rdls field to %s", rdlsfn);
			goto bailout;
        }
        rd_free_data(&rd);
        starxy_free_data(&xy);

		if (rdlist_fix_header(rdls)) {
			ERROR("Failed to fix rdls field header for %s", rdlsfn);
			goto bailout;
		}

        rdlist_next_field(rdls);
	}

	if (rdlist_fix_primary_header(rdls) ||
		rdlist_close(rdls)) {
		ERROR("Failed to fix header of RDLS file %s", rdlsfn);
		goto bailout;
	}
	rdls = NULL;

	if (xylist_close(xyls)) {
		ERROR("Failed to close XYLS file %s", xylsfn);
		goto bailout;
	}
	xyls = NULL;
	rtn = 0;

 bailout:
    if (alloced_fields)
        il_free(fields);
    if (rdls)
        rdlist_close(rdls);
    if (xyls)
        xylist_close(xyls);
	if (wcs)
		anwcs_free(wcs);
    return rtn;
}
예제 #5
0
int plot_xy_plot(const char* command, cairo_t* cairo,
				 plot_args_t* pargs, void* baton) {
	plotxy_t* args = (plotxy_t*)baton;
	// Plot it!
	xylist_t* xyls;
	starxy_t myxy;
	starxy_t* xy = NULL;
	starxy_t* freexy = NULL;
	int Nxy;
	int i;
#if 0
	double t0;
#endif

	plotstuff_builtin_apply(cairo, pargs);

	if (args->fn && dl_size(args->xyvals)) {
		ERROR("Can only plot one of xylist filename and xy_vals");
		return -1;
	}
	if (!args->fn && !dl_size(args->xyvals)) {
		ERROR("Neither xylist filename nor xy_vals given!");
		return -1;
	}

	if (args->fn) {
#if 0
		t0 = timenow();
#endif
		// Open xylist.
		xyls = xylist_open(args->fn);
		if (!xyls) {
			ERROR("Failed to open xylist from file \"%s\"", args->fn);
			return -1;
		}
		// we don't care about FLUX and BACKGROUND columns.
		xylist_set_include_flux(xyls, FALSE);
		xylist_set_include_background(xyls, FALSE);
		if (args->xcol)
			xylist_set_xname(xyls, args->xcol);
		if (args->ycol)
			xylist_set_yname(xyls, args->ycol);

		// Find number of entries in xylist.
		xy = xylist_read_field_num(xyls, args->ext, NULL);
		freexy = xy;
		xylist_close(xyls);
		if (!xy) {
			ERROR("Failed to read FITS extension %i from file %s.\n", args->ext, args->fn);
			return -1;
		}
		Nxy = starxy_n(xy);
		// If N is specified, apply it as a max.
		if (args->nobjs)
			Nxy = MIN(Nxy, args->nobjs);
		//logmsg("%g s to read xylist\n", timenow()-t0);
	} else {
		assert(dl_size(args->xyvals));
		starxy_from_dl(&myxy, args->xyvals, FALSE, FALSE);
		xy = &myxy;
		Nxy = starxy_n(xy);
	}

	// Transform through WCSes.
	if (args->wcs) {
		double ra, dec, x, y;
		assert(pargs->wcs);
		/*
		 // check for any overlap.
		 double pralo,prahi,pdeclo,pdechi;
		 double ralo,rahi,declo,dechi;
		 anwcs_get_radec_bounds(pargs->wcs, 100, &pralo, &prahi, &pdeclo, &pdechi);
		 anwcs_get_radec_bounds(args->wcs, 100, &ralo, &rahi, &declo, &dechi);
		 if (
		 */
		for (i=0; i<Nxy; i++) {
			anwcs_pixelxy2radec(args->wcs,
								// I used to add 1 here
								starxy_getx(xy, i), starxy_gety(xy, i),
								&ra, &dec);
			if (!plotstuff_radec2xy(pargs, ra, dec, &x, &y))
				continue;
			logverb("  xy (%g,%g) -> RA,Dec (%g,%g) -> plot xy (%g,%g)\n",
					starxy_getx(xy,i), starxy_gety(xy,i), ra, dec, x, y);

			// add shift and scale...
			// FIXME -- not clear that we want to do this here...
			/*
			 starxy_setx(xy, i, args->scale * (x - args->xoff));
			 starxy_sety(xy, i, args->scale * (y - args->yoff));
			 starxy_setx(xy, i, x-1);
			 starxy_sety(xy, i, y-1);
			 */

			// Output coords: FITS -> 0-indexed image
			starxy_setx(xy, i, x-1);
			starxy_sety(xy, i, y-1);
		}
	} else {
		// Shift and scale xylist entries.
		if (args->xoff != 0.0 || args->yoff != 0.0) {
			for (i=0; i<Nxy; i++) {
				starxy_setx(xy, i, starxy_getx(xy, i) - args->xoff);
				starxy_sety(xy, i, starxy_gety(xy, i) - args->yoff);
			}
		}
		if (args->scale != 1.0) {
			for (i=0; i<Nxy; i++) {
				starxy_setx(xy, i, args->scale * starxy_getx(xy, i));
				starxy_sety(xy, i, args->scale * starxy_gety(xy, i));
			}
		}
	}

	// Plot markers.
#if 0
	t0 = timenow();
#endif
	for (i=args->firstobj; i<Nxy; i++) {
		double x = starxy_getx(xy, i);
		double y = starxy_gety(xy, i);
		if (plotstuff_marker_in_bounds(pargs, x, y))
			plotstuff_stack_marker(pargs, x, y);
	}
	plotstuff_plot_stack(pargs, cairo);
	//logmsg("%g s to plot xylist\n", timenow()-t0);

	starxy_free(freexy);
	return 0;
}
예제 #6
0
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);
}
예제 #7
0
int main(int argc, char** args) {
	int c;
	char* rdlsfn = NULL;
	char* wcsfn = NULL;
	char* xylsfn = NULL;
	char* xcol = NULL;
	char* ycol = NULL;
	anbool forcetan = FALSE;
	anbool forcewcslib = FALSE;
	anbool forcewcstools = FALSE;
	anbool printhms = FALSE;
	il* fields;
	int ext = 0;
	double x, y;
	int loglvl = LOG_MSG;

	x = y = HUGE_VAL;
	fields = il_new(16);

    while ((c = getopt(argc, args, OPTIONS)) != -1) {
        switch (c) {
		case 'v':
			loglvl++;
			break;
		case 's':
			printhms = TRUE;
			break;
		case 'L':
			forcewcslib = TRUE;
			break;
		case 'T':
			forcewcstools = TRUE;
			break;
		case 'x':
			x = atof(optarg);
			break;
		case 'y':
			y = atof(optarg);
			break;
		case 'e':
			ext = atoi(optarg);
			break;
        case 'h':
			print_help(args[0]);
			exit(0);
		case 't':
			forcetan = TRUE;
			break;
		case 'o':
			rdlsfn = optarg;
			break;
		case 'i':
			xylsfn = optarg;
			break;
		case 'w':
			wcsfn = optarg;
			break;
		case 'f':
			il_append(fields, atoi(optarg));
			break;
		case 'X':
			xcol = optarg;
			break;
		case 'Y':
			ycol = optarg;
			break;
		}
	}

	log_init(loglvl);

	if (optind != argc) {
		print_help(args[0]);
		exit(-1);
	}

	if (!(wcsfn && ((rdlsfn && xylsfn) || ((x != HUGE_VAL) && (y != HUGE_VAL))))) {
		print_help(args[0]);
		exit(-1);
	}

	if (!xylsfn) {
		double ra,dec;
		anwcs_t* wcs = NULL;

		// read WCS.
		if (forcewcslib) {
			wcs = anwcs_open_wcslib(wcsfn, ext);
		} else if (forcewcstools) {
			wcs = anwcs_open_wcstools(wcsfn, ext);
		} else if (forcetan) {
			wcs = anwcs_open_tan(wcsfn, ext);
		} else {
			wcs = anwcs_open(wcsfn, ext);
		}
		if (!wcs) {
			ERROR("Failed to read WCS file \"%s\", extension %i", wcsfn, ext);
			exit(-1);
		}

		logverb("Read WCS:\n");
		if (log_get_level() >= LOG_VERB) {
			anwcs_print(wcs, log_get_fid());
		}

		// convert immediately.
		anwcs_pixelxy2radec(wcs, x, y, &ra, &dec);
		printf("Pixel (%.10f, %.10f) -> RA,Dec (%.10f, %.10f)\n", x, y, ra, dec);
		if (printhms) {
			char str[32];
			ra2hmsstring(ra, str);
			printf("                       RA,Dec (%20s, ", str);
			dec2dmsstring(dec, str);
			printf("%20s)\n", str);
		}
		anwcs_free(wcs);
		exit(0);
	}

    if (wcs_xy2rd(wcsfn, ext, xylsfn, rdlsfn, 
                  xcol, ycol, forcetan, forcewcslib, fields)) {
        ERROR("wcs-xy2rd failed");
        exit(-1);
    }

	return 0;
}
예제 #8
0
int main(int argc, char** args) {
  int ext = 0,c;
  double ra,dec;
  double sol[2];
  const gsl_multiroot_fsolver_type *T;
  gsl_multiroot_fsolver *s;
  int status;
  size_t iter=0;
  const size_t n=2;
  gsl_multiroot_function f={&fvec,n,NULL};
  gsl_vector *x = gsl_vector_alloc(n);
  char *wcsfn1=NULL, *wcsfn2=NULL;
  
  while ((c = getopt(argc, args, OPTIONS)) != -1) {
    switch(c) {
    case 'v':
      loglvl++;
      break;
    case 'h':
      print_help(args[0]);
      exit(0);
    case '1':
      wcsfn1 = optarg;
      break;
    case '2':
      wcsfn2 = optarg;
      break;
    }
  }
  log_init(loglvl);
  if (optind != argc) {
    print_help(args[0]);
    exit(-1);
  }
  if (!(wcsfn1) || !(wcsfn2)) {
    print_help(args[0]);
    exit(-1);
  }
  /* open the two wcs systems */
  wcs1 = anwcs_open(wcsfn1, ext);
  if (!wcs1) {
    ERROR("Failed to read WCS file");
    exit(-1);
  }
  logverb("Read WCS:\n");
  if (log_get_level() >= LOG_VERB) {
    anwcs_print(wcs1, log_get_fid());
  }
  wcs2 = anwcs_open(wcsfn2, ext);
  if (!wcs2) {
    ERROR("Failed to read WCS file");
    exit(-1);
  }
  logverb("Read WCS:\n");
  if (log_get_level() >= LOG_VERB) {
    anwcs_print(wcs2, log_get_fid());
  }
  
  /* setup the solver, start in the middle */

  gsl_vector_set(x,0,anwcs_imagew(wcs1)/2.0);
  gsl_vector_set(x,1,anwcs_imageh(wcs1)/2.0);
  T = gsl_multiroot_fsolver_hybrids;
  s = gsl_multiroot_fsolver_alloc (T,2);
  gsl_multiroot_fsolver_set(s,&f,x);
  print_state(iter,s);
  do {
    iter++;
    status = gsl_multiroot_fsolver_iterate(s);
    print_state(iter,s);
    if (status) break;
    status = gsl_multiroot_test_residual(s->f,1e-7);
  } while (status == GSL_CONTINUE && iter < 1000);
  sol[0]=gsl_vector_get(s->x,0);
  sol[1]=gsl_vector_get(s->x,1);


  /* write some diagnostics on stderr */
  /* transform to ra/dec */
  anwcs_pixelxy2radec(wcs1, sol[0], sol[1], &ra, &dec);
  if (loglvl > LOG_MSG)
    fprintf(stderr,"Pixel (%.10f, %.10f) -> RA,Dec (%.10f, %.10f)\n", 
	    sol[0], sol[1], ra, dec);
  /* transform to x/y with second wcs 
     center of rotation should stay the same x/y */
  anwcs_radec2pixelxy(wcs2, ra, dec, &sol[0], &sol[1]);
  if (loglvl > LOG_MSG)
    fprintf(stderr,"RA,Dec (%.10f, %.10f) -> Pixel (%.10f, %.10f) \n", 
	    ra, dec, sol[0], sol[1]);

  /* write the solution */
  fprintf(stdout,"%f\n",sol[0]); 
  fprintf(stdout,"%f\n",sol[1]);
  
  return(0);
}