37 #include "casu_mods.h"
38 #include "catalogue/casu_utils.h"
39 #include "catalogue/casu_wcsutils.h"
40 #include "casu_stats.h"
42 static int casu_plate6(
double *xpos,
double *ypos,
double *xi,
double *eta,
43 unsigned char *flag,
int nstds,
double *a,
double *b,
44 double *c,
double *d,
double *e,
double *f);
45 static int casu_plate4(
double *xpos,
double *ypos,
double *xi,
double *eta,
46 unsigned char *flag,
int nstds,
double *a,
double *b,
47 double *c,
double *d,
double *e,
double *f);
48 static void tidy(
double *work,
unsigned char *iwork);
127 cpl_table *matchedstds,
int nconst,
int shiftan,
129 int nstds,nc2,i,niter,nrej,ngood,nreq=6,xcol,ycol,nulval;
131 const char *fctid =
"casu_platesol";
133 double *xptr,*yptr,*xptr2,*yptr2,*ra,*dec,*xiptr,*etaptr,*wptr,averr;
134 double r1,r2,d1,d2,newcrval1,newcrval2,a,b,c,d,e,f,xifit,etafit,dxi,deta;
135 double crpix1,crpix2,xi,eta,rac_before,rac_after,decc_before,decc_after;
136 double xcen,ycen,oldtheta,scale,oldcrpix1,oldcrpix2,*dptr;
138 const double *crdata;
139 unsigned char *isbad,*wptr2,*iwork=NULL;
140 const char *reqcols[] = {
"X_coordinate",
"Y_coordinate",
"xpredict",
141 "ypredict",
"RA",
"Dec"};
145 const cpl_matrix *crm;
149 if (*status != CASU_OK)
154 if (nconst != 4 && nconst != 6) {
156 "Value of nconst = %" CPL_SIZE_FORMAT
" is unsupported",
163 nstds = (int)cpl_table_get_nrow(matchedstds);
167 "Too few standards (%" CPL_SIZE_FORMAT
") in table for %" CPL_SIZE_FORMAT
" coefficient fit",
168 (cpl_size)nstds,(cpl_size)nconst);
174 for (i = 0; i < nreq; i++) {
175 if (cpl_table_has_column(matchedstds,reqcols[i]) != 1) {
176 cpl_msg_error(fctid,
"Matched standards table missing column %s",
184 work = cpl_malloc(10*nstds*
sizeof(*work));
185 iwork = cpl_calloc(3*nstds,
sizeof(*isbad));
188 xptr2 = work + 2*nstds;
189 yptr2 = work + 3*nstds;
191 dec = work + 5*nstds;
192 xiptr = work + 6*nstds;
193 etaptr = work + 7*nstds;
194 wptr = work + 8*nstds;
196 wptr2 = iwork + nstds;
201 tptr = cpl_table_get_data_float(matchedstds,
"X_coordinate");
202 for (i = 0; i < nstds; i++)
203 xptr[i] = (
double)tptr[i];
204 tptr = cpl_table_get_data_float(matchedstds,
"Y_coordinate");
205 for (i = 0; i < nstds; i++)
206 yptr[i] = (
double)tptr[i];
207 tptr = cpl_table_get_data_float(matchedstds,
"xpredict");
208 for (i = 0; i < nstds; i++)
209 xptr2[i] = (
double)tptr[i];
210 tptr = cpl_table_get_data_float(matchedstds,
"ypredict");
211 for (i = 0; i < nstds; i++)
212 yptr2[i] = (
double)tptr[i];
213 if (cpl_table_get_column_type(matchedstds,
"RA") == CPL_TYPE_FLOAT) {
214 tptr = cpl_table_get_data_float(matchedstds,
"RA");
215 for (i = 0; i < nstds; i++)
216 ra[i] = (
double)tptr[i];
217 tptr = cpl_table_get_data_float(matchedstds,
"Dec");
218 for (i = 0; i < nstds; i++)
219 dec[i] = (
double)tptr[i];
221 dptr = cpl_table_get_data_double(matchedstds,
"RA");
222 for (i = 0; i < nstds; i++)
224 dptr = cpl_table_get_data_double(matchedstds,
"Dec");
225 for (i = 0; i < nstds; i++)
235 wcs = cpl_wcs_new_from_propertylist(plist);
236 cr = cpl_wcs_get_crval(wcs);
237 crdata = cpl_array_get_data_double_const(cr);
238 for (i = 0; i < nstds; i++) {
246 newcrval1 = crdata[0] + r1;
247 newcrval2 = crdata[1] + d1;
248 cpl_propertylist_update_double(plist,
"CRVAL1",newcrval1);
249 cpl_propertylist_update_double(plist,
"CRVAL2",newcrval2);
255 wcs = cpl_wcs_new_from_propertylist(plist);
256 n1 = cpl_array_get_int(cpl_wcs_get_image_dims(wcs),0,&nulval);
257 n2 = cpl_array_get_int(cpl_wcs_get_image_dims(wcs),1,&nulval);
258 cr = cpl_wcs_get_crpix(wcs);
259 crdata = cpl_array_get_data_double_const(cr);
260 oldcrpix1 = crdata[0];
261 oldcrpix2 = crdata[1];
262 xcen = 0.5*(double)n1;
263 ycen = 0.5*(double)n2;
268 for (i = 0; i < nstds; i++) {
284 *status = casu_plate6(xptr,yptr,xiptr,etaptr,isbad,nstds,&a,&b,
288 *status = casu_plate4(xptr,yptr,xiptr,etaptr,isbad,nstds,&a,&b,
292 *status = casu_plate6(xptr,yptr,xiptr,etaptr,isbad,nstds,&a,&b,
296 if (*status != CASU_OK) {
297 cpl_msg_error(fctid,
"Plate constant solution failed");
304 for (i = 0; i < nstds; i++) {
305 xifit = xptr[i]*a + yptr[i]*b + c;
306 etafit = xptr[i]*d + yptr[i]*e + f;
307 dxi = fabs(xifit - xiptr[i]);
308 deta = fabs(etafit - etaptr[i]);
311 wptr2[i*2] = isbad[i];
312 wptr2[i*2+1] = isbad[i];
320 for (i = 0; i < nstds; i++) {
321 if (!isbad[i] && (wptr[i*2] > 3.0*averr || wptr[i*2+1] > 3.0*averr)) nrej++;
326 if (nrej == 0 || ngood < nconst)
328 for (i = 0; i < nstds; i++) {
329 if (!isbad[i] && (wptr[i*2] > 3.0*averr || wptr[i*2+1] > 3.0*averr)) isbad[i] = 1;
336 crpix1 = (e*c - b*f)/(d*b - e*a);
337 crpix2 = (a*f - d*c)/(d*b - e*a);
346 for (i = 0; i < nstds; i++)
349 averr *= DEGRAD*3600.0;
355 cpl_propertylist_update_double(plist,
"CD1_1",a);
356 cpl_propertylist_update_double(plist,
"CD1_2",b);
357 cpl_propertylist_update_double(plist,
"CD2_1",d);
358 cpl_propertylist_update_double(plist,
"CD2_2",e);
359 cpl_propertylist_update_int(plist,
"ESO DRS NUMBRMS",ngood);
360 cpl_propertylist_set_comment(plist,
"ESO DRS NUMBRMS",
361 "Number of stars in WCS fit");
362 cpl_propertylist_update_double(plist,
"ESO DRS STDCRMS",averr);
363 cpl_propertylist_set_comment(plist,
"ESO DRS STDCRMS",
364 "[arcsec] Average error in WCS fit");
368 cpl_propertylist_update_string(plist,
"CUNIT1",
"deg");
369 cpl_propertylist_set_comment(plist,
"CUNIT1",
370 "Unit of coordinate transformation");
371 cpl_propertylist_update_string(plist,
"CUNIT2",
"deg");
372 cpl_propertylist_set_comment(plist,
"CUNIT2",
373 "Unit of coordinate transformation");
374 cpl_propertylist_update_double(plist,
"CSYER1",0.0);
375 cpl_propertylist_set_comment(plist,
"CSYER1",
"Systematic error axis 1");
376 cpl_propertylist_update_double(plist,
"CSYER2",0.0);
377 cpl_propertylist_set_comment(plist,
"CSYER2",
"Systematic error axis 2");
378 cpl_propertylist_update_double(plist,
"CRDER1",averr);
379 cpl_propertylist_set_comment(plist,
"CRDER1",
"Random error axis 1");
380 cpl_propertylist_update_double(plist,
"CRDER2",averr);
381 cpl_propertylist_set_comment(plist,
"CRDER2",
"Random error axis 2");
382 cpl_propertylist_update_string(plist,
"RADECSYS",
"ICRS");
383 cpl_propertylist_set_comment(plist,
"RADECSYS",
"Coordinate reference frame");
384 cpl_propertylist_update_double(plist,
"EQUINOX",2000.0);
385 cpl_propertylist_set_comment(plist,
"EQUINOX",
"Coordinate equinox");
389 crm = cpl_wcs_get_cd(wcs);
390 crdata = cpl_matrix_get_data_const(crm);
391 oldtheta = 0.5*(fabs(atan2(crdata[1],crdata[0])) +
392 fabs(atan2(crdata[2],crdata[3])));
394 wcs = cpl_wcs_new_from_propertylist(plist);
396 xcen = 3600.0*(rac_after - rac_before);
397 ycen = 3600.0*(decc_after - decc_before);
398 cpl_propertylist_update_double(plist,
"ESO DRS WCSRAOFF",xcen);
399 cpl_propertylist_set_comment(plist,
"ESO DRS WCSRAOFF",
400 "[arcsec] cenpix RA_after - RA_before)");
401 cpl_propertylist_update_double(plist,
"ESO DRS WCSDECOFF",ycen);
402 cpl_propertylist_set_comment(plist,
"ESO DRS WCSDECOFF",
403 "[arcsec] cenpix Dec_after - Dec_before)");
410 if (xcol != -1 && ycol != -1) {
411 snprintf(key,9,
"TCRPX%d",xcol);
413 snprintf(key,9,
"TCRPX%d",ycol);
416 snprintf(key,9,
"TCRVL%d",xcol);
417 cpl_propertylist_update_double(tlist,key,newcrval1);
418 snprintf(key,9,
"TCRVL%d",ycol);
419 cpl_propertylist_update_double(tlist,key,newcrval2);
421 snprintf(key,9,
"TC%d_%d",xcol,xcol);
422 cpl_propertylist_update_double(tlist,key,a);
423 snprintf(key,9,
"TC%d_%d",xcol,ycol);
424 cpl_propertylist_update_double(tlist,key,b);
425 snprintf(key,9,
"TC%d_%d",ycol,xcol);
426 cpl_propertylist_update_double(tlist,key,d);
427 snprintf(key,9,
"TC%d_%d",ycol,ycol);
428 cpl_propertylist_update_double(tlist,key,e);
429 cpl_propertylist_update_int(tlist,
"ESO DRS NUMBRMS",ngood);
430 cpl_propertylist_set_comment(tlist,
"ESO DRS NUMBRMS",
431 "Number of stars in WCS fit");
432 cpl_propertylist_update_double(tlist,
"ESO DRS STDCRMS",averr);
433 cpl_propertylist_set_comment(tlist,
"ESO DRS STDCRMS",
434 "[arcsec] Average error in WCS fit");
435 cpl_propertylist_update_double(tlist,
"ESO DRS WCSRAOFF",xcen);
436 cpl_propertylist_set_comment(tlist,
"ESO DRS WCSRAOFF",
437 "[arcsec] cenpix RA_after - RA_before)");
438 cpl_propertylist_update_double(tlist,
"ESO DRS WCSDECOFF",ycen);
439 cpl_propertylist_set_comment(tlist,
"ESO DRS WCSDECOFF",
440 "[arcsec] cenpix Dec_after - Dec_before)");
444 snprintf(key,9,
"TCUNI%d",xcol);
445 cpl_propertylist_update_string(tlist,key,
"deg");
446 cpl_propertylist_set_comment(tlist,key,
447 "Unit of coordinate transformation");
448 snprintf(key,9,
"TCUNI%d",ycol);
449 cpl_propertylist_update_string(tlist,key,
"deg");
450 cpl_propertylist_set_comment(tlist,key,
451 "Unit of coordinate transformation");
452 snprintf(key,9,
"TCSY%d",xcol);
453 cpl_propertylist_update_double(tlist,key,0.0);
454 cpl_propertylist_set_comment(tlist,key,
"Systematic error x-axis ");
455 snprintf(key,9,
"TCSY%d",ycol);
456 cpl_propertylist_update_double(tlist,key,0.0);
457 cpl_propertylist_set_comment(tlist,key,
"Systematic error y-axis");
458 snprintf(key,9,
"TCRD%d",xcol);
459 cpl_propertylist_update_double(tlist,key,averr);
460 cpl_propertylist_set_comment(tlist,key,
"Random error axis 1");
461 snprintf(key,9,
"TCRD%d",ycol);
462 cpl_propertylist_update_double(tlist,key,averr);
463 cpl_propertylist_set_comment(tlist,key,
"Random error axis 2");
464 cpl_propertylist_update_string(tlist,
"RADECSYS",
"ICRS");
465 cpl_propertylist_set_comment(tlist,
"RADECSYS",
466 "Coordinate reference frame");
467 cpl_propertylist_update_double(tlist,
"EQUINOX",2000.0);
468 cpl_propertylist_set_comment(tlist,
"EQUINOX",
"Coordinate equinox");
475 cr = cpl_wcs_get_crval(wcs);
476 crdata = cpl_array_get_data_double_const(cr);
478 rac_after -= crdata[0];
479 decc_after -= crdata[1];
480 cpl_propertylist_update_double(plist,
"ESO QC WCS_DCRVAL1",rac_after);
481 cpl_propertylist_set_comment(plist,
"ESO QC WCS_DCRVAL1",
482 "[deg] change in crval1");
483 cpl_propertylist_update_double(plist,
"ESO QC WCS_DCRVAL2",decc_after);
484 cpl_propertylist_set_comment(plist,
"ESO QC WCS_DCRVAL2",
485 "[deg] change in crval2");
489 crm = cpl_wcs_get_cd(wcs);
490 crdata = cpl_matrix_get_data_const(crm);
491 oldtheta = 0.5*(fabs(atan2(crdata[1],crdata[0])) +
492 fabs(atan2(crdata[2],crdata[3]))) - oldtheta;
494 cpl_propertylist_update_double(plist,
"ESO QC WCS_DTHETA",oldtheta);
495 cpl_propertylist_set_comment(plist,
"ESO QC WCS_DTHETA",
496 "[deg] change in rotation");
500 scale = 1800.0*(sqrt(pow(crdata[0],2.0) + pow(crdata[1],2.0)) +
501 sqrt(pow(crdata[2],2.0) + pow(crdata[3],2.0)));
502 cpl_propertylist_update_double(plist,
"ESO QC WCS_SCALE",scale);
503 cpl_propertylist_set_comment(plist,
"ESO QC WCS_SCALE",
504 "[arcsec] mean plate scale");
508 oldtheta = fabs(atan2(crdata[1],crdata[0])) -
509 fabs(atan2(crdata[2],crdata[3]));
510 cpl_propertylist_update_double(plist,
"ESO QC WCS_SHEAR",oldtheta);
511 cpl_propertylist_set_comment(plist,
"ESO QC WCS_SHEAR",
512 "[deg] abs(xrot) - abs(yrot)");
516 cpl_propertylist_update_double(plist,
"ESO QC WCS_RMS",averr);
517 cpl_propertylist_set_comment(plist,
"ESO QC WCS_RMS",
518 "[arcsec] Average error in WCS fit");
569 extern int casu_platexy(cpl_table *matchedxy,
int nconst, cpl_array **coefs,
571 const char *fctid =
"casu_platexy";
572 int nxy,nc2,i,niter,nrej,ngood;
573 unsigned char *isbad,*wptr2,*iwork=NULL;
574 double *xptr1,*xptr2,*yptr1,*yptr2,*wptr,a,b,c,d,e,f,*cc,xfit,yfit;
575 double dx,dy,averr,*work=NULL;
578 const char *reqcols[] = {
"X_coordinate_1",
"Y_coordinate_1",
"X_coordinate_2",
584 if (*status != CASU_OK)
589 if (nconst != 4 && nconst != 6) {
591 "Value of nconst = %" CPL_SIZE_FORMAT
" is unsupported",
598 nxy = (int)cpl_table_get_nrow(matchedxy);
602 "Too few objects (%" CPL_SIZE_FORMAT
") in table for %" CPL_SIZE_FORMAT
" coefficient fit",
603 (cpl_size)nxy,(cpl_size)nconst);
609 for (i = 0; i < nreq; i++) {
610 if (cpl_table_has_column(matchedxy,reqcols[i]) != 1) {
611 cpl_msg_error(fctid,
"Input table missing column %s\n",reqcols[i]);
618 work = cpl_malloc(6*nxy*
sizeof(*work));
619 iwork = cpl_calloc(3*nxy,
sizeof(*isbad));
622 xptr2 = work + 2*nxy;
623 yptr2 = work + 3*nxy;
631 tptr = cpl_table_get_data_float(matchedxy,
"X_coordinate_1");
632 for (i = 0; i < nxy; i++)
633 xptr1[i] = (
double)tptr[i];
634 tptr = cpl_table_get_data_float(matchedxy,
"Y_coordinate_1");
635 for (i = 0; i < nxy; i++)
636 yptr1[i] = (
double)tptr[i];
637 tptr = cpl_table_get_data_float(matchedxy,
"X_coordinate_2");
638 for (i = 0; i < nxy; i++)
639 xptr2[i] = (
double)tptr[i];
640 tptr = cpl_table_get_data_float(matchedxy,
"Y_coordinate_2");
641 for (i = 0; i < nxy; i++)
642 yptr2[i] = (
double)tptr[i];
654 *status = casu_plate6(xptr2,yptr2,xptr1,yptr1,isbad,nxy,&a,&b,
658 *status = casu_plate4(xptr2,yptr2,xptr1,yptr1,isbad,nxy,&a,&b,
662 *status = casu_plate6(xptr2,yptr2,xptr1,yptr1,isbad,nxy,&a,&b,
666 if (*status != CASU_OK) {
667 cpl_msg_error(fctid,
"Plate constant solution failed");
674 for (i = 0; i < nxy; i++) {
675 xfit = xptr2[i]*a + yptr2[i]*b + c;
676 yfit = xptr2[i]*d + yptr2[i]*e + f;
677 dx = fabs(xfit - xptr1[i]);
678 dy = fabs(yfit - yptr1[i]);
681 wptr2[i*2] = isbad[i];
682 wptr2[i*2+1] = isbad[i];
690 for (i = 0; i < nxy; i++) {
691 if (!isbad[i] && (wptr[i*2] > 3.0*averr || wptr[i*2+1] > 3.0*averr)) nrej++;
696 if (nrej == 0 || ngood < nconst)
698 for (i = 0; i < nxy; i++) {
699 if (!isbad[i] && (wptr[i*2] > 3.0*averr || wptr[i*2+1] > 3.0*averr)) isbad[i] = 1;
706 *coefs = cpl_array_new(6,CPL_TYPE_DOUBLE);
707 cc = cpl_array_get_data_double(*coefs);
757 static int casu_plate6(
double *xpos,
double *ypos,
double *xi,
double *eta,
758 unsigned char *flag,
int nstds,
double *a,
double *b,
759 double *c,
double *d,
double *e,
double *f) {
760 double sx1sq,sy1sq,sx1y1,sx1x2,sy1x2;
761 double sy1y2,sx1y2,xposmean,yposmean,ximean,etamean,xx1,yy1,xx2,yy2;
767 ngood = nstds - nbad;
790 for (i = 0; i < nstds; i++) {
792 xx1 = xpos[i] - xposmean;
793 yy1 = ypos[i] - yposmean;
794 xx2 = xi[i] - ximean;
795 yy2 = eta[i] - etamean;
808 *a = (sx1y1*sy1x2 - sx1x2*sy1sq)/(sx1y1*sx1y1 - sx1sq*sy1sq);
809 *b = (sx1x2*sx1y1 - sx1sq*sy1x2)/(sx1y1*sx1y1 - sx1sq*sy1sq);
810 *c = -xposmean*(*a) - yposmean*(*b) + ximean;
814 *d = (sx1y1*sx1y2 - sy1y2*sx1sq)/(sx1y1*sx1y1 - sy1sq*sx1sq);
815 *e = (sy1y2*sx1y1 - sy1sq*sx1y2)/(sx1y1*sx1y1 - sy1sq*sx1sq);
816 *f = -xposmean*(*e) - yposmean*(*d) + etamean;
859 static int casu_plate4(
double *xpos,
double *ypos,
double *xi,
double *eta,
860 unsigned char *flag,
int nstds,
double *a,
double *b,
861 double *c,
double *d,
double *e,
double *f) {
862 double sx1sq,sy1sq,sx1x2,sy1x2,sy1y2,sx1y2,xposmean,yposmean;
863 double ximean,etamean,xx1,yy1,xx2,yy2,det,num,denom,theta,mag;
864 double stheta,ctheta;
870 ngood = nstds - nbad;
892 for (i = 0; i < nstds; i++) {
894 xx1 = xpos[i] - xposmean;
895 yy1 = ypos[i] - yposmean;
896 xx2 = xi[i] - ximean;
897 yy2 = eta[i] - etamean;
909 det = sx1x2*sy1y2 - sy1x2*sx1y2;
912 denom = -sx1x2 + sy1y2;
915 denom = sx1x2 + sy1y2;
917 if (num == 0.0 && denom == 0.0) {
920 theta = atan2(num,denom);
922 theta += CPL_MATH_2PI;
929 num = denom*ctheta + num*stheta;
930 denom = sx1sq + sy1sq;
948 *c = -xposmean*(*a) - yposmean*(*b) + ximean;
949 *f = -xposmean*(*e) - yposmean*(*d) + etamean;
956 static void tidy(
double *work,
unsigned char *iwork) {
957 if(work != NULL) {cpl_free(work);}
959 if(iwork != NULL) {cpl_free(iwork);}
int casu_platexy(cpl_table *matchedxy, int nconst, cpl_array **coefs, int *status)
Work out an xy transformation for two coodinate lists.
int casu_platesol(cpl_propertylist *plist, cpl_propertylist *tlist, cpl_table *matchedstds, int nconst, int shiftan, int *status)
Work out a WCS for an image.
double casu_dmean(double *data, unsigned char *bpm, long npts)
double casu_dmed(double *data, unsigned char *bpm, long npts)
int casu_sumbpm(unsigned char *bpm, long npts, int *sumb)
int casu_findcol(cpl_propertylist *p, const char *col)
Find from catalogue header which are x,y columns.
void casu_propertylist_update_double(cpl_propertylist *plist, const char *name, double val)
void casu_radectoxieta(cpl_wcs *wcs, double ra, double dec, double *xi, double *eta)
void casu_xytoradec(cpl_wcs *wcs, double x, double y, double *ra, double *dec)