#include #include #include #include #include #include #include "qfits.h" #include #define LIM 256 #define D2R M_PI/180.0 #define R2D 180.0/M_PI #define NMAX 4096 struct catalog { int n; float x[NMAX],y[NMAX]; double ra[NMAX],de[NMAX]; double rx[NMAX],ry[NMAX]; float xres[NMAX],yres[NMAX],res[NMAX]; float xrms,yrms,rms; int usage[NMAX]; }; struct image { int naxis1,naxis2,nframes; float *zavg,*zstd,*zmax,*znum; double ra0,de0; float x0,y0; float a[2],b[2]; double mjd; float *dt; }; struct transformation { double ra0,de0; float a[3],b[3]; float x0,y0; }; int fgetline(FILE *file,char *s,int lim); void forward(double ra0,double de0,double ra,double de,double *x,double *y); void reverse(double ra0,double de0,double x,double y,double *ra,double *de); struct catalog read_catalog(char *filename); void lfit2d(float *x,float *y,float *z,int n,float *a); void add_fits_keywords(struct transformation t,char *filename); struct image read_fits(char *filename); // Modify FITS keywords void modify_fits_keywords(struct transformation t,char *filename) { char card[FITS_LINESZ+1]; char key[FITS_LINESZ+1]; char val[FITS_LINESZ+1]; char com[FITS_LINESZ+1]; sprintf(val,"%f",t.x0); keytuple2str(card,"CRPIX1",val,""); qfits_replace_card(filename,"CRPIX1",card); sprintf(val,"%f",t.y0); keytuple2str(card,"CRPIX2",val,""); qfits_replace_card(filename,"CRPIX2",card); sprintf(val,"%f",t.ra0); keytuple2str(card,"CRVAL1",val,""); qfits_replace_card(filename,"CRVAL1",card); sprintf(val,"%f",t.de0); keytuple2str(card,"CRVAL2",val,""); qfits_replace_card(filename,"CRVAL2",card); sprintf(val,"%e",t.a[1]/3600.0); keytuple2str(card,"CD1_1",val,""); qfits_replace_card(filename,"CD1_1",card); sprintf(val,"%e",t.a[2]/3600.0); keytuple2str(card,"CD1_2",val,""); qfits_replace_card(filename,"CD1_2",card); sprintf(val,"%e",t.b[1]/3600.0); keytuple2str(card,"CD2_1",val,""); qfits_replace_card(filename,"CD2_1",card); sprintf(val,"%e",t.b[2]/3600.0); keytuple2str(card,"CD2_2",val,""); qfits_replace_card(filename,"CD2_2",card); return; } int main(int argc,char *argv[]) { int i,j,k,l,m; struct catalog c; struct transformation t; double ra0,de0; float rmsmin; float x[NMAX],y[NMAX],rx[NMAX],ry[NMAX]; struct image img; char filename[128]; if (argc==1) strcpy(filename,"test.fits"); else if (argc==2) strcpy(filename,argv[1]); img=read_fits(filename); printf("files read\n"); c=read_catalog("out.dat"); printf("files read\n"); // Initial fit t.ra0=c.ra[0]; t.de0=c.de[0]; t.x0=(float) img.naxis1/2.0; t.y0=(float) img.naxis2/2.0; for (l=0;l<10;l++) { for (j=0;j<5;j++) { // Transform for (i=0;i2*c.rms) c.usage[i]=0; } } printf("%12.8lf %10.6lf %10.6lf %8.4f %8.4f %8.4f %8.4f\n",img.mjd,t.ra0,t.de0,t.a[1],t.a[2],t.b[1],t.b[2]); printf("%d/%d %f %f %f\n",m,c.n,c.xrms,c.yrms,c.rms); // add_fits_keywords(t,"test.fits"); modify_fits_keywords(t,filename); return 0; } // Read a line of maximum length int lim from file FILE into string s int fgetline(FILE *file,char *s,int lim) { int c,i=0; while (--lim > 0 && (c=fgetc(file)) != EOF && c != '\n') s[i++] = c; if (c == '\n') s[i++] = c; s[i] = '\0'; return i; } // Read catalog struct catalog read_catalog(char *filename) { int i=0; char line[LIM]; FILE *file; struct catalog c; file=fopen(filename,"r"); while (fgetline(file,line,LIM)>0) { sscanf(line,"%f %f %lf %lf",&c.x[i],&c.y[i],&c.ra[i],&c.de[i]); c.usage[i]=1; i++; } fclose(file); c.n=i; return c; } // Linear 2D fit void lfit2d(float *x,float *y,float *z,int n,float *a) { int i,j,m; double chisq; gsl_matrix *X,*cov; gsl_vector *yy,*w,*c; X=gsl_matrix_alloc(n,3); yy=gsl_vector_alloc(n); w=gsl_vector_alloc(n); c=gsl_vector_alloc(3); cov=gsl_matrix_alloc(3,3); // Fill matrices for(i=0;i