MUSE Pipeline Reference Manual  2.1.1
muse_wcs.c
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2005-2015 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20  */
21 
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25 
26 /*----------------------------------------------------------------------------*
27  * Includes *
28  *----------------------------------------------------------------------------*/
29 #include <cpl.h>
30 #include <string.h>
31 
32 #include "muse_wcs.h"
33 #include "muse_instrument.h"
34 
35 #include "muse_astro.h"
36 #include "muse_combine.h"
37 #include "muse_dar.h"
38 #include "muse_pfits.h"
39 #include "muse_quality.h"
40 #include "muse_resampling.h"
41 #include "muse_utils.h"
42 
43 /*----------------------------------------------------------------------------*
44  * Debugging Macros *
45  * Set these to 1 or higher for (lots of) debugging output *
46  *----------------------------------------------------------------------------*/
47 #define FAKE_POS_ROT 0 /* activate some fake positions+rotation for debugging */
48 
49 /*----------------------------------------------------------------------------*/
53 /*----------------------------------------------------------------------------*/
54 
57 /*----------------------------------------------------------------------------*/
69 /*----------------------------------------------------------------------------*/
72 {
73  muse_wcs_object *wcs = cpl_calloc(1, sizeof(muse_wcs_object));
74  return wcs;
75 }
76 
77 /*----------------------------------------------------------------------------*/
87 /*----------------------------------------------------------------------------*/
88 void
90 {
91  if (!aWCSObj) {
92  return;
93  }
94  muse_datacube_delete(aWCSObj->cube);
95  aWCSObj->cube = NULL;
96  cpl_table_delete(aWCSObj->detected);
97  aWCSObj->detected = NULL;
98  cpl_propertylist_delete(aWCSObj->wcs);
99  aWCSObj->wcs = NULL;
100  cpl_free(aWCSObj);
101 }
102 
103 /*---------------------------------------------------------------------------*/
115 /*---------------------------------------------------------------------------*/
117  { "XPOS", CPL_TYPE_DOUBLE, "pix", "%f", "horizontal position", CPL_TRUE },
118  { "YPOS", CPL_TYPE_DOUBLE, "pix", "%f", "vertical position", CPL_TRUE },
119  { "XERR", CPL_TYPE_DOUBLE, "pix", "%f",
120  "error estimate of the horizontal position", CPL_TRUE },
121  { "YERR", CPL_TYPE_DOUBLE, "pix", "%f",
122  "error estimate of the vertical position", CPL_TRUE },
123  { "FLUX", CPL_TYPE_DOUBLE, "count", "%e", "(relative) flux of the source", CPL_TRUE },
124  { "XFWHM", CPL_TYPE_DOUBLE, "pix", "%f", "horizontal FWHM", CPL_TRUE },
125  { "YFWHM", CPL_TYPE_DOUBLE, "pix", "%f", "vertical FWHM", CPL_TRUE },
126  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
127 };
128 
129 /*---------------------------------------------------------------------------*/
140 /*---------------------------------------------------------------------------*/
142  { "SourceID", CPL_TYPE_STRING, "", "%s", "the source identification", CPL_TRUE },
143  { "RA", CPL_TYPE_DOUBLE, "deg", "%f", "right ascension (decimal degrees)", CPL_TRUE },
144  { "DEC", CPL_TYPE_DOUBLE, "deg", "%f", "declination (decimal degrees)", CPL_TRUE },
145  { "filter", CPL_TYPE_STRING, "", "%s", "the filter used for the magnitude", CPL_TRUE },
146  { "mag", CPL_TYPE_DOUBLE, "mag", "%f", "the object (Vega) magnitude", CPL_TRUE },
147  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
148 };
149 
150 
151 /*----------------------------------------------------------------------------*/
176 /*----------------------------------------------------------------------------*/
177 cpl_table *
178 muse_wcs_centroid_stars(muse_image *aImage, float aSigma,
179  muse_wcs_centroid_type aCentroid)
180 {
181  cpl_ensure(aImage && aImage->data && aImage->stat, CPL_ERROR_NULL_INPUT,
182  NULL);
183  cpl_ensure(aSigma > 0., CPL_ERROR_ILLEGAL_INPUT, NULL);
184  switch (aCentroid) {
186  cpl_msg_info(__func__, "Gaussian profile fits for object centroiding");
187  break;
189  cpl_msg_info(__func__, "Moffat profile fits for object centroiding");
190  break;
192  cpl_msg_info(__func__, "Simple square box object centroiding");
193  break;
194  default:
195  cpl_msg_error(__func__, "Unknown centroiding method!");
196  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
197  return NULL;
198  }
199 
200  /* muse_image_duplicate() doesn't work here, because we don't want *
201  * to require the dq component, so do it "by hand"; then we can *
202  * convert the variance to an error = sqrt(variance) at the same time */
203  muse_image *image = muse_image_new();
204  image->data = cpl_image_duplicate(aImage->data);
205  image->stat = cpl_image_power_create(aImage->stat, 0.5);
206  if (aImage->header) {
207  image->header = cpl_propertylist_duplicate(aImage->header);
208  }
209  /* make sure to exclude bad pixels from the fits below, *
210  * either from an existing dq component or by using NANs */
211  if (aImage->dq) {
212  image->dq = cpl_image_duplicate(aImage->dq);
213  muse_quality_image_reject_using_dq(image->data, image->dq, image->stat);
214  } else {
215  cpl_image_reject_value(image->data, CPL_VALUE_NAN); /* mark NANs */
216  cpl_image_reject_value(image->stat, CPL_VALUE_NAN);
217  }
218  /* since cpl_image_get_fwhm() ignores the image masks, *
219  * treat the masked areas so that they don't contribute */
220  cpl_image_fill_rejected(image->data, cpl_image_get_min(image->data));
221  cpl_image_fill_rejected(image->stat, cpl_image_get_max(image->stat));
222 
223  /* do the source detection */
224  cpl_apertures *apertures = cpl_apertures_extract_sigma(image->data, aSigma);
225  int ndet = apertures ? cpl_apertures_get_size(apertures) : 0;
226  if (ndet < 1) {
227  muse_image_delete(image);
228  cpl_msg_error(__func__, "No sources for astrometric calibration found down "
229  "to %.3f sigma limit", aSigma);
230  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
231  return NULL;
232  }
233  double med_dist, median = cpl_image_get_median_dev(image->data, &med_dist),
234  threshold = median + aSigma * med_dist;
235  cpl_msg_debug(__func__, "The %.3f sigma threshold (%.3f, %.3f +/- %.3f) was "
236  "used to find %d sources", aSigma, threshold, median, med_dist,
237  ndet);
238 
239  int nx = cpl_image_get_size_x(image->data),
240  ny = cpl_image_get_size_y(image->data);
241  /* parameters for the fits */
242  const double bgguess = cpl_image_get_median(image->data),
243  beta = 2.5, /* moffat beta */
244  moffat_alpha_to_fwhm = 1 / (2 * sqrt(pow(2, 1/beta) - 1));
245 #if 0 /* unused idea to further robustify FWHM estimate */
246  double fwhmguess = (muse_pfits_get_fwhm_start(image->header)
247  + muse_pfits_get_fwhm_end(image->header)) / 2.;
248  if (fwhmguess < FLT_EPSILON) { /* in case FWHM headers are not there */
249  fwhmguess = 4.; /* [pix] assume median seeing in WFM */
250  } else { /* use nominal MUSE spaxel scale to get FWHM guess [pix] */
251  if (muse_pfits_get_mode(image->header) < MUSE_MODE_NFM_AO_N) {
252  fwhmguess /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeX_WFM) / 2;
253  } else {
254  fwhmguess /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeX_NFM) / 2;
255  }
256  }
257 #endif
258 
259  /* create the output table and loop over all apertures *
260  * to centroid the detected stars and fill the table */
261  cpl_table *detections = muse_cpltable_new(muse_wcs_detections_def, ndet);
262  cpl_table_unselect_all(detections);
263  int idet;
264  for (idet = 0; idet < ndet; idet++) {
265  double xc = cpl_apertures_get_centroid_x(apertures, idet+1),
266  yc = cpl_apertures_get_centroid_y(apertures, idet+1),
267  fluxaper = cpl_apertures_get_flux(apertures, idet+1);
268  double xwaper, ywaper;
269  cpl_image_get_fwhm(image->data, lround(xc), lround(yc), &xwaper, &ywaper);
270  /* Remove the detection if the FWHM could not be computed in both *
271  * axes, since this likely points to an object on the border or with *
272  * some other problem which should not be relied on for centroiding. */
273  if (xwaper < 0 || ywaper < 0) {
274  cpl_msg_debug(__func__, "FWHM computation unsuccessful at %f,%f, result "
275  "was %.3f,%.3f", xc, yc, xwaper, ywaper);
276  cpl_table_select_row(detections, idet);
277  continue;
278  }
279 
280  /* set the halfsize of the window to measure the centroids *
281  * with respect to the detection; +/-5 pix should suffice */
282  int x1 = floor(xc) - 5,
283  x2 = ceil(xc) + 5,
284  y1 = floor(yc) - 5,
285  y2 = ceil(yc) + 5;
286  /* force window to be inside the image */
287  if (x1 < 1) x1 = 1;
288  if (y1 < 1) y1 = 1;
289  if (x2 > nx) x2 = nx;
290  if (y2 > ny) y2 = ny;
291 
292  double xcen, ycen, xerr = 0., yerr = 0., xw, yw, flux;
293  cpl_errorstate state = cpl_errorstate_get();
294  switch (aCentroid) {
296  cpl_array *pgauss = cpl_array_new(7, CPL_TYPE_DOUBLE),
297  *pgerr = cpl_array_new(7, CPL_TYPE_DOUBLE),
298  *pgfit = cpl_array_new(7, CPL_TYPE_INT);
299  /* first, only fit the centroid */
300  cpl_array_set(pgfit, 3, 1); /* xc */
301  cpl_array_set(pgfit, 4, 1); /* yc */
302  cpl_array_set(pgfit, 0, 0);
303  cpl_array_set(pgfit, 1, 0);
304  cpl_array_set(pgfit, 2, 0);
305  cpl_array_set(pgfit, 5, 0);
306  cpl_array_set(pgfit, 6, 0);
307  cpl_array_set(pgauss, 3, xc);
308  cpl_array_set(pgauss, 4, yc);
309  cpl_array_set(pgauss, 0, bgguess);
310  cpl_array_set(pgauss, 1, fluxaper);
311  cpl_array_set(pgauss, 2, 0.); /* rho */
312  cpl_array_set(pgauss, 5, xwaper * CPL_MATH_SIG_FWHM);
313  cpl_array_set(pgauss, 6, ywaper * CPL_MATH_SIG_FWHM);
314  cpl_fit_image_gaussian(image->data, image->stat, lround(xc), lround(yc),
315  x2-x1+1, y2-y1+1, pgauss, pgerr, pgfit, NULL, NULL,
316  NULL, NULL, NULL, NULL, NULL);
317  xcen = cpl_array_get(pgauss, 3, NULL);
318  ycen = cpl_array_get(pgauss, 4, NULL);
319  xerr = cpl_array_get(pgerr, 3, NULL);
320  yerr = cpl_array_get(pgerr, 4, NULL);
321  /* second, keep center fixed, only fit flux and width */
322  cpl_array_set(pgfit, 3, 0); /* xc */
323  cpl_array_set(pgfit, 4, 0); /* yc */
324  cpl_array_set(pgfit, 1, 1); /* flux */
325  cpl_array_set(pgfit, 5, 1); /* sigma_x */
326  cpl_array_set(pgfit, 6, 1); /* sigma_y */
327  cpl_fit_image_gaussian(image->data, image->stat, lround(xc), lround(yc),
328  x2-x1+1, y2-y1+1, pgauss, pgerr, pgfit, NULL, NULL,
329  NULL, NULL, NULL, NULL, NULL);
330  xw = cpl_array_get(pgauss, 5, NULL) * CPL_MATH_FWHM_SIG;
331  yw = cpl_array_get(pgauss, 6, NULL) * CPL_MATH_FWHM_SIG;
332  flux = cpl_array_get(pgauss, 1, NULL);
333  cpl_array_delete(pgauss);
334  cpl_array_delete(pgerr);
335  cpl_array_delete(pgfit);
336  break;
337  }
339  cpl_size nmax = (x2-x1+1) * (y2-y1+1);
340  cpl_matrix *pos = cpl_matrix_new(nmax, 2);
341  cpl_vector *val = cpl_vector_new(nmax),
342  *err = cpl_vector_new(nmax);
343  int i, npix = 0;
344  for (i = x1; i <= x2; i++) {
345  int j;
346  for (j = y1; j <= y2; j++) {
347  int error;
348  double value = cpl_image_get(image->data, i, j, &error);
349  if (error != 0) {
350  continue;
351  }
352  cpl_matrix_set(pos, npix, 0, i);
353  cpl_matrix_set(pos, npix, 1, j);
354  cpl_vector_set(val, npix, value);
355  /* stat is already sigma! */
356  cpl_vector_set(err, npix, cpl_image_get(image->stat, i, j, &error));
357  npix++;
358  } /* for j */
359  } /* for i */
360  cpl_matrix_set_size(pos, npix, 2);
361  cpl_vector_set_size(val, npix);
362  cpl_vector_set_size(err, npix);
363  cpl_array *pmoffat = cpl_array_new(8, CPL_TYPE_DOUBLE),
364  *pmerror = cpl_array_new(8, CPL_TYPE_DOUBLE),
365  *pmfit = cpl_array_new(8, CPL_TYPE_INT);
366  /* first, only fit the centroid */
367  cpl_array_set(pmfit, 2, 1); /* xc */
368  cpl_array_set(pmfit, 3, 1); /* yc */
369  cpl_array_set(pmfit, 0, 0);
370  cpl_array_set(pmfit, 1, 0);
371  cpl_array_set(pmfit, 4, 0);
372  cpl_array_set(pmfit, 5, 0);
373  cpl_array_set(pmfit, 6, 0);
374  cpl_array_set(pmfit, 7, 0);
375  cpl_array_set(pmoffat, 2, xc);
376  cpl_array_set(pmoffat, 3, yc);
377  cpl_array_set(pmoffat, 0, bgguess);
378  cpl_array_set(pmoffat, 1, fluxaper);
379  cpl_array_set(pmoffat, 4, xwaper * moffat_alpha_to_fwhm);
380  cpl_array_set(pmoffat, 5, ywaper * moffat_alpha_to_fwhm);
381  cpl_array_set(pmoffat, 6, beta);
382  cpl_array_set(pmoffat, 7, 0.); /* rho */
383  muse_utils_fit_moffat_2d(pos, val, err, pmoffat, pmerror, pmfit, NULL, NULL);
384  xcen = cpl_array_get(pmoffat, 2, NULL);
385  ycen = cpl_array_get(pmoffat, 3, NULL);
386  xerr = cpl_array_get(pmerror, 2, NULL);
387  yerr = cpl_array_get(pmerror, 3, NULL);
388  /* second, keep center fixed, only fit flux and width */
389  cpl_array_set(pmfit, 2, 0); /* xc */
390  cpl_array_set(pmfit, 3, 0); /* yc */
391  cpl_array_set(pmfit, 1, 1); /* flux */
392  cpl_array_set(pmfit, 4, 1); /* xhwhm */
393  cpl_array_set(pmfit, 5, 1); /* yhwhm */
394  muse_utils_fit_moffat_2d(pos, val, err, pmoffat, pmerror, pmfit, NULL, NULL);
395  xw = cpl_array_get(pmoffat, 4, NULL) / moffat_alpha_to_fwhm;
396  yw = cpl_array_get(pmoffat, 5, NULL) / moffat_alpha_to_fwhm;
397  flux = cpl_array_get(pmoffat, 1, NULL);
398  cpl_array_delete(pmoffat);
399  cpl_array_delete(pmerror);
400  cpl_array_delete(pmfit);
401  cpl_matrix_delete(pos);
402  cpl_vector_delete(val);
403  cpl_vector_delete(err);
404  break;
405  }
406  default: /* MUSE_WCS_CENTROID_BOX is left */
407  muse_utils_image_get_centroid_window(image->data, x1, y1, x2, y2,
408  &xcen, &ycen,
410 #if 0
411  printf("%d apertures %f %f boxes %f %f deltas %f %f\n", idet+1, xc, yc,
412  xcen, ycen, xcen - xc, ycen - yc);
413  fflush(stdout);
414 #endif
415  /* set error to 0.15 pix which is the typical *
416  * stdev compared to Gaussian fits */
417  xerr = 0.15;
418  yerr = 0.15;
419  /* compute FWHM again with the revised central position */
420  cpl_image_get_fwhm(image->data, lround(xcen), lround(ycen), &xw, &yw);
421  flux = fluxaper; /* take the aperture flux */
422  }
423 
424  /* mandatory columns: position */
425  cpl_table_set(detections, "XPOS", idet, xcen);
426  cpl_table_set(detections, "YPOS", idet, ycen);
427  /* and the error on the position */
428  cpl_table_set(detections, "XERR", idet, xerr);
429  cpl_table_set(detections, "YERR", idet, yerr);
430  /* extra columns: FWHM... */
431  cpl_table_set(detections, "XFWHM", idet, xw);
432  cpl_table_set(detections, "YFWHM", idet, yw);
433  /* ... and flux */
434  cpl_table_set(detections, "FLUX", idet, flux);
435 
436  if (cpl_errorstate_is_equal(state) && xw > 0 && yw > 0 &&
437  xcen >= 1 && xcen <= nx && ycen >= 1 && ycen <= ny) {
438  continue;
439  }
440  /* some error occurred, so mark this entry for removal */
441  cpl_table_select_row(detections, idet);
442  } /* for idet (aperture index) */
443  cpl_table_erase_selected(detections);
444  cpl_apertures_delete(apertures);
445  muse_image_delete(image);
446  cpl_msg_debug(__func__, "%d of %d sources were recorded in the detections "
447  "table", (int)cpl_table_get_nrow(detections), ndet);
448 
449  return detections;
450 } /* muse_wcs_centroid_stars() */
451 
452 /*----------------------------------------------------------------------------*/
494 /*----------------------------------------------------------------------------*/
495 cpl_error_code
496 muse_wcs_locate_sources(muse_pixtable *aPixtable, float aSigma,
497  muse_wcs_centroid_type aCentroid,
498  muse_wcs_object *aWCSObj)
499 {
500  cpl_ensure_code(aPixtable && aPixtable->header && aWCSObj,
501  CPL_ERROR_NULL_INPUT);
502  cpl_ensure_code(aSigma > 0., CPL_ERROR_ILLEGAL_INPUT);
503  switch (aCentroid) {
507  break;
508  default:
509  return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
510  }
511  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) > 2) {
512  const char *fn = "wcs__pixtable.fits";
513  cpl_msg_info(__func__, "Saving pixel table as \"%s\"", fn);
514  muse_pixtable_save(aPixtable, fn);
515  }
516 
517  /* check that DAR correction was carried out in some way */
518  cpl_boolean darcorrected = CPL_FALSE;
519  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_NAME) &&
520  cpl_propertylist_get_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME) > 0.) {
521  darcorrected = CPL_TRUE;
522  } else if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
523  darcorrected = CPL_TRUE;
524  }
525  if (!darcorrected) {
526  const char *message = "Correction for differential atmospheric refraction "
527  "was not applied! Deriving astrometric correction "
528  "from such data is unlikely to give good results!";
529  cpl_msg_warning(__func__, "%s", message);
530  cpl_error_set_message(__func__, CPL_ERROR_UNSUPPORTED_MODE, "%s", message);
531  }
532 
533  muse_resampling_params *params =
535  params->pfx = 1.; /* large pixfrac to be sure to cover most gaps */
536  params->pfy = 1.;
537  params->pfl = 1.;
538  params->dlambda = 50.; /* we can integrate over lots of Angstrom here */
539  params->crtype = MUSE_RESAMPLING_CRSTATS_MEDIAN; /* median for clean cube */
540  params->crsigma = 25.; /* very moderate CR rejection */
541  muse_datacube *cube = muse_resampling_cube(aPixtable, params, NULL);
543  /* reset cosmic ray statuses in aPixtable, since the CR rejection *
544  * done here might not be appropriate for the final datacube */
545  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
546  if (!cube) {
547  return cpl_error_set_message(__func__, cpl_error_get_code(),
548  "Failure while creating cube!");
549  }
550  aWCSObj->cube = cube;
551  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
552  const char *fn = "wcs__cube.fits";
553  cpl_msg_info(__func__, "Saving cube as \"%s\"", fn);
554  muse_datacube_save(cube, fn);
555  }
556 
557  /* detect objects in the cube, using image planes around the central one */
558  int iplane, irefplane = cpl_imagelist_get_size(cube->data) / 2;
559  /* medianing three images removes all cosmic rays but continuum objects stay *
560  * (need to use muse_images and the muse_combine function because *
561  * cpl_imagelist_collapse_median_create() disregards all bad pixels) */
563  unsigned int ilist = 0;
564  for (iplane = irefplane - 1; iplane <= irefplane + 1; iplane++) {
565  muse_image *image = muse_image_new();
566  image->data = cpl_image_duplicate(cpl_imagelist_get(cube->data, iplane));
567  image->dq = cpl_image_duplicate(cpl_imagelist_get(cube->dq, iplane));
568  image->stat = cpl_image_duplicate(cpl_imagelist_get(cube->stat, iplane));
569  muse_imagelist_set(list, image, ilist++);
570  } /* for iplane (planes around ref. wavelength) */
571  muse_image *median = muse_combine_median_create(list);
572  /* fill existing empty header with the header of the *
573  * cube, but without the third dimension in the WCS */
574  cpl_propertylist_copy_property_regexp(median->header, cube->header,
575  "^C...*3$|^CD3_.$", 1);
576  muse_imagelist_delete(list);
577  if (!median) {
578  return cpl_error_set_message(__func__, cpl_error_get_code(),
579  "Failure while creating detection image!");
580  }
581  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
582  const char *fn = "wcs__image_median.fits";
583  cpl_msg_info(__func__, "Saving median detection image as \"%s\"", fn);
584  muse_image_save(median, fn);
585  }
586  cube->recimages = muse_imagelist_new();
587  cube->recnames = cpl_array_new(2, CPL_TYPE_STRING);
588  /* create a white-light image; this is unused here but should *
589  * be saved as first extension in the output cube */
590  muse_table *fwhite = muse_table_load_filter(NULL, "white");
591  muse_image *white = muse_datacube_collapse(cube, fwhite);
592  muse_table_delete(fwhite);
593  muse_imagelist_set(cube->recimages, white, 0);
594  cpl_array_set_string(cube->recnames, 0, "white");
595  /* append this median-combined detection image also as *
596  * part of the reconstructed images of the datacube */
597  muse_imagelist_set(cube->recimages, median, 1);
598  cpl_array_set_string(cube->recnames, 1, MUSE_WCS_DETIMAGE_EXTNAME);
599 
600  /* now do the centroiding of all detected sources */
601  cpl_table *detections = muse_wcs_centroid_stars(median, aSigma, aCentroid);
602  if (!detections || (detections && cpl_table_get_nrow(detections) < 1)) {
603  return cpl_error_get_code();
604  }
605  /* save pixel value and sky coordinates of center of the data */
606  aWCSObj->xcenter = cpl_image_get_size_x(median->data) / 2.,
607  aWCSObj->ycenter = cpl_image_get_size_y(median->data) / 2.;
608  aWCSObj->ra = muse_pfits_get_ra(median->header);
609  aWCSObj->dec = muse_pfits_get_dec(median->header);
610 #if 0
611  cpl_msg_debug(__func__, "image size: %d x %d --> center %f, %f",
612  (int)cpl_image_get_size_x(median->data),
613  (int)cpl_image_get_size_y(median->data),
614  aWCSObj->xcenter, aWCSObj->ycenter);
615 #endif
616 
617  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
618  const char *fn = "wcs__detections.fits";
619  cpl_msg_info(__func__, "Saving table of detections as \"%s\"", fn);
620  cpl_table_save(detections, NULL, NULL, fn, CPL_IO_CREATE);
621  }
622 
623  aWCSObj->detected = detections;
624  return CPL_ERROR_NONE;
625 } /* muse_wcs_locate_sources() */
626 
627 /*----------------------------------------------------------------------------*/
669 /*----------------------------------------------------------------------------*/
670 cpl_error_code
671 muse_wcs_solve(muse_wcs_object *aWCSObj, cpl_table *aReference,
672  float aRadius, float aFAccuracy, int aIter, float aThresh)
673 {
674  cpl_ensure_code(aWCSObj, CPL_ERROR_NULL_INPUT);
675  cpl_table *detected = aWCSObj->detected;
676  int ndet = cpl_table_get_nrow(detected),
677  nref = cpl_table_get_nrow(aReference);
678  cpl_ensure_code(ndet > 0 && nref > 0, CPL_ERROR_ILLEGAL_INPUT);
679  cpl_ensure_code(muse_cpltable_check(detected, muse_wcs_detections_def)
680  == CPL_ERROR_NONE &&
681  muse_cpltable_check(aReference, muse_wcs_reference_def)
682  == CPL_ERROR_NONE,
683  CPL_ERROR_BAD_FILE_FORMAT);
684  cpl_ensure_code(aRadius > 0.&& aFAccuracy > 0., CPL_ERROR_ILLEGAL_INPUT);
685 
686  /* sort tables by the source brightness */
687  cpl_propertylist *order = cpl_propertylist_new();
688  cpl_propertylist_append_bool(order, "FLUX", TRUE);
689  cpl_table_sort(detected, order);
690  cpl_propertylist_erase(order, "FLUX");
691  cpl_propertylist_append_bool(order, "mag", FALSE);
692  cpl_table_sort(aReference, order);
693  cpl_propertylist_delete(order);
694  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 2) {
695  FILE *fp = fopen("wcs__detections.ascii", "w");
696  fprintf(fp, "%s: table of detected sources (sorted by flux):\n", __func__);
697  cpl_table_dump(detected, 0, 1000, fp);
698  fclose(fp);
699  fp = fopen("wcs__references.ascii", "w");
700  fprintf(fp, "%s: table of reference objects (sorted by flux):\n", __func__);
701  cpl_table_dump(aReference, 0, 1000, fp);
702  fclose(fp);
703  }
704 
705  /* construct a basic input WCS */
706  cpl_propertylist *wcsin = muse_wcs_create_default();
707  cpl_propertylist_append_double(wcsin, "CRVAL1", aWCSObj->ra);
708  cpl_propertylist_append_double(wcsin, "CRVAL2", aWCSObj->dec);
709  cpl_propertylist_update_double(wcsin, "CRPIX1", aWCSObj->crpix1);
710  cpl_propertylist_update_double(wcsin, "CRPIX2", aWCSObj->crpix2);
711  /* add NAXIS to fool cpl_wcs_new_from_propertylist() */
712  cpl_propertylist_append_int(wcsin, "NAXIS", 2);
713  cpl_propertylist_append_int(wcsin, "NAXIS1", aWCSObj->xcenter * 2.);
714  cpl_propertylist_append_int(wcsin, "NAXIS2", aWCSObj->ycenter * 2.);
715 
716  /* convert input tables into matrices for the pattern-matching function */
717  cpl_matrix *data = cpl_matrix_new(2, ndet),
718  *patt = cpl_matrix_new(2, nref);
719  int i;
720  for (i = 0; i < ndet; i++) {
721  cpl_matrix_set(data, 0, i, cpl_table_get(detected, "XPOS", i, NULL));
722  cpl_matrix_set(data, 1, i, cpl_table_get(detected, "YPOS", i, NULL));
723  } /* for i (all detections) */
724 
725  /* use the basic WCS to transform input reference positions *
726  * to pixel positions, to take out deformations by gnomonic *
727  * projection before attempting pattern matching */
728  for (i = 0; i < nref; i++) {
729  double ra = cpl_table_get(aReference, "RA", i, NULL),
730  dec = cpl_table_get(aReference, "DEC", i, NULL),
731  x, y;
732  muse_wcs_pixel_from_celestial(wcsin, ra, dec, &x, &y);
733  cpl_matrix_set(patt, 0, i, x);
734  cpl_matrix_set(patt, 1, i, y);
735 #if 0
736  printf("conversion: %2d\t%.7f %.7f\t%6.2f %6.2f\n", i + 1, ra, dec, x, y);
737  fflush(stdout);
738 #endif
739  } /* for i (all reference points) */
740 #if 0
741  if (cpl_msg_get_level() == CPL_MSG_DEBUG) {
742  printf("%s: data matrix:\n", __func__);
743  cpl_matrix_dump(data, stdout);
744  printf("%s: patt matrix:\n", __func__);
745  cpl_matrix_dump(patt, stdout);
746  fflush(stdout);
747  }
748 #endif
749 
750  /* compute typical data error from the error columns, *
751  * to use this as input for the pattern matching */
752  double accuracy0 = sqrt(pow(cpl_table_get_column_mean(detected, "XERR"), 2) +
753  pow(cpl_table_get_column_mean(detected, "YERR"), 2));
754  /* start with somewhat worse accuracy to make *
755  * it more likely that a match is found at all */
756  double accuracy = accuracy0 * aFAccuracy,
757  radius = aRadius; /* start with (hopefully) large search radius */
758  int nid = INT_MAX; /* number of identified detections */
759  cpl_array *aid = NULL;
760  cpl_boolean dupes = CPL_FALSE;
761  double xscale, yscale;
762  do {
763  double scale, angle;
764  /* As recommended in the CPL docs, initially select the 20 first (brightest!) *
765  * detections and the 10 first (brightest) reference sources, so that *
766  * hopefully all reference sources are contained in the brightest detections. */
767 #define USE_DATA 15
768 #define USE_PATT 10
769  int ndata = ndet < USE_DATA ? ndet : USE_DATA,
770  npatt = nref < USE_PATT ? nref : USE_PATT;
771  cpl_array_delete(aid);
772  do {
773  cpl_msg_debug(__func__, "trying pattern matching with accuracy %g and "
774  "radius %g", accuracy, radius);
775  aid = cpl_ppm_match_points(data, ndata, accuracy,
776  patt, npatt, 1. /* [pix] */,
777  0.1 /* 10% */, radius, NULL, NULL,
778  &scale, &angle);
779  /* decrease accuracy in case pattern matching didn't succeed */
780  accuracy /= aid ? 1. : 1.5;
781  if (accuracy < FLT_EPSILON) {
782  break; /* doesn't make sense to go any lower */
783  }
784  } while (!aid); /* no matched points likely means to low accuracy */
785  cpl_errorstate state = cpl_errorstate_get();
786  nid = cpl_array_get_size(aid);
787  if (nid > 0) { /* subtract valid only if the array exists */
788  nid -= cpl_array_count_invalid(aid);
789  }
790 #if 0
791  printf("aid (valid=%d):\n", nid);
792  cpl_array_dump(aid, 0, cpl_array_get_size(aid), stdout);
793  fflush(stdout);
794 #endif
795  if (nid < 1) {
796  cpl_array_delete(aid);
797  cpl_matrix_delete(data);
798  cpl_matrix_delete(patt);
799  cpl_errorstate_set(state); /* swallow error about NULL cpl_array */
800  cpl_propertylist_delete(wcsin);
801  return cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND, "None of "
802  "the %d detections could be identified with "
803  "the %d reference positions with radius %.3f "
804  "pix (starting value %.3f) and data accuracy "
805  "%.3e pix (intrinsic mean error %.3e)", ndet,
806  nref, radius, aRadius, accuracy, accuracy0);
807  }
808  /* the first-guess scale computed here only makes sense *
809  * if multiplied by the pixel scale in the first-guess WCS */
810  muse_wcs_get_scales(wcsin, &xscale, &yscale);
811  dupes = muse_cplarray_has_duplicate(aid);
812  cpl_msg_debug(__func__, "%d %sidentified points [scale = %g, angle = %g; "
813  "used radius = %.3f pix (starting value %.3f), data accuracy "
814  "= %.3e pix (intrinsic mean error %.3e)]", nid,
815  dupes ? "(non-unique!) " : "unique ",
816  scale*1800.*(xscale+yscale), angle, radius, aRadius, accuracy,
817  accuracy0);
818  radius /= 1.25; /* next loop with smaller radius, if necessary */
819  } while (dupes);
820  cpl_matrix_delete(data);
821  cpl_matrix_delete(patt);
822 
823  /* create matrices again for cpl_wcs_platesol(), this *
824  * time with detected pixel positions and sky positions, *
825  * but only for the identified detections */
826  cpl_matrix *mpx = cpl_matrix_new(nid, 2),
827  *msky = cpl_matrix_new(nid, 2);
828  int iid = 0; /* index of identified object in matrices */
829  for (i = 0; i < cpl_array_get_size(aid); i++) { /* index in reference table */
830  if (!cpl_array_is_valid(aid, i)) {
831  continue;
832  }
833  int idata = cpl_array_get_int(aid, i, NULL); /* index in detections table */
834  cpl_matrix_set(mpx, iid, 0, cpl_table_get(detected, "XPOS", idata, NULL));
835  cpl_matrix_set(mpx, iid, 1, cpl_table_get(detected, "YPOS", idata, NULL));
836  cpl_matrix_set(msky, iid, 0, cpl_table_get(aReference, "RA", i, NULL));
837  cpl_matrix_set(msky, iid, 1, cpl_table_get(aReference, "DEC", i, NULL));
838 #if 0
839  printf("matrices: %2d\t%.7f %.7f\t%6.2f %6.2f\n", iid + 1,
840  cpl_table_get(aReference, "RA", i, NULL),
841  cpl_table_get(aReference, "DEC", i, NULL),
842  cpl_table_get(detected, "XPOS", idata, NULL),
843  cpl_table_get(detected, "YPOS", idata, NULL));
844 #endif
845  iid++;
846  }
847 #if 0
848  printf("mpx:\n");
849  cpl_matrix_dump(mpx, stdout);
850  printf("msky:\n");
851  cpl_matrix_dump(msky, stdout);
852  fflush(stdout);
853 #endif
854  cpl_array_delete(aid);
855 
856  /* compute the solution */
857  cpl_propertylist *wcsout = NULL;
858  cpl_error_code rc = cpl_wcs_platesol(wcsin, msky, mpx, aIter, aThresh,
859  CPL_WCS_PLATESOL_6, CPL_WCS_MV_CRVAL,
860  &wcsout);
861  if (aWCSObj->cube) {
862  cpl_propertylist_copy_property_regexp(wcsout, aWCSObj->cube->header,
863  CPL_WCS_REGEXP, 1);
864  } /* if cube */
865  cpl_matrix_delete(mpx);
866  cpl_matrix_delete(msky);
867  cpl_propertylist_delete(wcsin);
868  if (rc != CPL_ERROR_NONE) {
869  cpl_msg_warning(__func__, "Computing the WCS solution returned an error "
870  "(%d): %s", rc, cpl_error_get_message());
871  }
872 
873  /* Compute the rotation angle and scales */
874  double cd11 = muse_pfits_get_cd(wcsout, 1, 1),
875  cd22 = muse_pfits_get_cd(wcsout, 2, 2),
876  cd12 = muse_pfits_get_cd(wcsout, 1, 2),
877  cd21 = muse_pfits_get_cd(wcsout, 2, 1),
878  det = cd11 * cd22 - cd12 * cd21;
879  if (det < 0.) {
880  cd12 *= -1;
881  cd11 *= -1;
882  }
883  /* the values we want, by default for non-rotation */
884  double xang, yang;
885  muse_wcs_get_angles(wcsout, &xang, &yang);
886  muse_wcs_get_scales(wcsout, &xscale, &yscale);
887  xscale *= 3600.; /* scales in arc seconds */
888  yscale *= 3600.;
889  cpl_msg_info(__func__, "astrometric calibration results: scales %f/%f "
890  "arcsec/spaxel, rotation %g/%g deg", xscale, yscale, xang, yang);
891 
892 #if 0
893  printf("%s: output propertylist (rc = %d):\n", __func__, rc);
894  fflush(stdout);
895  cpl_propertylist_save(wcsout, "astrometry_wcsout.fits", CPL_IO_CREATE);
896  system("fold -w 80 astrometry_wcsout.fits | grep -v \"^ \"");
897  remove("astrometry_wcsout.fits");
898 #endif
899 
900  /* number of stars input to the astrometric fit as QC */
901  cpl_propertylist_update_int(wcsout, QC_ASTROMETRY_NSTARS, nid);
902  /* scales for QC in arcsec */
903  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_SCX, xscale);
904  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_SCY, yscale);
905  /* angles in degrees */
906  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_ANGX, xang);
907  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_ANGY, yang);
908  /* copy the "systematic error" as residuals for QC */
909  double medresx = cpl_propertylist_get_double(wcsout, "CSYER1") * 3600.,
910  medresy = cpl_propertylist_get_double(wcsout, "CSYER2") * 3600.;
911  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_RESX, medresx);
912  cpl_propertylist_update_float(wcsout, QC_ASTROMETRY_RESY, medresy);
913 
914  /* update header with the final values used, just for information */
915  cpl_propertylist_update_float(wcsout, MUSE_HDR_WCS_RADIUS, radius);
916  cpl_propertylist_set_comment(wcsout, MUSE_HDR_WCS_RADIUS,
917  MUSE_HDR_WCS_RADIUS_C);
918  cpl_propertylist_update_float(wcsout, MUSE_HDR_WCS_ACCURACY, accuracy);
919  cpl_propertylist_set_comment(wcsout, MUSE_HDR_WCS_ACCURACY,
920  MUSE_HDR_WCS_ACCURACY_C);
921  cpl_propertylist_update_float(wcsout, MUSE_HDR_WCS_FACCURACY,
922  accuracy / accuracy0);
923  cpl_propertylist_set_comment(wcsout, MUSE_HDR_WCS_FACCURACY,
924  MUSE_HDR_WCS_FACCURACY_C);
925 
926  aWCSObj->wcs = wcsout;
927  return rc;
928 } /* muse_wcs_solve() */
929 
930 /*----------------------------------------------------------------------------*/
986 /*----------------------------------------------------------------------------*/
987 cpl_error_code
988 muse_wcs_optimize_solution(muse_wcs_object *aWCSObj, float aDetSigma,
989  muse_wcs_centroid_type aCentroid,
990  cpl_table *aReference, float aRadius,
991  float aFAccuracy, int aIter, float aRejSigma)
992 {
993  cpl_ensure_code(aWCSObj && aWCSObj->cube, CPL_ERROR_NULL_INPUT);
994  /* check that the cube contains the detection image */
995  cpl_ensure_code(!strncmp(cpl_array_get_string(aWCSObj->cube->recnames, 1),
997  CPL_ERROR_DATA_NOT_FOUND);
998  cpl_ensure_code(aDetSigma < 0., CPL_ERROR_ILLEGAL_INPUT);
999  switch (aCentroid) {
1002  case MUSE_WCS_CENTROID_BOX:
1003  break;
1004  default:
1005  return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1006  }
1007  int nref = cpl_table_get_nrow(aReference);
1008  cpl_ensure_code(nref > 0, CPL_ERROR_ILLEGAL_INPUT);
1009  cpl_ensure_code(muse_cpltable_check(aReference, muse_wcs_reference_def)
1010  == CPL_ERROR_NONE, CPL_ERROR_BAD_FILE_FORMAT);
1011  cpl_ensure_code(aRadius > 0.&& aFAccuracy > 0., CPL_ERROR_ILLEGAL_INPUT);
1012 
1013  float detsigma = fabsf(aDetSigma),
1014  lolimit = 0.9999; /* 1.0 with a bit of margin */
1015  muse_image *detimage = muse_imagelist_get(aWCSObj->cube->recimages, 1);
1016 
1017  /* clear incoming detections and solutions (if any) */
1018  cpl_table_delete(aWCSObj->detected);
1019  aWCSObj->detected = NULL;
1020  cpl_propertylist_delete(aWCSObj->wcs);
1021  aWCSObj->wcs = NULL;
1022  /* create table that can be used to find the optimal solution afterwards */
1023  cpl_table *tres = cpl_table_new(lround((detsigma - lolimit) / 0.1) + 1);
1024  cpl_table_new_column(tres, "detsigma", CPL_TYPE_FLOAT);
1025  cpl_table_set_column_format(tres, "detsigma", "%.3f");
1026  cpl_table_new_column(tres, "ndet", CPL_TYPE_INT);
1027  cpl_table_new_column(tres, "nid", CPL_TYPE_INT);
1028  cpl_table_new_column(tres, "scalex", CPL_TYPE_FLOAT);
1029  cpl_table_set_column_format(tres, "scalex", "%.4f");
1030  cpl_table_new_column(tres, "scaley", CPL_TYPE_FLOAT);
1031  cpl_table_set_column_format(tres, "scaley", "%.4f");
1032  cpl_table_new_column(tres, "anglex", CPL_TYPE_FLOAT);
1033  cpl_table_set_column_format(tres, "anglex", "%.3f");
1034  cpl_table_new_column(tres, "angley", CPL_TYPE_FLOAT);
1035  cpl_table_set_column_format(tres, "angley", "%.3f");
1036  cpl_table_new_column(tres, "medresx", CPL_TYPE_FLOAT);
1037  cpl_table_set_column_format(tres, "medresx", "%.3f");
1038  cpl_table_new_column(tres, "medresy", CPL_TYPE_FLOAT);
1039  cpl_table_set_column_format(tres, "medresy", "%.3f");
1040 
1041  /* run detection and solution for all sigma levels between |detsigma| and 1.0 */
1042  cpl_error_code rc = CPL_ERROR_NONE;
1043  float ds;
1044  int irow;
1045  for (ds = detsigma, irow = 0; ds > lolimit; ds -= 0.1, irow++) {
1046  cpl_msg_debug(__func__, "searching for solution with detection sigma %.3f", ds);
1047  cpl_msg_indent_more();
1048  /* quieten output from cpl_ppm_match_points(),d muse_wcs_centroid_stars(), *
1049  * etc.; we'll get the necessary output at the end with the final run */
1050  cpl_msg_severity level = cpl_msg_get_level();
1051  cpl_msg_set_level(CPL_MSG_WARNING);
1052 
1053  /* first, detect stars at this level */
1054  aWCSObj->detected = muse_wcs_centroid_stars(detimage, ds, aCentroid);
1055  cpl_table_set_float(tres, "detsigma", irow, ds);
1056  cpl_table_set_int(tres, "ndet", irow, cpl_table_get_nrow(aWCSObj->detected));
1057 
1058  /* now see, if we identify stars and find a solution */
1059  cpl_errorstate state = cpl_errorstate_get();
1060  rc = muse_wcs_solve(aWCSObj, aReference, aRadius, aFAccuracy, aIter,
1061  aRejSigma);
1062  if (rc == CPL_ERROR_NONE && aWCSObj->wcs) {
1063  cpl_table_set_int(tres, "nid", irow,
1064  cpl_propertylist_get_int(aWCSObj->wcs, QC_ASTROMETRY_NSTARS));
1065  cpl_table_set_float(tres, "scalex", irow,
1066  cpl_propertylist_get_float(aWCSObj->wcs, QC_ASTROMETRY_SCX));
1067  cpl_table_set_float(tres, "scaley", irow,
1068  cpl_propertylist_get_float(aWCSObj->wcs, QC_ASTROMETRY_SCY));
1069  cpl_table_set_float(tres, "anglex", irow,
1070  cpl_propertylist_get_float(aWCSObj->wcs, QC_ASTROMETRY_ANGX));
1071  cpl_table_set_float(tres, "angley", irow,
1072  cpl_propertylist_get_float(aWCSObj->wcs, QC_ASTROMETRY_ANGY));
1073  cpl_table_set_float(tres, "medresx", irow,
1074  cpl_propertylist_get_float(aWCSObj->wcs, QC_ASTROMETRY_RESX));
1075  cpl_table_set_float(tres, "medresy", irow,
1076  cpl_propertylist_get_float(aWCSObj->wcs, QC_ASTROMETRY_RESY));
1077  cpl_propertylist_delete(aWCSObj->wcs);
1078  aWCSObj->wcs = NULL;
1079  } else {
1080  cpl_errorstate_set(state); /* ignore errors at this point */
1081  }
1082  cpl_table_delete(aWCSObj->detected);
1083  aWCSObj->detected = NULL;
1084  cpl_msg_set_level(level);
1085  cpl_msg_indent_less();
1086  } /* for ds */
1087 
1088  cpl_boolean debug = getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) >= 1;
1089  if (debug) {
1090  printf("%s: full table of results:\n", __func__);
1091  cpl_table_dump(tres, 0, 1000, stdout);
1092  fflush(stdout);
1093  }
1094 
1095  /* reject failed results and those far away from the expected scale */
1096  double scale = muse_pfits_get_mode(aWCSObj->cube->header) == MUSE_MODE_NFM_AO_N
1097  ? kMuseSpaxelSizeX_NFM : kMuseSpaxelSizeY_WFM;
1098  cpl_msg_debug(__func__, "pruning results +/-10%% away from expected spaxel "
1099  "scale of %.3f arcsec/pixel", scale);
1100  cpl_table_unselect_all(tres);
1101  cpl_table_or_selected_float(tres, "scalex", CPL_GREATER_THAN, scale * 1.1);
1102  cpl_table_or_selected_float(tres, "scalex", CPL_LESS_THAN, scale * 0.9);
1103  cpl_table_or_selected_float(tres, "scaley", CPL_GREATER_THAN, scale * 1.1);
1104  cpl_table_or_selected_float(tres, "scaley", CPL_LESS_THAN, scale * 0.9);
1105  cpl_table_or_selected_invalid(tres, "nid"); /* any of the result columns ... */
1106  cpl_table_erase_selected(tres);
1107  if (debug) {
1108  printf("%s: pruned table of results:\n", __func__);
1109  cpl_table_dump(tres, 0, 1000, stdout);
1110  fflush(stdout);
1111  }
1112  /* check if there are still table rows left */
1113  if (cpl_table_get_nrow(tres) < 1) {
1114  cpl_table_delete(tres);
1115  cpl_msg_error(__func__, "No valid solution found in the range %.3f .. %.3f "
1116  "sigma", detsigma, lolimit);
1117  return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_OUTPUT);
1118  }
1119 
1120  /* add table column with a weight to evaluate the goodness of each result; *
1121  * use *
1122  * NID / 50. * min(medresx) / medresx * min(medresy) / medresy *
1123  * as formula to fill it */
1124  cpl_table_cast_column(tres, "nid", "weight", CPL_TYPE_DOUBLE);
1125  cpl_table_set_column_format(tres, "weight", "%.5g");
1126  cpl_table_divide_scalar(tres, "weight", 50.);
1127  cpl_table_divide_columns(tres, "weight", "medresx");
1128  cpl_table_multiply_scalar(tres, "weight",
1129  cpl_table_get_column_min(tres, "medresx"));
1130  cpl_table_divide_columns(tres, "weight", "medresy");
1131  cpl_table_multiply_scalar(tres, "weight",
1132  cpl_table_get_column_min(tres, "medresy"));
1133  cpl_propertylist *sort = cpl_propertylist_new();
1134  cpl_propertylist_append_bool(sort, "weight", CPL_TRUE);
1135  cpl_propertylist_append_bool(sort, "nid", CPL_TRUE);
1136  cpl_table_sort(tres, sort);
1137  cpl_propertylist_delete(sort);
1138  /* the best detection sigma is now the one at the top of the table */
1139  detsigma = cpl_table_get_float(tres, "detsigma", 0, NULL);
1140  if (debug) {
1141  printf("%s: pruned and sorted table of results:\n", __func__);
1142  cpl_table_dump(tres, 0, 1000, stdout);
1143  printf("%s: ===> use the %.3f-sigma level\n", __func__, detsigma);
1144  fflush(stdout);
1145  }
1146  cpl_table_delete(tres);
1147 
1148  /* so, run detection and solution one last time */
1149  aWCSObj->detected = muse_wcs_centroid_stars(detimage, detsigma, aCentroid);
1150  rc = muse_wcs_solve(aWCSObj, aReference, aRadius, aFAccuracy, aIter, aRejSigma);
1151 
1152  /* update header with the final value used, for information */
1153  if (aWCSObj->wcs) {
1154  cpl_propertylist_update_float(aWCSObj->wcs, MUSE_HDR_WCS_DETSIGMA, detsigma);
1155  cpl_propertylist_set_comment(aWCSObj->wcs, MUSE_HDR_WCS_DETSIGMA,
1156  MUSE_HDR_WCS_DETSIGMA_C);
1157  }
1158 
1159  return rc;
1160 } /* muse_wcs_optimize_solution() */
1161 
1162 /*----------------------------------------------------------------------------*/
1172 /*----------------------------------------------------------------------------*/
1173 cpl_propertylist *
1175 {
1176  cpl_propertylist *wcs = cpl_propertylist_new();
1177 
1178  /* We only care about the spatial axes, pretend that we have a 2D image; *
1179  * set WCSAXES to make fitsverify happy when it finds the other ^C headers. */
1180  cpl_propertylist_append_int(wcs, "WCSAXES", 2);
1181 
1182  /* Check, if we are dealing with wcslib >= 4.23, which has the fix with the *
1183  * floating point numbers. If an older version is found, use a small number *
1184  * instead of zero for CRPIX, so that on loading CPL sees it as float value. */
1185  double smallvalue = FLT_MIN;
1186  const char *cpldesc = cpl_get_description(CPL_DESCRIPTION_DEFAULT); /* never fails */
1187  /* search for the WCSLIB version string */
1188  char *pcpldesc = strstr(cpldesc, "WCSLIB = ");
1189  if (pcpldesc) {
1190  pcpldesc += 8;
1191  double wcslibversion = atof(pcpldesc);
1192  if (wcslibversion >= 4.23) {
1193  smallvalue = 0.;
1194  } /* if newer wcslib */
1195  } /* if wcslib string found */
1196 
1197  /* axis 1 */
1198  /* CRPIX is the zero order correction to the astrometry, set zero here */
1199  cpl_propertylist_append_double(wcs, "CRPIX1", smallvalue);
1200  /* negative value, to signify that east is to the left */
1201  cpl_propertylist_append_double(wcs, "CD1_1", -kMuseSpaxelSizeX_WFM / 3600);
1202  cpl_propertylist_append_string(wcs, "CTYPE1", "RA---TAN");
1203  cpl_propertylist_append_string(wcs, "CUNIT1", "deg");
1204 
1205  /* axis 2 */
1206  cpl_propertylist_append_double(wcs, "CRPIX2", smallvalue);
1207  cpl_propertylist_append_double(wcs, "CD2_2", kMuseSpaxelSizeY_WFM / 3600);
1208  cpl_propertylist_append_string(wcs, "CTYPE2", "DEC--TAN");
1209  cpl_propertylist_append_string(wcs, "CUNIT2", "deg");
1210 
1211  /* cross-terms are assumed to be not present (no rotation known!) */
1212  cpl_propertylist_append_double(wcs, "CD1_2", 0.);
1213  cpl_propertylist_append_double(wcs, "CD2_1", 0.);
1214 
1215  /* leave out CRVAL1/2, as we use this only to do the gnomonic projection */
1216 
1217  return wcs;
1218 } /* muse_wcs_create_default() */
1219 
1220 /*----------------------------------------------------------------------------*/
1235 /*----------------------------------------------------------------------------*/
1236 cpl_propertylist *
1237 muse_wcs_apply_cd(const cpl_propertylist *aExpHeader,
1238  const cpl_propertylist *aCalHeader)
1239 {
1240  cpl_ensure(aCalHeader && aExpHeader, CPL_ERROR_NULL_INPUT, NULL);
1241 
1242  /* duplicate the input WCS because we want to adapt it to the current data */
1243  cpl_propertylist *header = cpl_propertylist_duplicate(aCalHeader);
1244 
1245  /* Find field rotation (in radians) and create a CD matrix that reflects *
1246  * the exposure position angle corrected by the astrometric solution. */
1247  double pa = muse_astro_posangle(aExpHeader) * CPL_MATH_RAD_DEG,
1248  xang, yang, xsc, ysc;
1249  muse_wcs_get_scales(header, &xsc, &ysc);
1250  muse_wcs_get_angles(header, &xang, &yang);
1251  cpl_msg_debug(__func__, "WCS solution: scales %f / %f arcsec, angles %f / %f "
1252  "deg", xsc * 3600., ysc * 3600., xang, yang);
1253  /* create PCi_j matrix from the CDi_j in the header and invert it */
1254  cpl_matrix *pc = cpl_matrix_new(2, 2);
1255  cpl_matrix_set(pc, 0, 0, muse_pfits_get_cd(header, 1, 1) / xsc);
1256  cpl_matrix_set(pc, 0, 1, muse_pfits_get_cd(header, 1, 2) / xsc);
1257  cpl_matrix_set(pc, 1, 0, muse_pfits_get_cd(header, 2, 1) / ysc);
1258  cpl_matrix_set(pc, 1, 1, muse_pfits_get_cd(header, 2, 2) / ysc);
1259  cpl_matrix *pcinv = cpl_matrix_invert_create(pc);
1260  cpl_matrix_delete(pc);
1261  /* now create corrective CDi_j elements from the inverted PCi_j matrix */
1262  double cd11cor, cd12cor, cd21cor, cd22cor;
1263  if (pcinv) {
1264  cd11cor = xsc * cpl_matrix_get(pcinv, 0, 0);
1265  cd12cor = xsc * cpl_matrix_get(pcinv, 0, 1);
1266  cd21cor = ysc * cpl_matrix_get(pcinv, 1, 0);
1267  cd22cor = ysc * cpl_matrix_get(pcinv, 1, 1);
1268  cpl_matrix_delete(pcinv);
1269  } else {
1270  cpl_msg_warning(__func__, "CD matrix of astrometric solution could not "
1271  "be inverted! Assuming negligible zeropoint rotation.");
1272  cd11cor = xsc * 1.;
1273  cd12cor = xsc * 0.;
1274  cd21cor = ysc * 0.;
1275  cd22cor = ysc * 1.;
1276  }
1277  /* now we can finally compute the effective CDi_j of the exposure */
1278  double cd11 = cos(pa) * cd11cor - sin(pa) * cd21cor,
1279  cd12 = cos(pa) * cd12cor - sin(pa) * cd22cor,
1280  cd21 = sin(pa) * cd11cor + cos(pa) * cd21cor,
1281  cd22 = sin(pa) * cd12cor + cos(pa) * cd22cor;
1282  cpl_propertylist_update_double(header, "CD1_1", cd11),
1283  cpl_propertylist_update_double(header, "CD1_2", cd12),
1284  cpl_propertylist_update_double(header, "CD2_1", cd21);
1285  cpl_propertylist_update_double(header, "CD2_2", cd22),
1286  muse_wcs_get_scales(header, &xsc, &ysc);
1287  muse_wcs_get_angles(header, &xang, &yang);
1288  cpl_msg_debug(__func__, "Updated CD matrix (%f deg field rotation): "
1289  "%e %e %e %e (scales %f / %f arcsec, angles %f / %f deg)",
1290  pa * CPL_MATH_DEG_RAD, cd11, cd12, cd21, cd22,
1291  xsc * 3600., ysc * 3600., xang, yang);
1292  return header;
1293 } /* muse_wcs_apply_cd() */
1294 
1295 /*----------------------------------------------------------------------------*/
1324 /*----------------------------------------------------------------------------*/
1325 cpl_error_code
1326 muse_wcs_project_tan(muse_pixtable *aPixtable, const cpl_propertylist *aWCS)
1327 {
1328  cpl_size nrow = muse_pixtable_get_nrow(aPixtable);
1329  /* nrow == 0 implies a NULL or broken pixel table */
1330  cpl_ensure_code(nrow > 0 && aPixtable->header && aWCS, CPL_ERROR_NULL_INPUT);
1331  cpl_ensure_code(muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_PIXEL,
1332  CPL_ERROR_INCOMPATIBLE_INPUT);
1333  /* ensure that the input WCS is for gnomonic projection */
1334  const char *type1 = muse_pfits_get_ctype(aWCS, 1),
1335  *type2 = muse_pfits_get_ctype(aWCS, 2);
1336  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1337  !strncmp(type2, "DEC--TAN", 9),
1338  CPL_ERROR_UNSUPPORTED_MODE);
1339 
1340  /* clean main WCS keys from input pixel table, in *
1341  * case there are keys we do not overwrite below */
1342  cpl_propertylist_erase_regexp(aPixtable->header, MUSE_WCS_KEYS, 0);
1343 
1344  /* apply the exposure rotation on top of the zeropoint *
1345  * rotation from the astrometric calibration */
1346  cpl_propertylist *header = muse_wcs_apply_cd(aPixtable->header, aWCS);
1347  /* don't want CRVAL or LON/LATPOLE in the output, *
1348  * because we create native coords here */
1349  cpl_propertylist_erase_regexp(header, "^CRVAL[0-9]+$|^L[OA][NT]POLE$", 0);
1350  double cd11 = muse_pfits_get_cd(header, 1, 1),
1351  cd12 = muse_pfits_get_cd(header, 1, 2),
1352  cd21 = muse_pfits_get_cd(header, 2, 1),
1353  cd22 = muse_pfits_get_cd(header, 2, 2);
1354 
1355  /* Compute reference pixel from center of the data plus the zero *
1356  * order correction of the corrective WCS (CRPIX1/2 of aWCS); use *
1357  * the spatial extents before the DAR correction, if possible. */
1358  cpl_errorstate prestate = cpl_errorstate_get();
1359  double xlo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO),
1360  xhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI),
1361  ylo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO),
1362  yhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI);
1363  if (!cpl_errorstate_is_equal(prestate)) {
1364  cpl_errorstate_set(prestate);
1365  /* try normal pixel table headers now */
1366  xlo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO);
1367  xhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI);
1368  ylo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO);
1369  yhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI);
1370  }
1371  double wcspix1 = muse_pfits_get_crpix(header, 1),
1372  wcspix2 = muse_pfits_get_crpix(header, 2),
1373  crpix1 = (xhi + xlo) / 2. + wcspix1,
1374  crpix2 = (yhi + ylo) / 2. + wcspix2;
1375  cpl_propertylist_update_double(header, "CRPIX1", crpix1),
1376  cpl_propertylist_update_double(header, "CRPIX2", crpix2),
1377  cpl_msg_debug(__func__, "Using reference pixel %f/%f (limits in pixel table "
1378  "%f..%f/%f..%f, WCS correction %f,%f)", crpix1, crpix2,
1379  xlo, xhi, ylo, yhi, wcspix1, wcspix2);
1380 
1381  /* delete the units of the x/y columns while we are working on them */
1382  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "");
1383  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "");
1384 
1385  /* replace coordinate values in the pixel table by the computed ones, *
1386  * carry out the _partial_ WCS coordinate transform for each pixel */
1387  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
1388  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS);
1389  cpl_size n;
1390  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1391  shared(cd11, cd12, cd21, cd22, crpix1, crpix2, nrow, xpos, ypos)
1392  for (n = 0; n < nrow; n++) {
1393  /* conversion from pixel coordinates to projection plane coordinates; *
1394  * x,y as in Calabretta & Greisen 2002 A&A 395, 1077 (Paper II), Fig. 1 */
1395  double x = cd11 * (xpos[n] - crpix1) + cd12 * (ypos[n] - crpix2),
1396  y = cd22 * (ypos[n] - crpix2) + cd21 * (xpos[n] - crpix1);
1397 
1398  /* conversion from projection plane to native spherical coordinates: *
1399  * phi,theta as in Calabretta & Greisen 2002 (Paper II), Fig. 1; *
1400  * these formulae are for the gnomomic (TAN) case, Sect. 5.1/5.1.3, *
1401  * i.e. Eq. (14), (15), and (55) *
1402  * As we only further use these in other parts of the pipeline, the *
1403  * values are not converted to degrees but stay in radians. */
1404  double phi = atan2(x, -y),
1405  theta = atan(CPL_MATH_DEG_RAD / sqrt(x*x + y*y));
1406  if (phi < 0) { /* phi should be between 0 and 2pi, to let tests pass */
1407  phi += CPL_MATH_2PI;
1408  }
1409 
1410  /* conversion from native spherical to celestial spherical coordinates *
1411  * is done later when combining exposures, see muse_xcombine_tables() */
1412  xpos[n] = phi;
1413  ypos[n] = theta - CPL_MATH_PI_2; /* subtract pi/2 for better accuracy */
1414  } /* for n (table/matrix rows) */
1415 
1416  /* Here we produce units in radians; use "rad" unit string to signify *
1417  * that these are native spherical coordinates but haven't gotten full *
1418  * WCS treatment to alpha,delta yet */
1419  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "rad");
1420  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "rad");
1421  muse_pixtable_compute_limits(aPixtable);
1422  /* copy spatial WCS to the pixel table */
1423  cpl_propertylist_copy_property_regexp(aPixtable->header, header,
1424  MUSE_WCS_KEYS, 0);
1425  cpl_propertylist_delete(header);
1426 
1427  /* add the status header */
1428  cpl_propertylist_update_string(aPixtable->header, MUSE_HDR_PT_WCS,
1430  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_WCS,
1431  MUSE_HDR_PT_WCS_COMMENT_PROJ);
1432  return CPL_ERROR_NONE;
1433 } /* muse_wcs_project_tan() */
1434 
1435 /*----------------------------------------------------------------------------*/
1464 /*----------------------------------------------------------------------------*/
1465 cpl_error_code
1466 muse_wcs_position_celestial(muse_pixtable *aPixtable, double aRA, double aDEC)
1467 {
1468  cpl_size nrow = muse_pixtable_get_nrow(aPixtable);
1469  /* nrow == 0 implies a NULL or broken pixel table */
1470  cpl_ensure_code(nrow > 0 && aPixtable->header, CPL_ERROR_NULL_INPUT);
1471  cpl_ensure_code(muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_NATSPH,
1472  CPL_ERROR_INCOMPATIBLE_INPUT);
1473  /* ensure that the input WCS is for gnomonic projection */
1474  const char *type1 = muse_pfits_get_ctype(aPixtable->header, 1),
1475  *type2 = muse_pfits_get_ctype(aPixtable->header, 2);
1476  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1477  !strncmp(type2, "DEC--TAN", 9),
1478  CPL_ERROR_UNSUPPORTED_MODE);
1479 
1480  cpl_msg_info(__func__, "Adapting WCS to RA/DEC=%f/%f deg", aRA, aDEC);
1481 
1482  /* delete the units of the x/y columns while we are working on them */
1483  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "");
1484  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "");
1485 
1486  /* replace coordinate values in the pixel table by the computed ones, *
1487  * carry out the spherical coordinate rotation for each pixel */
1488  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
1489  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS);
1490  double dp = aDEC / CPL_MATH_DEG_RAD; /* delta_p in Paper II (in radians) */
1491  cpl_size n;
1492  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1493  shared(aDEC, dp, nrow, xpos, ypos)
1494  for (n = 0; n < nrow; n++) {
1495  /* conversion from native spherical to celestial spherical coordinates *
1496  * alpha,delta as in Calabretta & Greisen 2002 A&A 395, 1077 (Paper II), *
1497  * Fig. 1; the formulae used are Eq. (2) in Sect. 2.3 in that paper and *
1498  * are generic but use that phi_p = 180 deg for zenithal projections *
1499  * (like TAN), i.e. cos(phi - phi_p) = -cos(phi) and similar for sin. */
1500  double phi = xpos[n],
1501  theta = ypos[n] + CPL_MATH_PI_2, /* add pi/2 again */
1502  ra = atan2(cos(theta) * sin(phi),
1503  sin(theta) * cos(dp) + cos(theta) * sin(dp) * cos(phi))
1504  * CPL_MATH_DEG_RAD,
1505  dec = asin(sin(theta) * sin(dp) - cos(theta) * cos(dp) * cos(phi))
1506  * CPL_MATH_DEG_RAD;
1507  /* The following should be *
1508  * xpos = aRA + ra; *
1509  * ypos = dec; *
1510  * but let's remove the zeropoint, to help accuracy despite *
1511  * the limited float precision (~7 digits, i.e. ~0.36 arcsec). */
1512  xpos[n] = ra;
1513  ypos[n] = dec - aDEC;
1514  } /* for n (table/matrix rows) */
1515 
1516  /* We produce output units in degrees (like wcslib); signify *
1517  * full WCS transformation by setting the "deg" output unit */
1518  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS, "deg");
1519  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_YPOS, "deg");
1520  /* append CRVAL now that we have a full WCS */
1521  cpl_propertylist_update_double(aPixtable->header, "CRVAL1", aRA);
1522  cpl_propertylist_update_double(aPixtable->header, "CRVAL2", aDEC);
1523  /* do the limit recomputation after setting the CRVALs */
1524  muse_pixtable_compute_limits(aPixtable);
1525 
1526  /* add the status header */
1527  cpl_propertylist_update_string(aPixtable->header, MUSE_HDR_PT_WCS,
1529  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_WCS,
1530  MUSE_HDR_PT_WCS_COMMENT_POSI);
1531  return CPL_ERROR_NONE;
1532 } /* muse_wcs_position_celestial() */
1533 
1534 /*----------------------------------------------------------------------------*/
1557 /*----------------------------------------------------------------------------*/
1558 cpl_error_code
1559 muse_wcs_celestial_from_pixel(cpl_propertylist *aHeader, double aX, double aY,
1560  double *aRA, double *aDEC)
1561 {
1562  cpl_ensure_code(aHeader && aRA && aDEC, CPL_ERROR_NULL_INPUT);
1563  const char *type1 = muse_pfits_get_ctype(aHeader, 1),
1564  *type2 = muse_pfits_get_ctype(aHeader, 2);
1565  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1566  !strncmp(type2, "DEC--TAN", 9),
1567  CPL_ERROR_UNSUPPORTED_MODE);
1568 
1569  muse_wcs *wcs = muse_wcs_new(aHeader);
1570  muse_wcs_celestial_from_pixel_fast(wcs, aX, aY, aRA, aDEC);
1571  cpl_free(wcs);
1572 
1573  return CPL_ERROR_NONE;
1574 } /* muse_wcs_celestial_from_pixel() */
1575 
1576 /*----------------------------------------------------------------------------*/
1601 /*----------------------------------------------------------------------------*/
1602 cpl_error_code
1603 muse_wcs_pixel_from_celestial(cpl_propertylist *aHeader, double aRA, double aDEC,
1604  double *aX, double *aY)
1605 {
1606  cpl_ensure_code(aHeader && aX && aY, CPL_ERROR_NULL_INPUT);
1607  /* make sure that the header represents TAN */
1608  const char *type1 = muse_pfits_get_ctype(aHeader, 1),
1609  *type2 = muse_pfits_get_ctype(aHeader, 2);
1610  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1611  !strncmp(type2, "DEC--TAN", 9),
1612  CPL_ERROR_UNSUPPORTED_MODE);
1613 
1614  muse_wcs *wcs = muse_wcs_new(aHeader);
1615  if (wcs->cddet == 0.) { /* that's important here */
1616  *aX = *aY = NAN;
1617  cpl_error_set(__func__, CPL_ERROR_SINGULAR_MATRIX);
1618  cpl_free(wcs);
1619  return CPL_ERROR_SINGULAR_MATRIX;
1620  }
1621  wcs->crval1 /= CPL_MATH_DEG_RAD; /* convert to radians */
1622  wcs->crval2 /= CPL_MATH_DEG_RAD;
1623  muse_wcs_pixel_from_celestial_fast(wcs, aRA / CPL_MATH_DEG_RAD,
1624  aDEC / CPL_MATH_DEG_RAD, aX, aY);
1625  cpl_free(wcs);
1626 
1627  return CPL_ERROR_NONE;
1628 } /* muse_wcs_pixel_from_celestial() */
1629 
1630 /*----------------------------------------------------------------------------*/
1656 /*----------------------------------------------------------------------------*/
1657 cpl_error_code
1658 muse_wcs_projplane_from_celestial(cpl_propertylist *aHeader, double aRA,
1659  double aDEC, double *aX, double *aY)
1660 {
1661  cpl_ensure_code(aHeader && aX && aY, CPL_ERROR_NULL_INPUT);
1662  /* make sure that the header represents TAN */
1663  const char *type1 = muse_pfits_get_ctype(aHeader, 1),
1664  *type2 = muse_pfits_get_ctype(aHeader, 2);
1665  cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
1666  !strncmp(type2, "DEC--TAN", 9),
1667  CPL_ERROR_UNSUPPORTED_MODE);
1668 
1669  /* spherical coordinate shift/translation */
1670  double a = aRA / CPL_MATH_DEG_RAD, /* RA in radians */
1671  d = aDEC / CPL_MATH_DEG_RAD, /* DEC in radians */
1672  /* alpha_p and delta_p in Paper II (in radians) */
1673  ap = muse_pfits_get_crval(aHeader, 1) / CPL_MATH_DEG_RAD,
1674  dp = muse_pfits_get_crval(aHeader, 2) / CPL_MATH_DEG_RAD,
1675  phi = atan2(-cos(d) * sin(a - ap),
1676  sin(d) * cos(dp) - cos(d) * sin(dp) * cos(a-ap))
1677  + 180 / CPL_MATH_DEG_RAD,
1678  theta = asin(sin(d) * sin(dp) + cos(d) * cos(dp) * cos(a-ap)),
1679  R_theta = CPL_MATH_DEG_RAD / tan(theta);
1680  /* spherical deprojection */
1681  *aX = R_theta * sin(phi),
1682  *aY = -R_theta * cos(phi);
1683 
1684  return CPL_ERROR_NONE;
1685 } /* muse_wcs_projplane_from_celestial() */
1686 
1687 /*----------------------------------------------------------------------------*/
1704 /*----------------------------------------------------------------------------*/
1705 cpl_error_code
1706 muse_wcs_projplane_from_pixel(cpl_propertylist *aHeader, double aX, double aY,
1707  double *aXOut, double *aYOut)
1708 {
1709  cpl_ensure_code(aHeader && aXOut && aYOut, CPL_ERROR_NULL_INPUT);
1710 
1711  muse_wcs *wcs = muse_wcs_new(aHeader);
1712  muse_wcs_projplane_from_pixel_fast(wcs, aX, aY, aXOut, aYOut);
1713  cpl_free(wcs);
1714 
1715  return CPL_ERROR_NONE;
1716 } /* muse_wcs_projplane_from_pixel() */
1717 
1718 /*----------------------------------------------------------------------------*/
1737 /*----------------------------------------------------------------------------*/
1738 cpl_error_code
1739 muse_wcs_pixel_from_projplane(cpl_propertylist *aHeader, double aX, double aY,
1740  double *aXOut, double *aYOut)
1741 {
1742  cpl_ensure_code(aHeader && aXOut && aYOut, CPL_ERROR_NULL_INPUT);
1743 
1744  muse_wcs *wcs = muse_wcs_new(aHeader);
1745  if (wcs->cddet == 0.) { /* that's important here */
1746  *aXOut = *aYOut = NAN;
1747  cpl_error_set(__func__, CPL_ERROR_SINGULAR_MATRIX);
1748  cpl_free(wcs);
1749  return CPL_ERROR_SINGULAR_MATRIX;
1750  }
1751  muse_wcs_pixel_from_projplane_fast(wcs, aX, aY, aXOut, aYOut);
1752  cpl_free(wcs);
1753 
1754  return CPL_ERROR_NONE;
1755 } /* muse_wcs_pixel_from_projplane() */
1756 
1757 /*----------------------------------------------------------------------------*/
1776 /*----------------------------------------------------------------------------*/
1777 cpl_error_code
1778 muse_wcs_get_angles(cpl_propertylist *aHeader, double *aXAngle, double *aYAngle)
1779 {
1780  cpl_ensure_code(aHeader && aXAngle && aYAngle, CPL_ERROR_NULL_INPUT);
1781 
1782  cpl_errorstate prestate = cpl_errorstate_get();
1783  double cd11 = muse_pfits_get_cd(aHeader, 1, 1),
1784  cd22 = muse_pfits_get_cd(aHeader, 2, 2),
1785  cd12 = muse_pfits_get_cd(aHeader, 1, 2),
1786  cd21 = muse_pfits_get_cd(aHeader, 2, 1),
1787  det = cd11 * cd22 - cd12 * cd21;
1788  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
1789  if (det < 0.) {
1790  cd12 *= -1;
1791  cd11 *= -1;
1792  }
1793  if (cd12 == 0. && cd21 == 0.) { /* matrix without rotation */
1794  *aXAngle = 0.;
1795  *aYAngle = 0.;
1796  return CPL_ERROR_NONE;
1797  }
1798  /* angles in degrees */
1799  *aXAngle = atan2(cd12, cd11) * CPL_MATH_DEG_RAD;
1800  *aYAngle = atan2(-cd21, cd22) * CPL_MATH_DEG_RAD;
1801  return CPL_ERROR_NONE;
1802 } /* muse_wcs_get_angles() */
1803 
1804 /*----------------------------------------------------------------------------*/
1823 /*----------------------------------------------------------------------------*/
1824 cpl_error_code
1825 muse_wcs_get_scales(cpl_propertylist *aHeader, double *aXScale, double *aYScale)
1826 {
1827  cpl_ensure_code(aHeader && aXScale && aYScale, CPL_ERROR_NULL_INPUT);
1828 
1829  cpl_errorstate prestate = cpl_errorstate_get();
1830  double cd11 = muse_pfits_get_cd(aHeader, 1, 1),
1831  cd22 = muse_pfits_get_cd(aHeader, 2, 2),
1832  cd12 = muse_pfits_get_cd(aHeader, 1, 2),
1833  cd21 = muse_pfits_get_cd(aHeader, 2, 1),
1834  det = cd11 * cd22 - cd12 * cd21;
1835  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
1836 
1837  if (det < 0.) {
1838  cd12 *= -1;
1839  cd11 *= -1;
1840  }
1841  if (cd12 == 0. && cd21 == 0.) { /* matrix without rotation */
1842  *aXScale = cd11;
1843  *aYScale = cd22;
1844  return CPL_ERROR_NONE;
1845  }
1846  *aXScale = sqrt(cd11*cd11 + cd12*cd12); /* only the absolute value */
1847  *aYScale = sqrt(cd22*cd22 + cd21*cd21);
1848  return CPL_ERROR_NONE;
1849 } /* muse_wcs_get_scales() */
1850 
1851 /*----------------------------------------------------------------------------*/
1868 /*----------------------------------------------------------------------------*/
1869 muse_wcs *
1870 muse_wcs_new(cpl_propertylist *aHeader)
1871 {
1872  muse_wcs *wcs = cpl_calloc(1, sizeof(muse_wcs));
1873  if (!aHeader) {
1874  wcs->cd11 = wcs->cd22 = wcs->cddet = 1.; /* see below */
1875  return wcs;
1876  }
1877 
1878  cpl_errorstate prestate = cpl_errorstate_get();
1879  wcs->crpix1 = muse_pfits_get_crpix(aHeader, 1);
1880  wcs->crpix2 = muse_pfits_get_crpix(aHeader, 2);
1881  wcs->crval1 = muse_pfits_get_crval(aHeader, 1);
1882  wcs->crval2 = muse_pfits_get_crval(aHeader, 2);
1883  if (!cpl_errorstate_is_equal(prestate)) {
1884  /* all these headers default to 0.0 following the FITS *
1885  * Standard, so we can ignore any errors set here */
1886  cpl_errorstate_set(prestate);
1887  }
1888 
1889  prestate = cpl_errorstate_get();
1890  wcs->cd11 = muse_pfits_get_cd(aHeader, 1, 1);
1891  wcs->cd22 = muse_pfits_get_cd(aHeader, 2, 2);
1892  wcs->cd12 = muse_pfits_get_cd(aHeader, 1, 2);
1893  wcs->cd21 = muse_pfits_get_cd(aHeader, 2, 1);
1894  if (!cpl_errorstate_is_equal(prestate) &&
1895  wcs->cd11 == 0. && wcs->cd12 == 0. && wcs->cd21 == 0. && wcs->cd22 == 0.) {
1896  /* FITS Standard says to handle the CD matrix like the PC *
1897  * matrix in this case, with 1 for the diagonal elements */
1898  wcs->cd11 = wcs->cd22 = wcs->cddet = 1.;
1899  cpl_errorstate_set(prestate); /* not a real error */
1900  }
1901  wcs->cddet = wcs->cd11 * wcs->cd22 - wcs->cd12 * wcs->cd21;
1902  if (wcs->cddet == 0.) {
1903  cpl_error_set(__func__, CPL_ERROR_SINGULAR_MATRIX);
1904  }
1905 
1906  /* wcs->iscelsph defaults to 0 = CPL_FALSE, leave it at that */
1907  if (getenv("MUSE_DEBUG_WCS") && atoi(getenv("MUSE_DEBUG_WCS")) > 0) {
1908  cpl_msg_debug(__func__, "wcs: axis1 = { %f %f %e }, axis2 = { %f %f %e }, "
1909  "crossterms = { %e %e }, det = %e",
1910  wcs->crpix1, wcs->crval1, wcs->cd11,
1911  wcs->crpix2, wcs->crval2, wcs->cd22,
1912  wcs->cd12, wcs->cd21, wcs->cddet);
1913  }
1914  return wcs;
1915 } /* muse_wcs_new() */
1916 
double crpix2
Definition: muse_wcs.h:106
#define MUSE_PIXTABLE_XPOS
Definition: muse_pixtable.h:51
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
#define MUSE_HDR_PT_PREDAR_YLO
Structure definition for a collection of muse_images.
cpl_error_code muse_wcs_get_scales(cpl_propertylist *aHeader, double *aXScale, double *aYScale)
Compute the spatial scales (in degrees) from the FITS header WCS.
Definition: muse_wcs.c:1825
cpl_error_code muse_wcs_pixel_from_projplane(cpl_propertylist *aHeader, double aX, double aY, double *aXOut, double *aYOut)
Convert from projection plane coordinates to pixel coordinates.
Definition: muse_wcs.c:1739
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:85
double muse_pfits_get_cd(const cpl_propertylist *aHeaders, unsigned int aAxisI, unsigned int aAxisJ)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:446
#define MUSE_HDR_PT_XLO
FITS header keyword contains the lower limit of the data in x-direction.
double muse_pfits_get_crval(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:423
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:105
double muse_pfits_get_ra(const cpl_propertylist *aHeaders)
find out the right ascension
Definition: muse_pfits.c:267
double cd22
Definition: muse_wcs.h:108
const muse_cpltable_def muse_wcs_detections_def[]
Definition of the table structure for the astrometric field detections.
Definition: muse_wcs.c:116
cpl_image * data
the data extension
Definition: muse_image.h:46
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
#define MUSE_WCS_DETIMAGE_EXTNAME
name for special astrometry detection image
Definition: muse_wcs.h:56
muse_image * muse_combine_median_create(muse_imagelist *aImages)
Median combine a list of input images.
Definition: muse_combine.c:316
static void muse_wcs_pixel_from_projplane_fast(muse_wcs *aWCS, double aX, double aY, double *aXOut, double *aYOut)
Convert from projection plane coordinates to pixel coordinates.
Definition: muse_wcs.h:245
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aGrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
cpl_image * stat
the statistics extension
Definition: muse_image.h:64
cpl_error_code muse_wcs_projplane_from_celestial(cpl_propertylist *aHeader, double aRA, double aDEC, double *aX, double *aY)
Convert from celestial spherical coordinates to projection plane coordinates.
Definition: muse_wcs.c:1658
#define MUSE_HDR_PT_WCS_PROJ
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
void muse_imagelist_delete(muse_imagelist *aList)
Free the memory of the MUSE image list.
WCS object to store data needed while computing the astrometric solution.
Definition: muse_wcs.h:70
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
cpl_array * recnames
the reconstructed image filter names
Definition: muse_datacube.h:71
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
cpl_boolean muse_cplarray_has_duplicate(const cpl_array *aArray)
Check, if an array contains duplicate values.
cpl_error_code muse_pixtable_reset_dq(muse_pixtable *aPixtable, unsigned int aDQ)
Reset a given bad pixel status (DQ flag) for all pixels in the table.
cpl_error_code muse_cpltable_check(const cpl_table *aTable, const muse_cpltable_def *aDef)
Check whether the table contains the fields of the definition.
cpl_error_code muse_datacube_save(muse_datacube *aCube, const char *aFilename)
Save the three cube extensions and the FITS headers of a MUSE datacube to a file. ...
cpl_error_code muse_wcs_project_tan(muse_pixtable *aPixtable, const cpl_propertylist *aWCS)
Carry out a gnomonic projection of a pixel table into native spherical coordinates.
Definition: muse_wcs.c:1326
cpl_table * muse_wcs_centroid_stars(muse_image *aImage, float aSigma, muse_wcs_centroid_type aCentroid)
Detect and centroid stars on an image of an astrometric exposure.
Definition: muse_wcs.c:178
muse_resampling_crstats_type crtype
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
double muse_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
Definition: muse_pfits.c:401
const muse_cpltable_def muse_wcs_reference_def[]
Definition of the table structure for the astrometric reference sources.
Definition: muse_wcs.c:141
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
Definition: muse_wcs.c:1870
double muse_pfits_get_fwhm_end(const cpl_propertylist *aHeaders)
find out the ambient seeing at end of exposure (in arcsec)
Definition: muse_pfits.c:1097
muse_wcs_centroid_type
Type of centroiding algorithm to use.
Definition: muse_wcs.h:90
#define MUSE_HDR_PT_PREDAR_YHI
Structure definition of MUSE pixel table.
#define MUSE_WCS_KEYS
regular expression for WCS properties
Definition: muse_wcs.h:48
muse_image * muse_imagelist_get(muse_imagelist *aList, unsigned int aIdx)
Get the muse_image of given list index.
muse_pixtable_wcs muse_pixtable_wcs_check(muse_pixtable *aPixtable)
Check the state of the world coordinate system of a pixel table.
cpl_error_code muse_wcs_celestial_from_pixel(cpl_propertylist *aHeader, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
Definition: muse_wcs.c:1559
#define MUSE_HDR_PT_YLO
FITS header keyword contains the lower limit of the data in y-direction.
cpl_error_code muse_utils_fit_moffat_2d(const cpl_matrix *aPositions, const cpl_vector *aValues, const cpl_vector *aErrors, cpl_array *aParams, cpl_array *aPErrors, const cpl_array *aPFlags, double *aRMS, double *aRedChisq)
Fit a 2D Moffat function to a given set of data.
Definition: muse_utils.c:1878
static void muse_wcs_pixel_from_celestial_fast(muse_wcs *aWCS, double aRA, double aDEC, double *aX, double *aY)
Convert from celestial spherical coordinates to pixel coordinates.
Definition: muse_wcs.h:176
static void muse_wcs_celestial_from_pixel_fast(muse_wcs *aWCS, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
Definition: muse_wcs.h:138
double cddet
Definition: muse_wcs.h:109
muse_resampling_params * muse_resampling_params_new(muse_resampling_type aMethod)
Create the resampling parameters structure.
cpl_error_code muse_wcs_solve(muse_wcs_object *aWCSObj, cpl_table *aReference, float aRadius, float aFAccuracy, int aIter, float aThresh)
Find the world coordinate solution for a given field using coordinates of detected sources and coordi...
Definition: muse_wcs.c:671
muse_wcs_object * muse_wcs_object_new(void)
Allocate memory for a new muse_wcs_object object.
Definition: muse_wcs.c:71
#define MUSE_HDR_PT_PREDAR_XLO
cpl_error_code muse_wcs_projplane_from_pixel(cpl_propertylist *aHeader, double aX, double aY, double *aXOut, double *aYOut)
Convert from pixel coordinates to projection plane coordinates.
Definition: muse_wcs.c:1706
double muse_pfits_get_fwhm_start(const cpl_propertylist *aHeaders)
find out the ambient seeing at start of exposure (in arcsec)
Definition: muse_pfits.c:1078
muse_table * muse_table_load_filter(muse_processing *aProcessing, const char *aFilterName)
Load a table for a given filter name.
Definition: muse_utils.c:799
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:76
cpl_error_code muse_wcs_pixel_from_celestial(cpl_propertylist *aHeader, double aRA, double aDEC, double *aX, double *aY)
Convert from celestial spherical coordinates to pixel coordinates.
Definition: muse_wcs.c:1603
double muse_pfits_get_dec(const cpl_propertylist *aHeaders)
find out the declination
Definition: muse_pfits.c:285
#define MUSE_HDR_PT_PREDAR_XHI
Structure to store a table together with a property list.
Definition: muse_table.h:43
void muse_table_delete(muse_table *aTable)
Deallocate memory associated to a muse_table object.
Definition: muse_table.c:80
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:81
const char * muse_pfits_get_ctype(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis type string
Definition: muse_pfits.c:469
cpl_error_code muse_pixtable_save(muse_pixtable *aPixtable, const char *aFilename)
Save a MUSE pixel table to a file on disk.
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
double muse_astro_posangle(const cpl_propertylist *aHeader)
Derive the position angle of an observation from information in a FITS header.
Definition: muse_astro.c:413
cpl_error_code muse_wcs_optimize_solution(muse_wcs_object *aWCSObj, float aDetSigma, muse_wcs_centroid_type aCentroid, cpl_table *aReference, float aRadius, float aFAccuracy, int aIter, float aRejSigma)
Find an optimal astrometry solution given a detection image.
Definition: muse_wcs.c:988
#define MUSE_HDR_PT_WCS_POSI
cpl_error_code muse_image_save(muse_image *aImage, const char *aFilename)
Save the three image extensions and the FITS headers of a MUSE image to a file.
Definition: muse_image.c:405
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
static void muse_wcs_projplane_from_pixel_fast(muse_wcs *aWCS, double aX, double aY, double *aXOut, double *aYOut)
Convert from pixel coordinates to projection plane coordinates.
Definition: muse_wcs.h:216
#define MUSE_HDR_PT_DAR_NAME
cpl_propertylist * muse_wcs_apply_cd(const cpl_propertylist *aExpHeader, const cpl_propertylist *aCalHeader)
Apply the CDi_j matrix of an astrometric solution to an observation.
Definition: muse_wcs.c:1237
cpl_error_code muse_utils_image_get_centroid_window(cpl_image *aImage, int aX1, int aY1, int aX2, int aY2, double *aX, double *aY, muse_utils_centroid_type aBgType)
Compute centroid over an image window, optionally marginalizing over the background.
Definition: muse_utils.c:1234
Resampling parameters.
Definition of a cpl table structure.
void muse_resampling_params_delete(muse_resampling_params *aParams)
Delete a resampling parameters structure.
cpl_error_code muse_wcs_get_angles(cpl_propertylist *aHeader, double *aXAngle, double *aYAngle)
Compute the rotation angles (in degrees) from the FITS header WCS.
Definition: muse_wcs.c:1778
muse_image * muse_datacube_collapse(muse_datacube *aCube, const muse_table *aFilter)
Integrate a FITS NAXIS=3 datacube along the wavelength direction.
void muse_wcs_object_delete(muse_wcs_object *aWCSObj)
Deallocate memory associated to a muse_wcs_object.
Definition: muse_wcs.c:89
#define MUSE_HDR_PT_DAR_CORR
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:66
cpl_error_code muse_wcs_position_celestial(muse_pixtable *aPixtable, double aRA, double aDEC)
Convert native to celestial spherical coordinates in a pixel table.
Definition: muse_wcs.c:1466
double crval2
Definition: muse_wcs.h:107
#define MUSE_HDR_PT_WCS
#define MUSE_PIXTABLE_YPOS
Definition: muse_pixtable.h:52
#define MUSE_HDR_PT_YHI
FITS header keyword contains the upper limit of the data in y-direction.
cpl_error_code muse_imagelist_set(muse_imagelist *aList, muse_image *aImage, unsigned int aIdx)
Set the muse_image of given list index.
cpl_error_code muse_pixtable_compute_limits(muse_pixtable *aPixtable)
(Re-)Compute the limits of the coordinate columns of a pixel table.
int muse_quality_image_reject_using_dq(cpl_image *aData, cpl_image *aDQ, cpl_image *aStat)
Reject pixels of one or two images on a DQ image.
Definition: muse_quality.c:628
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1352
muse_imagelist * recimages
the reconstructed image data
Definition: muse_datacube.h:64
cpl_error_code muse_wcs_locate_sources(muse_pixtable *aPixtable, float aSigma, muse_wcs_centroid_type aCentroid, muse_wcs_object *aWCSObj)
Determine the centers of stars on an astrometric exposure.
Definition: muse_wcs.c:496
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:86
cpl_propertylist * header
The FITS header.
#define MUSE_HDR_PT_XHI
FITS header keyword contains the upper limit of the data in x-ion.
cpl_propertylist * muse_wcs_create_default(void)
Create FITS headers containing a default (relative) WCS.
Definition: muse_wcs.c:1174